EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
minimal_model_comms_v3.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 ! minimal_model_comms just sets up communicators for empire for
3 ! testing purposes
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 
33  include 'mpif.h'
34  integer :: mpi_err
35  integer :: world_rank
36  integer :: world_size
37  integer, parameter :: mdl_num_proc=5
38  integer :: mdl_mpi_comm
39  integer :: temp_mdls_comm
40  integer :: temp_mdls_rank
41  integer :: mdl_rank
42  integer :: da
43  integer, parameter :: rk = kind(1.0d0)
44  real(kind=rk), allocatable, dimension(:) :: state_vector
45  integer :: state_dim
46  integer :: i
47  integer :: cpl_root
48  integer :: particle_rank
49  integer :: nda
50  integer :: nens
51  integer :: temp_mdls_size
52  integer :: temp_cpl_comm
53  integer :: first_ptcl
54  integer :: final_ptcl
55  integer :: cpl_mpi_comm
56  integer :: null_mpi_comm
57  integer :: tmp_colour_2
58  integer :: tmp_cpl_comm2
59  integer :: tmp_cpl_rank
60 
61  print*,'RUNNING MINIMAL_MODEL_COMMS_V3'
62 
63  call mpi_init(mpi_err)
64  if(mpi_err .eq. 0) then
65  print*,'mpi_init successful'
66  else
67  print*,'mpi_init unsuccessful'
68  end if
69 
70 
71  call mpi_comm_rank(mpi_comm_world,world_rank,mpi_err)
72  if(mpi_err .eq. 0) then
73  print*,'mpi_comm_rank successful'
74  print*,'world_rank = ',world_rank
75  else
76  print*,'mpi_comm_rank unsuccessful'
77  end if
78 
79 
80  call mpi_comm_size(mpi_comm_world,world_size,mpi_err)
81  if(mpi_err .eq. 0) then
82  print*,'mpi_comm_size successful'
83  print*,'world_size = ',world_size
84  else
85  print*,'mpi_comm_size unsuccessful'
86  end if
87 
88 
89  cpl_root = world_size-1
90  print*,'rank = ',world_rank,' on mpi_comm_world which has size ',world_size
91  da = 0
92 
93 
94  call mpi_allreduce(mdl_num_proc,i,1,mpi_integer,mpi_max&
95  &,mpi_comm_world,mpi_err)
96  if(mpi_err .eq. 0) then
97  print*,'mpi_allreduce successful'
98  print*,'i = ',i
99  else
100  print*,'mpi_allreduce unsuccessful'
101  end if
102 !!!==================================================!!!
103 !!!================ BRUCE LEE =======================!!!
104 !!!==================================================!!!
105 
106 
107  call mpi_comm_split(mpi_comm_world,da,world_rank,temp_mdls_comm&
108  &,mpi_err)
109  if(mpi_err .eq. 0) then
110  print*,'mpi_comm_split successful: temp_mdls_comm created'
111  else
112  print*,'mpi_comm_split unsuccessful: temp_mdls_comm not created'
113  end if
114 
115  call mpi_comm_size(temp_mdls_comm,temp_mdls_size,mpi_err)
116  if(mpi_err .eq. 0) then
117  print*,'mpi_comm_size successful'
118  print*,'temp_mdls_size = ',temp_mdls_size
119  else
120  print*,'mpi_comm_size unsuccessful'
121  end if
122 
123 
124  if(mod(temp_mdls_size,mdl_num_proc) .ne. 0) then
125  print*,'MINIMAL MODEL LAUNCH ERROR.'
126  print*,'MUST LAUNCH A MULTIPLE OF ',mdl_num_proc,' copies of the &
127  &model'
128  stop
129  end if
130 
131 
132  nda = (world_size-temp_mdls_size)/mdl_num_proc
133  if(nda .lt. 1) then
134  print*,'MINIMAL MODEL COMMS v3 ERROR: nda is less than 1.'
135  print*,'Make sure you launch with a DA CODE'
136  stop
137  end if
138 
139 !!!==================================================!!!
140 !!!================ JEAN CLAUDE VAN DAMME ===========!!!
141 !!!==================================================!!!
142 
143 
144  nens = temp_mdls_size/mdl_num_proc
145  call mpi_comm_rank(temp_mdls_comm,temp_mdls_rank,mpi_err)
146  if(mpi_err .eq. 0) then
147  print*,'mpi_comm_rank successful'
148  print*,'temp_mdls_rank = ',temp_mdls_rank
149  else
150  print*,'mpi_comm_rank unsuccessful'
151  end if
152 
153 
154  particle_rank = temp_mdls_rank/mdl_num_proc
155 
156  call mpi_comm_split(temp_mdls_comm,particle_rank,temp_mdls_rank&
157  &,mdl_mpi_comm,mpi_err)
158  if(mpi_err .eq. 0) then
159  print*,'mpi_comm_split successful: mdl_mpi_comm created'
160  else
161  print*,'mpi_comm_split unsuccessful: mdl_mpi_comm not created'
162  end if
163 
164 
165 
166  call mpi_comm_rank(mdl_mpi_comm,mdl_rank,mpi_err)
167  if(mpi_err .eq. 0) then
168  print*,'mpi_comm_rank successful'
169  print*,'mdl_rank = ',mdl_rank
170  else
171  print*,'mpi_comm_rank unsuccessful'
172  end if
173 
174 !!!==================================================!!!
175 !!!============== STEVEN SEGAL ======================!!!
176 !!!==================================================!!!
177 
178 
179 
180  cpl_root = nda*particle_rank/nens
181  print*,'cpl_root = ',cpl_root
182 
183  if(cpl_root .lt. 0) then
184  print*,'MINIMAL MODEL COMMS v3 ERROR: cpl_root is less than 0.'
185  print*,'Make sure you launch with a DA CODE'
186  stop
187  end if
188 
189  call mpi_comm_split(mpi_comm_world,cpl_root,temp_mdls_rank,temp_cpl_comm,mpi_err)
190  if(mpi_err .eq. 0) then
191  print*,'mpi_comm_split successful: temp_cpl_comm created'
192  else
193  print*,'mpi_comm_split unsuccessful: temp_cpl_comm not created'
194  end if
195 
196 !!!==================================================!!!
197 !!!============== CHUCK NORRIS ======================!!!
198 !!!==================================================!!!
199 
200  call mpi_comm_rank(temp_cpl_comm,tmp_cpl_rank,mpi_err)
201  if(mpi_err .eq. 0) then
202  print*,'mpi_comm_rank successful: tmp_cpl_rank = ',tmp_cpl_rank
203  else
204  print*,'mpl_comm_rank unsuccessful: tmp_cpl_rank not detected'
205  end if
206 
207  tmp_colour_2 = mod(tmp_cpl_rank,mdl_num_proc)
208  print*,'tmp_colour_2 = ',tmp_colour_2
209  call mpi_comm_split(temp_cpl_comm,tmp_colour_2,tmp_cpl_rank&
210  &,tmp_cpl_comm2,mpi_err)
211  if(mpi_err .eq. 0) then
212  print*,'mpi_comm_split successful: tmp_cpl_comm2 created'
213  else
214  print*,'mpi_comm_split unsuccessful: tmp_cpl_comm2 not created'
215  end if
216 
217 !!!==================================================!!!
218 !!!=============== JACKIE CHAN ======================!!!
219 !!!==================================================!!!
220 
221 
222  first_ptcl = -((-cpl_root)*nens/nda)
223  final_ptcl = -((-cpl_root-1)*nens/nda)-1
224 
225  first_ptcl = ceiling(real(cpl_root)*real(nens)/real(nda))
226  final_ptcl = ceiling(real(cpl_root+1)*real(nens)/real(nda))-1
227 
228 
229  print*,'range of particles = ',first_ptcl,final_ptcl
230 
231 
232 
233  do i = first_ptcl,final_ptcl
234  print*,'i = ',i,' particle_rank = ',particle_rank
235  if(i .eq. particle_rank) then
236  call mpi_comm_split(tmp_cpl_comm2,1,temp_mdls_rank&
237  &,cpl_mpi_comm,mpi_err)
238  print*,'created cpl_mpi_comm',cpl_mpi_comm
239  else
240  print*,'doing null splitting'
241  call mpi_comm_split(tmp_cpl_comm2,0,temp_mdls_rank&
242  &,null_mpi_comm,mpi_err)
243  print*,'created mpi_comm_null'
244  call mpi_comm_free(null_mpi_comm,mpi_err)
245  print*,'freed up null_mpi_comm'
246  end if
247 
248 
249  end do
250 
251  cpl_root = 1
252 
253 
254 !!!==================================================!!!
255 !!!=============== MICHELLE YEOH ====================!!!
256 !!!==================================================!!!
257 
258  !free up the temporary communicators
259  call mpi_comm_free(temp_cpl_comm ,mpi_err)
260  call mpi_comm_free(tmp_cpl_comm2,mpi_err)
261 
262 
263 
264  select case(mdl_rank)
265  case(0)
266  state_dim = 1
267  case(1)
268  state_dim = 3
269  case(2)
270  state_dim = 2
271  case(3)
272  state_dim = 5
273  case(4)
274  state_dim = 1
275  case default
276  print*,'it was at this point, model realised, he fucked up'
277  stop
278  end select
279  allocate(state_vector(state_dim))
280 
281 
282 
283 
284  state_vector = 10*mdl_rank + (/ (real(i,rk), i = 1,state_dim) /)
285 
286  print*,'state_vector = '
287  print*,state_vector
288 
289  print*,'doing a send on cpl_mpi_comm'
290 
291  call mpi_send(state_dim,1,mpi_integer,cpl_root,1,cpl_mpi_comm,mpi_err)
292  print*,'finished the send on cpl_mpi_comm'
293 !!!==================================================!!!
294 !!!=============== SCARLETT JOHANSSON ===============!!!
295 !!!==================================================!!!
296 
297 
298 
299  call mpi_finalize(mpi_err)
300  print*,'MINIMAL_MODEL_COMMS_v3 finished nicely'
301 end program minimal_model_comms_v3
program minimal_model_comms_v3