EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
Lorenz96_hidden_empire.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 ! Implements Lorenz 1996 model with hidden fast modes as in
3 ! Bergemann 2010 http://doi.wiley.com/10.1002/qj.672
4 !
5 !The MIT License (MIT)
6 !
7 !Copyright (c) 2015 Philip A. Browne
8 !
9 !Permission is hereby granted, free of charge, to any person obtaining a copy
10 !of this software and associated documentation files (the "Software"), to deal
11 !in the Software without restriction, including without limitation the rights
12 !to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
13 !copies of the Software, and to permit persons to whom the Software is
14 !furnished to do so, subject to the following conditions:
15 !
16 !The above copyright notice and this permission notice shall be included in all
17 !copies or substantial portions of the Software.
18 !
19 !THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
20 !IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
21 !FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
22 !AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
23 !LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
24 !OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
25 !SOFTWARE.
26 !
27 !Email: p.browne@reading.ac.uk
28 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
29 
31  implicit none
32  include 'mpif.h'
33  real(kind=kind(1.0D0)) :: dt=2.5d-3
34  real(kind=kind(1.0D0)), allocatable, dimension(:,:) :: x,k1,k2,k3,k4
35  real(kind=kind(1.0d0)) :: F=8.0d0
36  real(kind=kind(1.0d0)) :: alpha=0.5d0
37  real(kind=kind(1.0d0)) :: delta=0.5d0
38  real(kind=kind(1.0d0)) :: epsilon=2.5d-3
39  real(kind=kind(1.0d0)) :: gamma=0.1d0
40  integer :: N=40
41  integer :: total_timesteps=1000
42  integer :: t
43  integer :: mpi_err,mdl_id,cpl_root,cpl_mpi_comm
44  integer, dimension(MPI_STATUS_SIZE) :: mpi_status
45  integer :: ios
46  logical :: l96_exists
47 
48  namelist/l96/n,&
49  &total_timesteps,&
50  &f,&
51  &dt,&
52  &delta,&
53  &epsilon,&
54  &alpha,&
55  &gamma
56 
57  inquire(file='l96.nml', exist=l96_exists)
58  if(l96_exists) then
59  open(32,file='l96.nml',iostat=ios,action='read'&
60  &,status='old')
61  if(ios .ne. 0) stop 'Cannot open l96.nml'
62  read(32,nml=l96)
63  close(32)
64  end if
65 
66 
67 
68 
69  call initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
70 
71 
72  allocate(x(n,3),k1(n,3),k2(n,3),k3(n,3),k4(n,3))
73  x = f
74  x(:,1) = f
75  x(:,2) = f
76  x(:,3) = 0.
77  x(n/2,1) = f+0.05d0
78  print*,x(:,1)
79 
80  if(mdl_id .eq. 0) then
81  call mpi_send(x(:,1),n,mpi_double_precision,cpl_root&
82  &,1,cpl_mpi_comm,mpi_err)
83  call mpi_recv(x(:,1),n,mpi_double_precision,cpl_root&
84  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
85  end if
86 2 continue
87 
88 
89  do t = 1,total_timesteps
90  k1 = g(x , n , f ,alpha,delta,epsilon,gamma)
91  k2 = g(x +0.5d0 * dt * k1 , n , f ,alpha,delta,epsilon,gamma)
92  k3 = g(x +0.5d0 * dt * k2 , n , f ,alpha,delta,epsilon,gamma)
93  k4 = g(x + dt * k3 , n , f ,alpha,delta,epsilon,gamma)
94  x = x + dt *( k1 + 2.0d0 *( k2 + k3 ) + k4 )/6.0d0
95 
96  if(mdl_id .eq. 0) then
97  call mpi_send(x(:,1),n,mpi_double_precision,cpl_root&
98  &,1,cpl_mpi_comm,mpi_err)
99  call mpi_recv(x(:,1),n,mpi_double_precision,cpl_root&
100  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
101  end if
102  print*,x(:,1)
103 
104  select case(mpi_status(mpi_tag))
105  case(2)
106  go to 2
107  case(3)
108  go to 3
109  case default
110  end select
111 
112  end do
113 3 continue
114 ! call mpi_finalize(mpi_err)
115 
116 contains
117  function g (X , N, F ,alpha,delta,epsilon,gamma)
118  implicit none
119  real(kind=kind(1.0D0)),intent(in),dimension(N,3) :: X
120  integer, intent(in) :: N
121  real(kind=kind(1.0D0)),dimension(N,3) :: g
122  real(kind=kind(1.0D0)),intent(in) :: F
123  real(kind=kind(1.0D0)),intent(in) :: alpha
124  real(kind=kind(1.0D0)),intent(in) :: delta
125  real(kind=kind(1.0D0)),intent(in) :: epsilon
126  real(kind=kind(1.0D0)),intent(in) :: gamma
127  integer :: j
128 
129  g(1,1) = (1.0d0-delta)*(x(2,1)-x(n-1,1) )*x(n,1) + delta*(x(n,1)*x(2,2)-x(n-1,1)*x(n,2)) - x(1,1) + f
130  g(1,2) = x(1,3)
131  g(1,3) = (-x(1,2) + alpha**2*(x(2,2)-2*x(1,2)+x(n,2)) + x(1,1) - gamma*epsilon**2*x(1,3))/(epsilon**2)
132 
133  g(2,1) = (1.0d0-delta)*(x(3,1)-x(n,1) )*x(1,1) + delta*(x(1,1)*x(3,2)-x(n,1)*x(1,2)) - x(2,1) + f
134  g(2,2) = x(2,3)
135  g(2,3) = (-x(2,2) + alpha**2*(x(3,2)-2*x(2,2)+x(1,2)) + x(2,1) - gamma*epsilon**2*x(2,3))/(epsilon**2)
136 
137 
138  do j = 3,n-1
139  g(j,1) = (1.0d0-delta)*(x(j+1,1)-x(j-2,1))*x(j-1,1) + delta*(x(j-1,1)*x(j+1,2)-x(j-2,1)*x(j-1,2)) - x(j,1) + f
140  g(j,2) = x(j,3)
141  g(j,3) = (-x(j,2) + alpha**2*(x(j+1,2)-2*x(j,2)+x(j-1,2)) + x(j,1) - gamma*epsilon**2*x(j,3))/(epsilon**2)
142  end do
143 
144  g(n,1) = (1.0d0-delta)*(x(1,1)-x(n-2,1) )*x(n-1,1) + delta*(x(n-1,1)*x(1,2)-x(n-2,1)*x(n-1,2)) - x(n,1) + f
145  g(n,2) = x(n,3)
146  g(n,3) = (-x(n,2) + alpha**2*(x(1,2)-2*x(n,2)+x(n-1,2)) + x(n,1) - gamma*epsilon**2*x(n,3))/(epsilon**2)
147 
148 
149 
150  end function g
151  subroutine initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
152  implicit none
153  include 'mpif.h'
154  integer, intent(out) :: mdl_id,cpl_root,cpl_mpi_comm
155  integer :: mdl_num_proc=1
156  integer :: mpi_err,world_size,world_id
157  integer :: cpl_colour
158  integer :: particle_id,nens, da, nda
159  integer :: mdl_mpi_comm,mdlcolour
160  integer :: tmp_mdls_comm,models_id,models_size
161  call mpi_init(mpi_err)
162  da = 0
163  call mpi_comm_rank(mpi_comm_world,world_id,mpi_err)
164  call mpi_comm_size(mpi_comm_world,world_size,mpi_err)
165  call mpi_comm_split(mpi_comm_world,da,world_id,tmp_mdls_comm,mpi_err)
166  call mpi_comm_size(tmp_mdls_comm,models_size,mpi_err)
167  call mpi_comm_rank(tmp_mdls_comm,models_id, mpi_err)
168  mdlcolour = models_id/mdl_num_proc
169  call mpi_comm_split(tmp_mdls_comm,mdlcolour,models_id,mdl_mpi_comm,mpi_err)
170  call mpi_comm_rank(mdl_mpi_comm,mdl_id,mpi_err)
171  if(mdl_id .eq. 0) then
172  cpl_colour = 9999
173  else
174  cpl_colour = mpi_undefined
175  end if
176  call mpi_comm_split(mpi_comm_world,cpl_colour,mdlcolour,cpl_mpi_comm,mpi_err)
177  if(mdl_id .eq. 0) then
178  call mpi_comm_size(cpl_mpi_comm,nens,mpi_err)
179  call mpi_comm_rank(cpl_mpi_comm,particle_id,mpi_err)
180  nda = world_size-models_size;nens = nens - nda
181  cpl_root = ((nda*particle_id)/nens)+nens
182  else
183  cpl_root = -1
184  end if
185  end subroutine initialise_mpi
186 end program lorenz96_hidden
subroutine initialise_mpi(mdl_id, cpl_root, cpl_mpi_comm)
real(kind=kind(1.0d0)) function, dimension(n, 3) g(X, N, F, alpha, delta, epsilon, gamma)
program lorenz96_hidden