25 module ensemble_manager_mod
31 use mpp_mod,
only : input_nml_file
32 use fms2_io_mod,
only : fms2_io_set_filename_appendix=>set_filename_appendix
38 integer,
parameter :: MAX_ENSEMBLE_SIZE = 100
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
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
77 integer :: i, io_status, npes, ierr
79 namelist /ensemble_nml/ ensemble_size
81 read (input_nml_file, ensemble_nml, iostat=io_status)
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')
92 if (npes < ensemble_size)
then
93 call mpp_error(fatal,
'npes must be >= ensemble_size')
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')
127 integer,
intent(inout) :: pelist(:,:)
128 character(len=*),
intent(in),
optional :: name
130 if (
size(pelist,1) < ensemble_size) &
131 call mpp_error(fatal,
'get_ensemble_pelist: size of pelist 1st index < ensemble_size')
133 if(
present(name))
then
136 if (
size(pelist,2) < ocean_npes_pm)&
137 call mpp_error(fatal,
'get_ensemble_pelist: size of pelist 2nd index < ocean_npes_pm')
139 pelist(1:ensemble_size,1:ocean_npes_pm) = &
140 ensemble_pelist_ocean(1:ensemble_size,1:ocean_npes_pm)
143 if (
size(pelist,2) < atmos_npes_pm)&
144 call mpp_error(fatal,
'get_ensemble_pelist: size of pelist 2nd index < atmos_npes_pm')
146 pelist(1:ensemble_size,1:atmos_npes_pm) = &
147 ensemble_pelist_atmos(1:ensemble_size,1:atmos_npes_pm)
150 if (
size(pelist,2) < land_npes_pm)&
151 call mpp_error(fatal,
'get_ensemble_pelist: size of pelist 2nd index < land_npes_pm')
153 pelist(1:ensemble_size,1:land_npes_pm) = &
154 ensemble_pelist_land(1:ensemble_size,1:land_npes_pm)
157 if (
size(pelist,2) < ice_npes_pm)&
158 call mpp_error(fatal,
'get_ensemble_pelist: size of pelist 2nd index < ice_npes_pm')
160 pelist(1:ensemble_size,1:ice_npes_pm) = &
161 ensemble_pelist_ice(1:ensemble_size,1:ice_npes_pm)
164 call mpp_error(fatal,
'get_ensemble_pelist: unknown argument name='//name)
167 if (
size(pelist,2) < total_npes_pm)&
168 call mpp_error(fatal,
'get_ensemble_pelist: size of pelist 2nd index < total_npes_pm')
170 pelist(1:ensemble_size,1:total_npes_pm) = &
171 ensemble_pelist(1:ensemble_size,1:total_npes_pm)
186 integer,
intent(inout) :: pelist(:)
187 character(len=*),
intent(in) :: name
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')
194 pelist(1:ensemble_size*ocean_npes_pm) = &
195 ensemble_pelist_ocean_filter(1:ensemble_size*ocean_npes_pm)
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')
201 pelist(1:ensemble_size*atmos_npes_pm) = &
202 ensemble_pelist_atmos_filter(1:ensemble_size*atmos_npes_pm)
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')
208 pelist(1:ensemble_size*land_npes_pm) = &
209 ensemble_pelist_land_filter(1:ensemble_size*land_npes_pm)
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')
215 pelist(1:ensemble_size*ice_npes_pm) = &
216 ensemble_pelist_ice_filter(1:ensemble_size*ice_npes_pm)
219 call mpp_error(fatal,
'get_ensemble_filter_pelist: unknown argument name='//name)
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
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')
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) )
260 if( concurrent .OR. atmos_npes+ocean_npes == npes )
then
261 ocean_pe_start = ensemble_size*atmos_npes
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
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
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,:)
338 write(pelist_name,
'(a,i2.2)')
'atm_ens',n
340 write(pelist_name,
'(a,i2.2)')
'ocn_ens',n
342 write(pelist_name,
'(a,i2.2)')
'lnd_ens',n
344 write(pelist_name,
'(a,i2.2)')
'ice_ens',n
348 write(pelist_name,
'(a,i2.2)')
'atm_ens',ensemble_id
350 write(pelist_name,
'(a,i2.2)')
'ocn_ens',ensemble_id
352 write(pelist_name,
'(a,i2.2)')
'lnd_ens',ensemble_id
354 write(pelist_name,
'(a,i2.2)')
'ice_ens',ensemble_id
358 ocean_npes_pm = ocean_npes
359 atmos_npes_pm = atmos_npes
360 land_npes_pm = land_npes
361 ice_npes_pm = ice_npes
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))
370 i=(n-1)*ocean_npes_pm + m
371 ensemble_pelist_ocean_filter(i) = ensemble_pelist_ocean(n,m)
374 i=(n-1)*atmos_npes_pm + m
375 ensemble_pelist_atmos_filter(i) = ensemble_pelist_atmos(n,m)
378 i=(n-1)*land_npes_pm + m
379 ensemble_pelist_land_filter(i) = ensemble_pelist_land(n,m)
382 i=(n-1)*ice_npes_pm + m
383 ensemble_pelist_ice_filter(i) = ensemble_pelist_ice(n,m)
387 write(pelist_name,
'(a)')
'ocn_filter'
390 write(pelist_name,
'(a)')
'atm_filter'
393 write(pelist_name,
'(a)')
'lnd_filter'
396 write(pelist_name,
'(a)')
'ice_filter'
404 if (ensemble_size > 1)
then
405 write( text,
'(a,i2.2)' )
'ens_', ensemble_id
407 call fms2_io_set_filename_appendix(trim(text))
413 end module ensemble_manager_mod
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...
integer function stdout()
This function returns the current standard fortran unit numbers for output.
subroutine mpp_declare_pelist(pelist, name, commID)
Declare a pelist.
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
integer function mpp_npes()
Returns processor count for current pelist.
integer function mpp_pe()
Returns processor ID.