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