FMS  2025.03
Flexible Modeling System
ensemble_manager.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup ensemble_manager_mod ensemble_manager_mod
20 !> @ingroup coupler
21 !> @brief Routines for setting up and managing ensembles and ensemble pe lists.
22 
23 !> @addtogroup ensemble_manager_mod
24 !> @{
25 module ensemble_manager_mod
26 
27 
28  use fms_mod, only : check_nml_error
29  use mpp_mod, only : mpp_npes, stdout, stdlog, mpp_error, fatal
30  use mpp_mod, only : mpp_pe, mpp_declare_pelist
31  use mpp_mod, only : input_nml_file
32  use fms2_io_mod, only : fms2_io_set_filename_appendix=>set_filename_appendix
33 
34  IMPLICIT NONE
35 
36  private
37 
38  integer, parameter :: MAX_ENSEMBLE_SIZE = 100
39 
40 
41  integer, allocatable, dimension(:,:) :: ensemble_pelist
42  integer, allocatable, dimension(:,:) :: ensemble_pelist_ocean
43  integer, allocatable, dimension(:,:) :: ensemble_pelist_atmos
44  integer, allocatable, dimension(:,:) :: ensemble_pelist_land
45  integer, allocatable, dimension(:,:) :: ensemble_pelist_ice
46  integer, allocatable, dimension(:) :: ensemble_pelist_ocean_filter
47  integer, allocatable, dimension(:) :: ensemble_pelist_atmos_filter
48  integer, allocatable, dimension(:) :: ensemble_pelist_land_filter
49  integer, allocatable, dimension(:) :: ensemble_pelist_ice_filter
50 
51  integer :: ensemble_size = 1
52  integer :: ensemble_id = 1
53  integer :: pe, total_npes_pm=0,ocean_npes_pm=0,atmos_npes_pm=0
54  integer :: land_npes_pm=0,ice_npes_pm=0
55 
57  public :: ensemble_pelist_setup
59 contains
60 
61 !> @brief Initializes @ref ensemble_manager_mod
62 !!
63 !> @throw FATAL, "ensemble_manager_mod: ensemble_nml variable ensemble_size must be a positive integer"
64 !! @throw FATAL, "ensemble_manager_mod: ensemble_nml variable ensemble_size should be no larger
65 !! than MAX_ENSEMBLE_SIZE, change ensemble_size or increase MAX_ENSEMBLE_SIZE"
66 !! @throw FATAL, "ensemble_size must be divis by npes"
67 !! @throw FATAL, "get_ensemble_pelist: size of pelist 1st index < ensemble_size"
68 !! @throw FATAL, "get_ensemble_pelist: size of pelist 2nd index < ocean_npes_pm"
69 !! @throw FATAL, "get_ensemble_pelist: size of pelist 2nd index < atmos_npes_pm"
70 !! @throw FATAL, "get_ensemble_pelist: size of pelist 2nd index < land_npes_pm"
71 !! @throw FATAL, "get_ensemble_pelist: size of pelist 2nd index < ice_npes_pm"
72 !! @throw FATAL, "get_ensemble_pelist: unknown argument name=[name]"
73 !! @throw FATAL, "get_ensemble_pelist: size of pelist 2nd index < total_npes_pm"
74  subroutine ensemble_manager_init()
75 
76 
77  integer :: i, io_status, npes, ierr
78 
79  namelist /ensemble_nml/ ensemble_size
80 
81  read (input_nml_file, ensemble_nml, iostat=io_status)
82  ierr = check_nml_error(io_status, 'ensemble_nml')
83 
84  if(ensemble_size < 1) call mpp_error(fatal, &
85  'ensemble_manager_mod: ensemble_nml variable ensemble_size must be a positive integer')
86  if(ensemble_size > max_ensemble_size) call mpp_error(fatal, &
87  'ensemble_manager_mod: ensemble_nml variable ensemble_size should be no larger than MAX_ENSEMBLE_SIZE, '// &
88  'change ensemble_size or increase MAX_ENSEMBLE_SIZE')
89 
90  pe = mpp_pe()
91  npes = mpp_npes()
92  if (npes < ensemble_size) then
93  call mpp_error(fatal,'npes must be >= ensemble_size')
94  endif
95  total_npes_pm = npes/ensemble_size
96  if (mod(npes, total_npes_pm) /= 0) call mpp_error(fatal,'ensemble_size must be divis by npes')
97 
98  call mpp_declare_pelist((/(i,i=0,npes-1)/),'_ens0') ! for ensemble driver
99 
100  end subroutine ensemble_manager_init
101 
102  !> @brief Getter function for ensemble_id
103  !! @return integer of ensemble id
104  function get_ensemble_id()
105  integer :: get_ensemble_id
106  get_ensemble_id = ensemble_id
107  end function get_ensemble_id
108 
109  !> @brief Returns ensemble size integer array
110  !! @return integer array of sizes
111  function get_ensemble_size()
112 
113  integer, dimension(6) :: get_ensemble_size
114 
115  get_ensemble_size(1) = ensemble_size
116  get_ensemble_size(2) = total_npes_pm
117  get_ensemble_size(3) = ocean_npes_pm
118  get_ensemble_size(4) = atmos_npes_pm
119  get_ensemble_size(5) = land_npes_pm
120  get_ensemble_size(6) = ice_npes_pm
121 
122  end function get_ensemble_size
123 
124  !> @brief Gets pe list for current ensemble or a given ensemble component.
125  subroutine get_ensemble_pelist(pelist, name)
126 
127  integer, intent(inout) :: pelist(:,:) !< Ensemble pelist
128  character(len=*), intent(in), optional :: name !< Component name.
129 
130  if (size(pelist,1) < ensemble_size) &
131  call mpp_error(fatal,'get_ensemble_pelist: size of pelist 1st index < ensemble_size')
132 
133  if(present(name)) then
134  select case(name)
135  case('ocean')
136  if (size(pelist,2) < ocean_npes_pm)&
137  call mpp_error(fatal,'get_ensemble_pelist: size of pelist 2nd index < ocean_npes_pm')
138  pelist = 0
139  pelist(1:ensemble_size,1:ocean_npes_pm) = &
140  ensemble_pelist_ocean(1:ensemble_size,1:ocean_npes_pm)
141 
142  case('atmos')
143  if (size(pelist,2) < atmos_npes_pm)&
144  call mpp_error(fatal,'get_ensemble_pelist: size of pelist 2nd index < atmos_npes_pm')
145  pelist = 0
146  pelist(1:ensemble_size,1:atmos_npes_pm) = &
147  ensemble_pelist_atmos(1:ensemble_size,1:atmos_npes_pm)
148 
149  case('land')
150  if (size(pelist,2) < land_npes_pm)&
151  call mpp_error(fatal,'get_ensemble_pelist: size of pelist 2nd index < land_npes_pm')
152  pelist = 0
153  pelist(1:ensemble_size,1:land_npes_pm) = &
154  ensemble_pelist_land(1:ensemble_size,1:land_npes_pm)
155 
156  case('ice')
157  if (size(pelist,2) < ice_npes_pm)&
158  call mpp_error(fatal,'get_ensemble_pelist: size of pelist 2nd index < ice_npes_pm')
159  pelist = 0
160  pelist(1:ensemble_size,1:ice_npes_pm) = &
161  ensemble_pelist_ice(1:ensemble_size,1:ice_npes_pm)
162 
163  case default
164  call mpp_error(fatal,'get_ensemble_pelist: unknown argument name='//name)
165  end select
166  else
167  if (size(pelist,2) < total_npes_pm)&
168  call mpp_error(fatal,'get_ensemble_pelist: size of pelist 2nd index < total_npes_pm')
169  pelist = 0
170  pelist(1:ensemble_size,1:total_npes_pm) = &
171  ensemble_pelist(1:ensemble_size,1:total_npes_pm)
172  endif
173 
174  return
175  end subroutine get_ensemble_pelist
176 
177 !> @brief Gets filter pelist for a given ensemble component.
178 !!
179 !! @throw FATAL, "get_ensemble_filter_pelist: size of pelist argument < ensemble_size * ocean_npes_pm"
180 !! @throw FATAL, "get_ensemble_filter_pelist: size of pelist argument < ensemble_size * atmos_npes_pm"
181 !! @throw FATAL, "get_ensemble_filter_pelist: size of pelist argument < ensemble_size * land_npes_pm"
182 !! @throw FATAL, "get_ensemble_filter_pelist: size of pelist argument < ensemble_size * ice_npes_pm"
183 !! @throw FATAL, "get_ensemble_filter_pelist: unknown argument name=[name]"
184  subroutine get_ensemble_filter_pelist(pelist, name)
185 
186  integer, intent(inout) :: pelist(:) !< Returned filter pe list
187  character(len=*), intent(in) :: name !< Ensemble component name
188 
189  select case(name)
190  case('ocean')
191  if (size(pelist) < ensemble_size * ocean_npes_pm)&
192  call mpp_error(fatal,'get_ensemble_filter_pelist: size of pelist argument < ensemble_size * ocean_npes_pm')
193  pelist = 0
194  pelist(1:ensemble_size*ocean_npes_pm) = &
195  ensemble_pelist_ocean_filter(1:ensemble_size*ocean_npes_pm)
196 
197  case('atmos')
198  if (size(pelist) < ensemble_size * atmos_npes_pm)&
199  call mpp_error(fatal,'get_ensemble_filter_pelist: size of pelist argument < ensemble_size * atmos_npes_pm')
200  pelist = 0
201  pelist(1:ensemble_size*atmos_npes_pm) = &
202  ensemble_pelist_atmos_filter(1:ensemble_size*atmos_npes_pm)
203 
204  case('land')
205  if (size(pelist) < ensemble_size * land_npes_pm)&
206  call mpp_error(fatal,'get_ensemble_filter_pelist: size of pelist argument < ensemble_size * land_npes_pm')
207  pelist = 0
208  pelist(1:ensemble_size*land_npes_pm) = &
209  ensemble_pelist_land_filter(1:ensemble_size*land_npes_pm)
210 
211  case('ice')
212  if (size(pelist) < ensemble_size * ice_npes_pm)&
213  call mpp_error(fatal,'get_ensemble_filter_pelist: size of pelist argument < ensemble_size * ice_npes_pm')
214  pelist = 0
215  pelist(1:ensemble_size*ice_npes_pm) = &
216  ensemble_pelist_ice_filter(1:ensemble_size*ice_npes_pm)
217 
218  case default
219  call mpp_error(fatal,'get_ensemble_filter_pelist: unknown argument name='//name)
220  end select
221 
222 
223  return
224  end subroutine get_ensemble_filter_pelist
225 
226 !nnz: I think the following block of code should be contained in a subroutine
227 ! to consolidate and ensure the consistency of declaring the various pelists.
228 !> @brief Sets up pe list for an ensemble.
229 !!
230 !! @throw FATAL, "ensemble_manager_mod: land_npes > atmos_npes"
231 !! @throw FATAL, "ensemble_manager_mod: ice_npes > atmos_npes"
232  subroutine ensemble_pelist_setup(concurrent, atmos_npes, ocean_npes, land_npes, ice_npes, &
233  Atm_pelist, Ocean_pelist, Land_pelist, Ice_pelist)
234  logical, intent(in) :: concurrent
235  integer, intent(in) :: atmos_npes, ocean_npes
236  integer, intent(in) :: land_npes, ice_npes
237  integer, dimension(:), intent(inout) :: atm_pelist, ocean_pelist
238  integer, dimension(:), intent(inout) :: land_pelist, ice_pelist
239  integer :: atmos_pe_start, atmos_pe_end, ocean_pe_start, ocean_pe_end
240  integer :: land_pe_start, land_pe_end, ice_pe_start, ice_pe_end
241  character(len=10) :: pelist_name, text
242  integer :: npes, n, m ,i
243 
244  npes = total_npes_pm
245 
246  ! make sure land_npes and ice_npes are not greater than atmos_npes
247  if(land_npes > atmos_npes) call mpp_error(fatal, 'ensemble_manager_mod: land_npes > atmos_npes')
248  if(ice_npes > atmos_npes) call mpp_error(fatal, 'ensemble_manager_mod: ice_npes > atmos_npes')
249 
250  allocate(ensemble_pelist(ensemble_size,total_npes_pm))
251  allocate(ensemble_pelist_ocean(1:ensemble_size, 1:ocean_npes) )
252  allocate(ensemble_pelist_atmos(1:ensemble_size, 1:atmos_npes) )
253  allocate(ensemble_pelist_land(1:ensemble_size, 1:land_npes) )
254  allocate(ensemble_pelist_ice(1:ensemble_size, 1:ice_npes) )
255 
256  atmos_pe_start = 0
257  ocean_pe_start = 0
258  land_pe_start = 0
259  ice_pe_start = 0
260  if( concurrent .OR. atmos_npes+ocean_npes == npes )then
261  ocean_pe_start = ensemble_size*atmos_npes
262  endif
263  do n=1,ensemble_size
264  atmos_pe_end = atmos_pe_start + atmos_npes - 1
265  ocean_pe_end = ocean_pe_start + ocean_npes - 1
266  land_pe_end = land_pe_start + land_npes - 1
267  ice_pe_end = ice_pe_start + ice_npes - 1
268  ensemble_pelist_atmos(n, 1:atmos_npes) = (/(i,i=atmos_pe_start,atmos_pe_end)/)
269  ensemble_pelist_ocean(n, 1:ocean_npes) = (/(i,i=ocean_pe_start,ocean_pe_end)/)
270  ensemble_pelist_land(n, 1:land_npes) = (/(i,i=land_pe_start, land_pe_end)/)
271  ensemble_pelist_ice(n, 1:ice_npes) = (/(i,i=ice_pe_start, ice_pe_end)/)
272  ensemble_pelist(n, 1:atmos_npes) = ensemble_pelist_atmos(n, 1:atmos_npes)
273  if( concurrent .OR. atmos_npes+ocean_npes == npes ) &
274  ensemble_pelist(n, atmos_npes+1:npes) = ensemble_pelist_ocean(n, 1:ocean_npes)
275  if(any(ensemble_pelist(n,:) == pe)) ensemble_id = n
276  write(pelist_name,'(a,i2.2)') '_ens',n
277  call mpp_declare_pelist(ensemble_pelist(n,:), trim(pelist_name))
278  atmos_pe_start = atmos_pe_end + 1
279  ocean_pe_start = ocean_pe_end + 1
280  land_pe_start = atmos_pe_start
281  ice_pe_start = atmos_pe_start
282  enddo
283 
284  atm_pelist(:) = ensemble_pelist_atmos(ensemble_id,:)
285  ocean_pelist(:) = ensemble_pelist_ocean(ensemble_id,:)
286  land_pelist(:) = ensemble_pelist_land(ensemble_id,:)
287  ice_pelist(:) = ensemble_pelist_ice(ensemble_id,:)
288 
289  ! write(pelist_name,'(a,i2.2)') '_ocn_ens',ensemble_id
290  ! call mpp_declare_pelist(Ocean%pelist , trim(pelist_name) )
291 
292  ! write(pelist_name,'(a,i2.2)') '_atm_ens',ensemble_id
293  ! call mpp_declare_pelist(Atm%pelist , trim(pelist_name) )
294  !
295  !nnz: The above is sufficient for non-concurrent mode.
296  ! BUT
297  ! For atmosphere_init to work in ensemble, concurrent mode
298  ! the following Atm_pelist should be declared (per ensemble member)
299  ! instead of the above Atm%pelist!
300  !
301  ! allocate( Atm_pelist(1:ensemble_size, 1:atmos_npes) )
302  ! do n=1,ensemble_size
303  ! do i=1, atmos_npes
304  ! Atm_pelist(n, i) = ensemble_pelist(n, i)
305  ! enddo
306  ! write(pelist_name,'(a,i2.2)') '_atm_ens',n
307  ! call mpp_declare_pelist(Atm_pelist(n,:) , trim(pelist_name) )
308  ! enddo
309  !
310  ! The way I understand this with the help of Totalview is:
311  ! With mpp_declare_pelist(Atm%pelist)
312  ! When we are in fv_arrays_init when mp_init(comID) is called
313  ! comID is the same for the atmospheric PE's for both ensemble members
314  ! since peset(5)%id is the same (7) for those PE's, so the PE count is double what it should be inside
315  ! mp_init().
316  ! It is also true that for Ocean PE's, peset(4)%id is the same (6) for Ocean PE's in both ensemble members
317  ! but for Ocean it is not a problem because Ocean is not trying to create new communicators
318  ! from this peset whereas ATM does (vis mp_init).
319  !
320  ! Who sets peset(i)%id ? Can it be modified to assign different %id for the two subsets.
321  ! peset(i)%id = 0 for Ocean PE's on ATM pesets and for ATM PE's on Ocean pesets.
322  !
323  ! With mpp_declare_pelist(Atm_pelist(n,:)) n=1,...,ensemble_size
324  ! we get separate pesets for each ATM ensemble member and each with a different %id and mp_init is cured.
325  !
326  ! There is also a matter of precedence. If we have both calls
327  ! call mpp_declare_pelist(Atm%pelist , trim(pelist_name) )
328  ! and
329  ! call mpp_declare_pelist(Atm_pelist(n,:) , trim(pelist_name) )
330  ! then concurrent run fails because with call mpp_set_current_pelist( Atm%pelist )
331  ! peset(i) is searched for i=1,2,... and the first pelist that matches argument, its peset is set as current.
332  !
333  ! To be consistent with ATM and OCEAN we can do the following
334  ! (eventhough mpp_declare_pelist(Ocean%pelist) is adequate right now.)
335 
336  if( concurrent )then
337  do n=1,ensemble_size
338  write(pelist_name,'(a,i2.2)') 'atm_ens',n
339  call mpp_declare_pelist(ensemble_pelist_atmos(n,:) , trim(pelist_name) )
340  write(pelist_name,'(a,i2.2)') 'ocn_ens',n
341  call mpp_declare_pelist(ensemble_pelist_ocean(n,:) , trim(pelist_name) )
342  write(pelist_name,'(a,i2.2)') 'lnd_ens',n
343  call mpp_declare_pelist(ensemble_pelist_land(n,:) , trim(pelist_name) )
344  write(pelist_name,'(a,i2.2)') 'ice_ens',n
345  call mpp_declare_pelist(ensemble_pelist_ice(n,:) , trim(pelist_name) )
346  enddo
347  else
348  write(pelist_name,'(a,i2.2)') 'atm_ens',ensemble_id
349  call mpp_declare_pelist(atm_pelist , trim(pelist_name) )
350  write(pelist_name,'(a,i2.2)') 'ocn_ens',ensemble_id
351  call mpp_declare_pelist(ocean_pelist , trim(pelist_name) )
352  write(pelist_name,'(a,i2.2)') 'lnd_ens',ensemble_id
353  call mpp_declare_pelist(land_pelist , trim(pelist_name) )
354  write(pelist_name,'(a,i2.2)') 'ice_ens',ensemble_id
355  call mpp_declare_pelist(ice_pelist , trim(pelist_name) )
356  endif
357 
358  ocean_npes_pm = ocean_npes
359  atmos_npes_pm = atmos_npes
360  land_npes_pm = land_npes
361  ice_npes_pm = ice_npes
362 
363  !Declare pelist of all Ocean, Atmos, Land and Ice pes across all ensembles ( filters )
364  allocate(ensemble_pelist_ocean_filter(ensemble_size * ocean_npes_pm))
365  allocate(ensemble_pelist_atmos_filter(ensemble_size * atmos_npes_pm))
366  allocate(ensemble_pelist_land_filter(ensemble_size * land_npes_pm))
367  allocate(ensemble_pelist_ice_filter(ensemble_size * ice_npes_pm))
368  do n=1,ensemble_size
369  do m=1,ocean_npes_pm
370  i=(n-1)*ocean_npes_pm + m
371  ensemble_pelist_ocean_filter(i) = ensemble_pelist_ocean(n,m)
372  enddo
373  do m=1,atmos_npes_pm
374  i=(n-1)*atmos_npes_pm + m
375  ensemble_pelist_atmos_filter(i) = ensemble_pelist_atmos(n,m)
376  enddo
377  do m=1,land_npes_pm
378  i=(n-1)*land_npes_pm + m
379  ensemble_pelist_land_filter(i) = ensemble_pelist_land(n,m)
380  enddo
381  do m=1,ice_npes_pm
382  i=(n-1)*ice_npes_pm + m
383  ensemble_pelist_ice_filter(i) = ensemble_pelist_ice(n,m)
384  enddo
385  enddo
386 
387  write(pelist_name,'(a)') 'ocn_filter'
388  call mpp_declare_pelist(ensemble_pelist_ocean_filter, trim(pelist_name) )
389 
390  write(pelist_name,'(a)') 'atm_filter'
391  call mpp_declare_pelist(ensemble_pelist_atmos_filter, trim(pelist_name) )
392 
393  write(pelist_name,'(a)') 'lnd_filter'
394  call mpp_declare_pelist(ensemble_pelist_land_filter, trim(pelist_name) )
395 
396  write(pelist_name,'(a)') 'ice_filter'
397  call mpp_declare_pelist(ensemble_pelist_ice_filter, trim(pelist_name) )
398 
399  !
400  !Rename output files to identify the ensemble
401  !If ensemble_size=1 do not rename files so that the same coupler
402  !can be used for non-ensemble experiments
403  !
404  if (ensemble_size > 1) then
405  write( text,'(a,i2.2)' ) 'ens_', ensemble_id
406  !Append ensemble_id to the restart filenames
407  call fms2_io_set_filename_appendix(trim(text))
408  endif
409 
410  end subroutine ensemble_pelist_setup
411 
412 
413 end module ensemble_manager_mod
414 !> @}
415 ! close documentation grouping
integer function, public get_ensemble_id()
Getter function for ensemble_id.
integer function, dimension(6), public get_ensemble_size()
Returns ensemble size integer array.
subroutine, public ensemble_pelist_setup(concurrent, atmos_npes, ocean_npes, land_npes, ice_npes, Atm_pelist, Ocean_pelist, Land_pelist, Ice_pelist)
Sets up pe list for an ensemble.
subroutine, public ensemble_manager_init()
Initializes ensemble_manager_mod.
subroutine, public get_ensemble_pelist(pelist, name)
Gets pe list for current ensemble or a given ensemble component.
subroutine, public get_ensemble_filter_pelist(pelist, name)
Gets filter pelist for a given ensemble component.
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
Definition: fms.F90:524
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:43
subroutine mpp_declare_pelist(pelist, name, commID)
Declare a pelist.
Definition: mpp_util.inc:470
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
Definition: mpp_util.inc:59
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:421
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
Error handler.
Definition: mpp.F90:382