EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
minimal_model.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 ! minimal_model_comms sets up communicators for empire for
3 ! testing purposes and sends the state back and forward
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  integer :: mpi_err,mdl_id,cpl_root,cpl_mpi_comm
34  integer :: total_timesteps,i
35  integer, dimension(MPI_STATUS_SIZE) :: mpi_status
36  real(kind=kind(1.0d0)), dimension(3) :: x
37  call initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
38 
39  print*,'Reading total_timesteps from file timesteps: '
40  open(11,file='timesteps',action='read',status='old')
41  read(11,*) total_timesteps
42  close(11)
43 
44 
45 
46  !send the state to da code with mpi_send
47  if(mdl_id .eq. 0) then
48  call mpi_send(x,3,mpi_double_precision,cpl_root&
49  &,1,cpl_mpi_comm,mpi_err)
50 
51  !get the state back from da code with mpi_recv
52  call mpi_recv(x,3,mpi_double_precision,cpl_root&
53  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
54  print*,'Received tag = ',mpi_status(mpi_tag)
55  end if
56 
57  do i = 1,total_timesteps
58  print*,'Timestep = ',i
59 
60  if(mdl_id .eq. 0) then
61  !send the state to da code with mpi_send
62  call mpi_send(x,3,mpi_double_precision,cpl_root&
63  &,1,cpl_mpi_comm,mpi_err)
64 
65  !get the state back from da code with mpi_recv
66  call mpi_recv(x,3,mpi_double_precision,cpl_root&
67  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
68  print*,'Received tag = ',mpi_status(mpi_tag)
69  end if
70  end do
71 
72 
73 
74 
75 
76  call mpi_finalize(mpi_err)
77 contains
78  subroutine initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
79  implicit none
80  include 'mpif.h'
81  integer, intent(out) :: mdl_id,cpl_root,cpl_mpi_comm
82  integer :: mdl_num_proc=1
83  integer :: mpi_err,world_size,world_id
84  integer :: cpl_colour
85  integer :: particle_id,nens, da, nda
86  integer :: mdl_mpi_comm,mdlcolour
87  integer :: tmp_mdls_comm,models_id,models_size
88  call mpi_init(mpi_err)
89  if(mpi_err .eq. 0) then
90  print*,'mpi_init successful'
91  else
92  print*,'mpi_init unsuccessful'
93  end if
94 
95  da = 0
96  call mpi_comm_rank(mpi_comm_world,world_id,mpi_err)
97  if(mpi_err .eq. 0) then
98  print*,'mpi_comm_rank successful'
99  print*,'world_id = ',world_id
100  else
101  print*,'mpi_comm_rank unsuccessful'
102  end if
103 
104 
105  call mpi_comm_size(mpi_comm_world,world_size,mpi_err)
106  if(mpi_err .eq. 0) then
107  print*,'mpi_comm_size successful'
108  print*,'world_size = ',world_size
109  else
110  print*,'mpi_comm_size unsuccessful'
111  end if
112 
113 
114  call mpi_comm_split(mpi_comm_world,da,world_id,tmp_mdls_comm,mpi_err)
115  if(mpi_err .eq. 0) then
116  print*,'mpi_comm_split successful'
117  else
118  print*,'mpi_comm_split unsuccessful'
119  end if
120 
121 
122  call mpi_comm_size(tmp_mdls_comm,models_size,mpi_err)
123  if(mpi_err .eq. 0) then
124  print*,'mpi_comm_size successful'
125  print*,'models_size = ',models_size
126  else
127  print*,'mpi_comm_size unsuccessful'
128  end if
129 
130 
131  call mpi_comm_rank(tmp_mdls_comm,models_id, mpi_err)
132  if(mpi_err .eq. 0) then
133  print*,'mpi_comm_rank successful'
134  print*,'models_id = ',models_id
135  else
136  print*,'mpi_comm_rank unsuccessful'
137  end if
138 
139  mdlcolour = models_id/mdl_num_proc
140  call mpi_comm_split(tmp_mdls_comm,mdlcolour,models_id,mdl_mpi_comm,mpi_err)
141  if(mpi_err .eq. 0) then
142  print*,'mpi_comm_split successful'
143  else
144  print*,'mpi_comm_split unsuccessful'
145  end if
146 
147 
148  call mpi_comm_rank(mdl_mpi_comm,mdl_id,mpi_err)
149  if(mpi_err .eq. 0) then
150  print*,'mpi_comm_rank successful'
151  print*,'mdl_id = ',mdl_id
152  else
153  print*,'mpi_comm_rank unsuccessful'
154  end if
155 
156  if(mdl_id .eq. 0) then
157  cpl_colour = 9999
158  else
159  cpl_colour = mpi_undefined
160  end if
161  call mpi_comm_split(mpi_comm_world,cpl_colour,mdlcolour,cpl_mpi_comm,mpi_err)
162  if(mpi_err .eq. 0) then
163  print*,'mpi_comm_split successful'
164  else
165  print*,'mpi_comm_split unsuccessful'
166  end if
167 
168  if(mdl_id .eq. 0) then
169  call mpi_comm_size(cpl_mpi_comm,nens,mpi_err)
170  if(mpi_err .eq. 0) then
171  print*,'mpi_comm_size successful'
172  print*,'nens = ',nens
173  else
174  print*,'mpi_comm_size unsuccessful'
175  end if
176 
177  call mpi_comm_rank(cpl_mpi_comm,particle_id,mpi_err)
178  if(mpi_err .eq. 0) then
179  print*,'mpi_comm_rank successful'
180  print*,'particle_id = ',particle_id
181  else
182  print*,'mpi_comm_rank unsuccessful'
183  end if
184 
185  nda = world_size-models_size;nens = nens - nda
186  cpl_root = ((nda*particle_id)/nens)+nens
187  else
188  cpl_root = -1
189  end if
190  print*,'cpl_root = ',cpl_root
191  end subroutine initialise_mpi
192 end program minimal_model_comms
program minimal_model_comms
subroutine initialise_mpi(mdl_id, cpl_root, cpl_mpi_comm)