EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
Lorenz96_empire_v2.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 ! Lorenz96_empire.f90 Implements Lorenz 1996 model with EMPIRE coupling
3 !
4 !The MIT License (MIT)
5 !
6 !Copyright (c) 2015 Philip A. Browne
7 !
8 !Permission is hereby granted, free of charge, to any person obtaining a copy
9 !of this software and associated documentation files (the "Software"), to deal
10 !in the Software without restriction, including without limitation the rights
11 !to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
12 !copies of the Software, and to permit persons to whom the Software is
13 !furnished to do so, subject to the following conditions:
14 !
15 !The above copyright notice and this permission notice shall be included in all
16 !copies or substantial portions of the Software.
17 !
18 !THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
19 !IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
20 !FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
21 !AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
22 !LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
23 !OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
24 !SOFTWARE.
25 !
26 !Email: p.browne@reading.ac.uk
27 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
28 
29 program lorenz96_v2
30  implicit none
31  include 'mpif.h'
32  real(kind=kind(1.0D0)) :: dt=1.0d-2
33  real(kind=kind(1.0D0)), allocatable, dimension(:) :: x,k1,k2,k3,k4
34  real(kind=kind(1.0d0)) :: F=8.0d0
35  integer :: N=40
36  integer :: total_timesteps=100
37  integer :: t
38  integer :: mpi_err,mdl_rank,cpl_root,cpl_mpi_comm
39  real(kind=kind(1.0d0)), dimension(0) :: send_null
40  integer :: tag
41  integer :: ios
42  logical :: l96_exists
43 
44  namelist/l96/n,&
45  &total_timesteps,&
46  &f,&
47  &dt
48 
49  inquire(file='l96.nml', exist=l96_exists)
50  if(l96_exists) then
51  open(32,file='l96.nml',iostat=ios,action='read'&
52  &,status='old')
53  if(ios .ne. 0) stop 'Cannot open l96.nml'
54  read(32,nml=l96)
55  close(32)
56  end if
57 
58 
59  call initialise_mpi_v2(mdl_rank,cpl_root,cpl_mpi_comm)
60 
61  call empire_process_dimensions(n,cpl_root,cpl_mpi_comm)
62  allocate(x(n),k1(n),k2(n),k3(n),k4(n))
63  x = f
64  x(n/2) = f+0.05
65  print*,x
66 
67 
68  !send the state to da code with mpi_gatherv
69  call mpi_gatherv(x,n,mpi_double_precision,x&
70  &,n,n,mpi_double_precision,cpl_root&
71  &,cpl_mpi_comm,mpi_err)
72 
73 
74  !get the state back from da code with mpi_gatherv
75  call mpi_scatterv(send_null,0,0,mpi_double_precision,x&
76  &,n,mpi_double_precision,cpl_root,cpl_mpi_comm,mpi_err)
77  !get the tag from the da code
78  call mpi_bcast(tag,1,mpi_integer,cpl_root,cpl_mpi_comm,mpi_err)
79 2 continue
80 
81 
82  do t = 1,total_timesteps
83  k1 = g(x , n , f )
84  k2 = g(x +0.5d0 * dt * k1 , n , f )
85  k3 = g(x +0.5d0 * dt * k2 , n , f )
86  k4 = g(x + dt * k3 , n , f )
87  x = x + dt *( k1 + 2.0d0 *( k2 + k3 ) + k4 )/6.0d0
88 
89 
90  !send the state to da code with mpi_gatherv
91  call mpi_gatherv(x,n,mpi_double_precision,x&
92  &,n,n,mpi_double_precision,cpl_root&
93  &,cpl_mpi_comm,mpi_err)
94 
95  !get the state back from da code with mpi_gatherv
96  call mpi_scatterv(send_null,0,0,mpi_double_precision,x&
97  &,n,mpi_double_precision,cpl_root,cpl_mpi_comm,mpi_err)
98  !get the tag from the da code
99  call mpi_bcast(tag,1,mpi_integer,cpl_root,cpl_mpi_comm,mpi_err)
100 
101  print*,x
102 
103  select case(tag)
104  case(2)
105  go to 2
106  case(3)
107  go to 3
108  case default
109  end select
110 
111 
112 
113  end do
114 3 continue
115  call mpi_finalize(mpi_err)
116 
117 contains
118  function g (x , N, F )
119  implicit none
120  real(kind=kind(1.0D0)),intent(in),dimension(0:N-1) :: x
121  integer, intent(in) :: N
122  real(kind=kind(1.0D0)),dimension(0:N-1) :: g
123  real(kind=kind(1.0D0)),intent(in) :: F
124  integer :: j
125  do j = 0,n-1
126  g(j) = ( x(modulo(j+1,n)) -x( modulo(j-2,n)) )*&
127  &x(modulo(j-1,n)) - x(j) + f
128  end do
129  end function g
130  subroutine initialise_mpi_v2(mdl_rank,cpl_root,cpl_mpi_comm)
131 
132 
133  implicit none
134 
135  include 'mpif.h'
136  integer, intent(out) :: mdl_rank
137  integer, intent(out) :: cpl_root
138  integer, intent(out) :: cpl_mpi_comm
139 
140  integer, parameter :: mdl_num_proc=1
141  integer :: mdl_mpi_comm
142 
143  integer :: mpi_err
144  integer :: world_rank
145  integer :: world_size
146 
147  integer :: temp_mdls_size
148  integer :: temp_cpl_comm
149  integer :: temp_mdls_comm
150  integer :: temp_mdls_rank
151  integer :: da
152  integer :: i
153  integer :: particle_rank
154  integer :: nda
155  integer :: nens
156  integer :: first_ptcl
157  integer :: final_ptcl
158  integer :: null_mpi_comm
159 
160  logical :: msg=.true.
161 
162 
163  call mpi_init(mpi_err)
164  if(mpi_err .eq. 0) then
165  if(msg) print*,'mpi_init successful'
166  else
167  print*,'mpi_init unsuccessful'
168  end if
169 
170 
171  call mpi_comm_rank(mpi_comm_world,world_rank,mpi_err)
172  if(mpi_err .eq. 0) then
173  if(msg) print*,'mpi_comm_rank successful'
174  if(msg) print*,'world_rank = ',world_rank
175  else
176  print*,'mpi_comm_rank unsuccessful'
177  end if
178 
179 
180  call mpi_comm_size(mpi_comm_world,world_size,mpi_err)
181  if(mpi_err .eq. 0) then
182  if(msg) print*,'mpi_comm_size successful'
183  if(msg) print*,'world_size = ',world_size
184  else
185  print*,'mpi_comm_size unsuccessful'
186  end if
187 
188 
189  cpl_root = world_size-1
190  if(msg) then
191  print*,'rank = ',world_rank,' on mpi_comm_world which has &
192  &size ',world_size
193  end if
194 
195  da = 0
196 
197 
198  call mpi_allreduce(mdl_num_proc,i,1,mpi_integer,mpi_max&
199  &,mpi_comm_world,mpi_err)
200  if(mpi_err .eq. 0) then
201  if(msg) print*,'mpi_allreduce successful'
202  if(msg) print*,'i = ',i
203  else
204  print*,'mpi_allreduce unsuccessful'
205  end if
206 
207 
208 
209  call mpi_comm_split(mpi_comm_world,da,world_rank,temp_mdls_comm&
210  &,mpi_err)
211  if(mpi_err .eq. 0) then
212  if(msg) print*,'mpi_comm_split successful: temp_mdls_comm created'
213  else
214  print*,'mpi_comm_split unsuccessful: temp_mdls_comm not created'
215  end if
216 
217  call mpi_comm_size(temp_mdls_comm,temp_mdls_size,mpi_err)
218  if(mpi_err .eq. 0) then
219  if(msg) print*,'mpi_comm_size successful'
220  if(msg) print*,'temp_mdls_size = ',temp_mdls_size
221  else
222  print*,'mpi_comm_size unsuccessful'
223  end if
224 
225 
226  if(mod(temp_mdls_size,mdl_num_proc) .ne. 0) then
227  print*,'MINIMAL MODEL LAUNCH ERROR.'
228  print*,'MUST LAUNCH A MULTIPLE OF ',mdl_num_proc,' copies of the &
229  &model'
230  stop
231  end if
232 
233 
234  nda = world_size-temp_mdls_size
235  if(nda .lt. 1) then
236  print*,'MINIMAL MODEL COMMS v2 ERROR: nda is less than 1.'
237  print*,'Make sure you launch with a DA CODE'
238  stop
239  end if
240 
241 
242 
243  nens = temp_mdls_size/mdl_num_proc
244  call mpi_comm_rank(temp_mdls_comm,temp_mdls_rank,mpi_err)
245  if(mpi_err .eq. 0) then
246  if(msg) print*,'mpi_comm_rank successful'
247  if(msg) print*,'temp_mdls_rank = ',temp_mdls_rank
248  else
249  print*,'mpi_comm_rank unsuccessful'
250  end if
251 
252 
253  particle_rank = temp_mdls_rank/mdl_num_proc
254 
255  call mpi_comm_split(temp_mdls_comm,particle_rank,temp_mdls_rank&
256  &,mdl_mpi_comm,mpi_err)
257  if(mpi_err .eq. 0) then
258  if(msg) print*,'mpi_comm_split successful: mdl_mpi_comm created'
259  else
260  print*,'mpi_comm_split unsuccessful: mdl_mpi_comm not created'
261  end if
262 
263 
264 
265  call mpi_comm_rank(mdl_mpi_comm,mdl_rank,mpi_err)
266  if(mpi_err .eq. 0) then
267  if(msg) print*,'mpi_comm_rank successful'
268  if(msg) print*,'mdl_rank = ',mdl_rank
269  else
270  print*,'mpi_comm_rank unsuccessful'
271  end if
272 
273 
274  cpl_root = nda*particle_rank/nens
275  if(msg) print*,'cpl_root = ',cpl_root
276 
277  if(cpl_root .lt. 0) then
278  print*,'MINIMAL MODEL COMMS v2 ERROR: cpl_root is less than 0.'
279  print*,'Make sure you launch with a DA CODE'
280  stop
281  end if
282 
283  call mpi_comm_split(mpi_comm_world,cpl_root,temp_mdls_rank,temp_cpl_comm,mpi_err)
284  if(mpi_err .eq. 0) then
285  if(msg) print*,'mpi_comm_split successful: temp_cpl_comm created'
286  else
287  print*,'mpi_comm_split unsuccessful: temp_cpl_comm not created'
288  end if
289 
290 
291 
292 
293 
294 
295  first_ptcl = ceiling(real(cpl_root)*real(nens)/real(nda))
296  final_ptcl = ceiling(real(cpl_root+1)*real(nens)/real(nda))-1
297 
298 
299  if(msg) print*,'range of particles = ',first_ptcl,final_ptcl
300 
301 
302 
303  do i = first_ptcl,final_ptcl
304  if(msg) print*,'i = ',i,' particle_rank = ',particle_rank
305  if(i .eq. particle_rank) then
306  call mpi_comm_split(temp_cpl_comm,1,temp_mdls_rank&
307  &,cpl_mpi_comm,mpi_err)
308  if(msg) print*,'created cpl_mpi_comm'
309  else
310  if(msg) print*,'doing null splitting'
311  call mpi_comm_split(temp_cpl_comm,0,temp_mdls_rank&
312  &,null_mpi_comm,mpi_err)
313  if(msg) print*,'created mpi_comm_null'
314  call mpi_comm_free(null_mpi_comm,mpi_err)
315  if(msg) print*,'freed up null_mpi_comm'
316  end if
317 
318 
319  end do
320 
321  cpl_root = mdl_num_proc
322 
323 
324 
325 
326  end subroutine initialise_mpi_v2
327 
328  subroutine empire_process_dimensions(N,cpl_root,cpl_mpi_comm)
329  implicit none
330  include 'mpif.h'
331  integer, intent(in) :: N
332  integer, intent(in) :: cpl_root
333  integer, intent(in) :: cpl_mpi_comm
334  integer :: mpi_err
335  logical :: msg = .true.
336 
337  if(msg) print*,'called empire_process_dimensions'
338  call mpi_gather(n,1,mpi_integer,n&
339  &,1,mpi_integer,cpl_root,cpl_mpi_comm,mpi_err)
340  if(msg) print*,'finished the gather on cpl_mpi_comm for empire_process_dimensions'
341  end subroutine empire_process_dimensions
342 
343 
344 
345 end program lorenz96_v2
subroutine empire_process_dimensions(N, cpl_root, cpl_mpi_comm)
real(kind=kind(1.0d0)) function, dimension(n, 3) g(X, N, F, alpha, delta, epsilon, gamma)
subroutine initialise_mpi_v2(mdl_rank, cpl_root, cpl_mpi_comm)
program lorenz96_v2