Line data Source code
1 : module subcol_utils
2 : !---------------------------------------------------------------------------
3 : ! Purpose:
4 : !
5 : ! Provides utilities to support subcolumns
6 : !
7 : !---------------------------------------------------------------------------
8 :
9 : use shr_kind_mod, only: r8=>shr_kind_r8, r4=>shr_kind_r4, i4=>shr_kind_i4
10 : use infnan, only: nan, assignment(=)
11 : use physics_types, only: physics_state, physics_ptend, physics_tend, physics_tend_alloc, physics_state_alloc
12 : use ppgrid, only: pcols, psubcols, pver
13 : use constituents, only: pcnst
14 : use cam_abortutils, only: endrun
15 : use pio, only: var_desc_t
16 : use subcol_pack_mod, only: subcol_get_nsubcol, subcol_set_nsubcol, subcol_get_indcol
17 :
18 : implicit none
19 :
20 : private
21 : save
22 :
23 : integer, target, allocatable :: filter2d(:,:)
24 : real(r8),target, allocatable :: weight2d(:,:)
25 : logical :: weight_set, filter_set
26 :
27 : !! The active subcolumn scheme
28 : character(len=16) :: subcol_scheme = 'off'
29 :
30 : !! Public interface functions which do not depend on the subcolumn scheme
31 :
32 : public :: subcol_utils_init ! Initialize module data (e.g., nsubcol2d)
33 : public :: subcol_get_filter ! return the filter values
34 : public :: subcol_set_filter ! set the filter values
35 : public :: subcol_get_weight ! return the weight values
36 : public :: subcol_set_weight ! set the weight values
37 :
38 : public :: subcol_field_copy ! copy a physics buffer field into one with subcolumn dimensions
39 : public :: subcol_ptend_copy ! copy a physics_ptend object into one with subcolumn dimensions
40 : public :: subcol_set_subcols ! set nsubcols and copy state & tend objects into one with subcolumn dimensions
41 :
42 : public :: subcol_field_avg_shr ! Average subcol fields back into GBA fields
43 : public :: subcol_ptend_avg_shr ! average subcolumn ptend to grid ptend
44 :
45 : public :: subcol_field_get_firstsubcol ! Retrieve the first subcolumn and assign to grid
46 : public :: subcol_ptend_get_firstsubcol ! retrieve the first subcolumn from the ptend fields and assign to grid ptend
47 :
48 : public :: subcol_utils_init_restart ! Initialize restart with subcolumn specific fields
49 : public :: subcol_utils_read_restart ! Read subcolumn specific fields from restart
50 : public :: subcol_utils_write_restart ! Write subcolumn specific fields for restart
51 : public :: is_filter_set ! True if filters for averaging have been set
52 : public :: is_weight_set ! True if weights for averaging have been set
53 : public :: is_subcol_on ! true is any subcol_scheme other than "off" is set
54 : public :: subcol_get_scheme ! Return the active subcolumn scheme name
55 : public :: subcol_utils_readnl ! Set the active scheme based on namelist
56 :
57 : interface subcol_field_avg_shr
58 : ! TYPE int,double,real
59 : ! DIMS 1,2
60 : module procedure subcol_field_avg_shr_{DIMS}d{TYPE}
61 : end interface
62 :
63 : interface subcol_avg_inter
64 : ! TYPE int,double,real
65 : module procedure subcol_avg_inter_{TYPE}
66 : end interface
67 :
68 : interface subcol_avg
69 : ! TYPE int,double,real
70 : module procedure subcol_avg_{TYPE}
71 : end interface
72 :
73 : interface subcol_field_get_firstsubcol
74 : ! TYPE int,double,real
75 : ! DIMS 1,2
76 : module procedure subcol_field_get_firstsubcol_{DIMS}d{TYPE}
77 : end interface
78 :
79 : interface subcol_state_field_copy
80 : ! TYPE double
81 : ! DIMS 1,2,3
82 : module procedure subcol_state_field_copy_{DIMS}d{TYPE}
83 : end interface
84 :
85 : interface subcol_field_copy
86 : ! TYPE int,double,real
87 : ! DIMS 1,2,3,4,5
88 : module procedure subcol_field_copy_{DIMS}d{TYPE}
89 : end interface
90 :
91 : type(var_desc_t) :: weight2d_desc, filter2d_desc
92 : integer :: subcol_dimid = -1 ! subcol dimension for restart
93 :
94 : integer :: ret_nan_int
95 : real(r8):: ret_nan_double
96 : real(r4):: ret_nan_real
97 :
98 : integer :: fillval_int
99 : real(r8):: fillval_double
100 : real(r4):: fillval_real
101 :
102 : contains
103 :
104 768 : subroutine subcol_allocate_internal()
105 : use ppgrid, only: begchunk, endchunk
106 : use subcol_pack_mod, only: subcol_pack_allocate
107 :
108 768 : call subcol_pack_allocate()
109 :
110 768 : if (allocated(filter2d)) then
111 0 : deallocate(filter2d)
112 : end if
113 2304 : allocate(filter2d(pcols*psubcols, begchunk:endchunk))
114 53400 : filter2d = 0
115 :
116 768 : if (allocated(weight2d)) then
117 0 : deallocate(weight2d)
118 : end if
119 2304 : allocate(weight2d(pcols*psubcols, begchunk:endchunk))
120 53400 : weight2d = 0._r8
121 :
122 768 : end subroutine subcol_allocate_internal
123 :
124 768 : subroutine subcol_utils_init(subcol_scheme_init)
125 :
126 : character(len=*), optional, intent(in) :: subcol_scheme_init ! Name of subcolumn generator
127 :
128 768 : call subcol_allocate_internal()
129 :
130 768 : ret_nan_int = 0
131 768 : ret_nan_double = nan
132 768 : ret_nan_real = nan
133 :
134 768 : fillval_int = 0
135 768 : fillval_double = 0._r8
136 768 : fillval_real = 0._r4
137 :
138 768 : weight_set = .false.
139 768 : filter_set = .false.
140 :
141 768 : end subroutine subcol_utils_init
142 :
143 0 : subroutine subcol_get_filter(lchnk, filter)
144 : !-----------------------------------------------------------------------
145 : ! Retrieve the filter module variable
146 : !-----------------------------------------------------------------------
147 :
148 : integer, intent(in) :: lchnk
149 : integer, intent(out) :: filter(:)
150 :
151 0 : filter(:) = filter2d(:,lchnk)
152 0 : end subroutine subcol_get_filter
153 :
154 0 : subroutine subcol_get_weight(lchnk, weight)
155 : !-----------------------------------------------------------------------
156 : ! Retrieve the weight module variable
157 : !-----------------------------------------------------------------------
158 :
159 : integer, intent(in) :: lchnk
160 : real(r8), intent(out) :: weight(:)
161 :
162 0 : weight(:) = weight2d(:,lchnk)
163 0 : end subroutine subcol_get_weight
164 :
165 0 : integer function subcol_get_ncol(lchnk) result(ncol)
166 : !-----------------------------------------------------------------------
167 : ! Compute the number of (sub)columns for a chunk
168 : ! NB: This is considered an internal function so it can use nsubcol2d
169 : !-----------------------------------------------------------------------
170 : integer, intent(in) :: lchnk
171 :
172 : integer :: nsubcol(pcols)
173 :
174 0 : call subcol_get_nsubcol(lchnk, nsubcol)
175 0 : ncol = SUM(nsubcol)
176 :
177 0 : end function subcol_get_ncol
178 :
179 0 : subroutine subcol_set_filter(lchnk, filter)
180 : !-----------------------------------------------------------------------
181 : ! Set the filter module variable
182 : !-----------------------------------------------------------------------
183 :
184 : integer, intent(in) :: lchnk
185 : integer, intent(in) :: filter(:)
186 :
187 0 : filter2d(:,lchnk) = filter(:)
188 0 : filter_set = .true.
189 0 : end subroutine subcol_set_filter
190 :
191 0 : subroutine subcol_set_weight(lchnk, weight)
192 : !-----------------------------------------------------------------------
193 : ! Set the weight module variable
194 : !-----------------------------------------------------------------------
195 :
196 : integer, intent(in) :: lchnk
197 : real(r8), intent(in) :: weight(:)
198 :
199 0 : weight2d(:,lchnk) = weight(:)
200 0 : weight_set = .true.
201 0 : end subroutine subcol_set_weight
202 :
203 0 : logical function is_weight_set()
204 0 : is_weight_set=weight_set
205 0 : end function is_weight_set
206 :
207 0 : logical function is_filter_set()
208 0 : is_filter_set=filter_set
209 0 : end function is_filter_set
210 :
211 10434240 : logical function is_subcol_on()
212 10434240 : is_subcol_on = (trim(subcol_scheme) /= 'off')
213 10434240 : end function is_subcol_on
214 :
215 9216 : character(len=16) function subcol_get_scheme()
216 9216 : subcol_get_scheme = trim(subcol_scheme)
217 9216 : end function subcol_get_scheme
218 :
219 1536 : subroutine subcol_utils_readnl(nlfile)
220 : use namelist_utils, only: find_group_name
221 : use units, only: getunit, freeunit
222 : use spmd_utils, only: masterproc, mpi_character, masterprocid, mpicom
223 :
224 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
225 :
226 : ! Local variables
227 : integer :: unitn, ierr
228 :
229 : namelist /subcol_nl/ subcol_scheme
230 :
231 1536 : if (masterproc) then
232 2 : unitn = getunit()
233 2 : open( unitn, file=trim(nlfile), status='old' )
234 2 : call find_group_name(unitn, 'subcol_nl', status=ierr)
235 2 : if (ierr == 0) then
236 0 : read(unitn, subcol_nl, iostat=ierr)
237 0 : if (ierr /= 0) then
238 0 : call endrun('subcol_readnl: ERROR reading namelist')
239 : end if
240 : end if
241 2 : close(unitn)
242 2 : call freeunit(unitn)
243 : end if
244 :
245 : #ifdef SPMD
246 : ! Broadcast namelist variables
247 1536 : call mpi_bcast(subcol_scheme, len(subcol_scheme), mpi_character, masterprocid, mpicom, ierr)
248 : #endif
249 :
250 1536 : end subroutine subcol_utils_readnl
251 :
252 : ! TYPE int,double,real
253 : ! DIMS 1,2,3,4,5
254 0 : subroutine subcol_field_copy_{DIMS}d{TYPE} (field, lchnk, field_sc)
255 : !-----------------------------------------------------------------------
256 : ! Copy a pbuf field dimensioned pcols to one with pcols*psubcols and fill appropriately
257 : !-----------------------------------------------------------------------
258 :
259 : {VTYPE}, intent(in) :: field{DIMSTR}
260 : integer, intent(in) :: lchnk
261 : {VTYPE}, intent(out) :: field_sc{DIMSTR} ! intent out
262 :
263 : ! Local Variables
264 : integer :: ncol
265 : integer :: indcol(pcols*psubcols)
266 :
267 0 : if (size(field,dim=1) .ne. pcols) then
268 0 : call endrun('subcol_field_copy error: only fields with first dimension pcols may use this routine')
269 : end if
270 :
271 0 : field_sc = fillval_{TYPE}
272 0 : ncol = subcol_get_ncol(lchnk)
273 0 : call subcol_get_indcol(lchnk, indcol)
274 :
275 : #if ({DIMS} == 1)
276 0 : field_sc(:ncol) = field(indcol(:ncol))
277 : #endif
278 :
279 : #if ({DIMS} == 2)
280 0 : field_sc(:ncol,:) = field(indcol(:ncol),:)
281 : #endif
282 :
283 : #if ({DIMS} == 3)
284 0 : field_sc(:ncol,:,:) = field(indcol(:ncol),:,:)
285 : #endif
286 :
287 : #if ({DIMS} == 4)
288 0 : field_sc(:ncol,:,:,:) = field(indcol(:ncol),:,:,:)
289 : #endif
290 :
291 : #if ({DIMS} == 5)
292 0 : field_sc(:ncol,:,:,:,:) = field(indcol(:ncol),:,:,:,:)
293 : #endif
294 :
295 0 : end subroutine subcol_field_copy_{DIMS}d{TYPE}
296 :
297 : ! TYPE double
298 : ! DIMS 1,2,3
299 0 : subroutine subcol_state_field_copy_{DIMS}d{TYPE} (field, state, field_sc)
300 : !-----------------------------------------------------------------------
301 : ! Copy a state field dimensioned with pcols to one with pcols*psubcols and fill appropriately
302 : !-----------------------------------------------------------------------
303 :
304 : {VTYPE}, intent(in) :: field{DIMSTR}
305 : type(physics_state), intent(in) :: state
306 : {VTYPE}, allocatable :: field_sc{DIMSTR} ! intent out
307 :
308 : integer :: dim2, dim3
309 : integer :: indcol(pcols*psubcols)
310 :
311 0 : if (size(field,dim=1) .ne. pcols) then
312 0 : call endrun('subcol_state_field_copy error: only fields with first dimension pcols may use this routine')
313 : end if
314 :
315 0 : call subcol_get_indcol(state%lchnk, indcol)
316 :
317 : #if ({DIMS} == 1)
318 0 : if (.not. allocated(field_sc)) then
319 0 : allocate(field_sc(pcols*psubcols))
320 : end if
321 :
322 0 : field_sc = 0._r8
323 0 : field_sc(:state%ncol) = field(indcol(:state%ncol))
324 : #endif
325 :
326 : #if ({DIMS} == 2)
327 0 : if (.not. allocated(field_sc)) then
328 0 : dim2 = size(field,dim=2)
329 0 : allocate(field_sc(pcols*psubcols,dim2))
330 : end if
331 :
332 0 : field_sc = 0._r8
333 0 : field_sc(:state%ncol,:) = field(indcol(:state%ncol),:)
334 : #endif
335 :
336 : #if ({DIMS} == 3)
337 0 : if (.not. allocated(field_sc)) then
338 0 : dim2 = size(field,dim=2)
339 0 : dim3 = size(field,dim=3)
340 0 : allocate(field_sc(pcols*psubcols,dim2,dim3))
341 : end if
342 :
343 0 : field_sc = 0._r8
344 0 : field_sc(:state%ncol,:,:) = field(indcol(:state%ncol),:,:)
345 : #endif
346 :
347 0 : end subroutine subcol_state_field_copy_{DIMS}d{TYPE}
348 :
349 0 : subroutine subcol_tend_copy(tend, state_sc, tend_sc)
350 : !-----------------------------------------------------------------------
351 : ! Copy all of tend to subcolumns in tend_sc, allocating tend_sc if necessary
352 : !-----------------------------------------------------------------------
353 : type(physics_tend), intent(inout) :: tend
354 : type(physics_state), intent(in) :: state_sc ! subcolumn state
355 : type(physics_tend), intent(inout) :: tend_sc ! subcolumn tend
356 :
357 :
358 0 : if (.not. allocated(tend%dtdt)) then
359 0 : call physics_tend_alloc(tend_sc, state_sc%psetcols)
360 : end if
361 :
362 0 : tend_sc%psetcols = pcols*psubcols
363 0 : call subcol_state_field_copy(tend%dtdt, state_sc, tend_sc%dtdt)
364 0 : call subcol_state_field_copy(tend%dudt, state_sc, tend_sc%dudt)
365 0 : call subcol_state_field_copy(tend%dvdt, state_sc, tend_sc%dvdt)
366 0 : call subcol_state_field_copy(tend%flx_net, state_sc, tend_sc%flx_net)
367 0 : call subcol_state_field_copy(tend%te_tnd, state_sc, tend_sc%te_tnd)
368 0 : call subcol_state_field_copy(tend%tw_tnd, state_sc, tend_sc%tw_tnd)
369 :
370 0 : end subroutine subcol_tend_copy
371 :
372 0 : subroutine subcol_state_copy(state, state_sc)
373 : !-----------------------------------------------------------------------
374 : ! Copy all of state to subcolumns in state_sc, allocating state_sc if necessary
375 : !-----------------------------------------------------------------------
376 : type(physics_state), intent(in) :: state
377 : type(physics_state), intent(inout) :: state_sc ! subcolumn state
378 :
379 : !
380 : ! Local variables
381 : !
382 : integer :: ngrdcol
383 :
384 0 : if (.not. allocated(state_sc%lat)) then
385 0 : call endrun('subcol_state_copy: user must allocate state_sc prior to calling this routine')
386 : end if
387 :
388 0 : ngrdcol = state%ngrdcol
389 :
390 0 : call subcol_state_hdrinit(state, state_sc)
391 :
392 0 : call subcol_state_field_copy(state%lat, state_sc, state_sc%lat)
393 0 : call subcol_state_field_copy(state%lon, state_sc, state_sc%lon)
394 0 : call subcol_state_field_copy(state%ps, state_sc, state_sc%ps)
395 0 : call subcol_state_field_copy(state%psdry, state_sc, state_sc%psdry)
396 0 : call subcol_state_field_copy(state%phis, state_sc, state_sc%phis)
397 0 : call subcol_state_field_copy(state%te_ini, state_sc, state_sc%te_ini)
398 0 : call subcol_state_field_copy(state%te_cur, state_sc, state_sc%te_cur)
399 0 : call subcol_state_field_copy(state%tw_ini, state_sc, state_sc%tw_ini)
400 0 : call subcol_state_field_copy(state%tw_cur, state_sc, state_sc%tw_cur)
401 0 : call subcol_state_field_copy(state%t, state_sc, state_sc%t)
402 0 : call subcol_state_field_copy(state%u, state_sc, state_sc%u)
403 0 : call subcol_state_field_copy(state%v, state_sc, state_sc%v)
404 0 : call subcol_state_field_copy(state%s, state_sc, state_sc%s)
405 0 : call subcol_state_field_copy(state%omega, state_sc, state_sc%omega)
406 0 : call subcol_state_field_copy(state%pmid, state_sc, state_sc%pmid)
407 0 : call subcol_state_field_copy(state%pdel, state_sc, state_sc%pdel)
408 0 : call subcol_state_field_copy(state%rpdel, state_sc, state_sc%rpdel)
409 0 : call subcol_state_field_copy(state%lnpmid, state_sc, state_sc%lnpmid)
410 0 : call subcol_state_field_copy(state%exner, state_sc, state_sc%exner)
411 0 : call subcol_state_field_copy(state%zm, state_sc, state_sc%zm)
412 0 : call subcol_state_field_copy(state%pint, state_sc, state_sc%pint)
413 0 : call subcol_state_field_copy(state%lnpint, state_sc, state_sc%lnpint)
414 0 : call subcol_state_field_copy(state%zi, state_sc, state_sc%zi)
415 0 : call subcol_state_field_copy(state%lnpmiddry, state_sc, state_sc%lnpmiddry)
416 0 : call subcol_state_field_copy(state%pmiddry, state_sc, state_sc%pmiddry)
417 0 : call subcol_state_field_copy(state%pdeldry, state_sc, state_sc%pdeldry)
418 0 : call subcol_state_field_copy(state%rpdeldry, state_sc, state_sc%rpdeldry)
419 0 : call subcol_state_field_copy(state%pintdry, state_sc, state_sc%pintdry)
420 0 : call subcol_state_field_copy(state%lnpintdry, state_sc, state_sc%lnpintdry)
421 0 : call subcol_state_field_copy(state%q, state_sc, state_sc%q)
422 :
423 0 : end subroutine subcol_state_copy
424 :
425 0 : subroutine subcol_set_subcols(state, tend, nsubcol, state_sc, tend_sc)
426 : !-----------------------------------------------------------------------
427 : ! Propogate nsubcol information to common areas such as state, tend,
428 : ! nsubcol2d, and indcol2d
429 : !-----------------------------------------------------------------------
430 : type(physics_state), intent(in) :: state
431 : type(physics_tend), intent(inout) :: tend
432 : integer, intent(in) :: nsubcol(pcols)
433 : type(physics_state), intent(inout) :: state_sc ! subcolumn state
434 : type(physics_tend), intent(inout) :: tend_sc ! subcolumn tend
435 :
436 0 : call subcol_set_nsubcol(state%lchnk, state%ngrdcol, nsubcol)
437 0 : call subcol_state_copy(state, state_sc)
438 0 : call subcol_tend_copy(tend, state_sc, tend_sc)
439 0 : end subroutine subcol_set_subcols
440 :
441 :
442 0 : subroutine subcol_ptend_copy(ptend, state, ptend_cp)
443 : !-----------------------------------------------------------------------
444 : ! Copy a physics_ptend object into one which has subcolumns
445 : !-----------------------------------------------------------------------
446 : use physics_types, only: physics_ptend_init
447 :
448 : type(physics_ptend), intent(in) :: ptend ! ptend source, dimensioned with grid columns
449 : type(physics_state), intent(in) :: state ! state with subcolumns
450 : type(physics_ptend), intent(out) :: ptend_cp ! copy of ptend, dimensioned with subcolumns
451 :
452 : ! Local Variables
453 : integer :: ncol
454 : integer :: indcol(pcols*psubcols)
455 :
456 : !-----------------------------------------------------------------------
457 :
458 0 : if (ptend%psetcols .ne. pcols) then
459 0 : call endrun('subcol_ptend_copy: ptend must be dimensioned pcols to use this routine')
460 : end if
461 :
462 : call physics_ptend_init(ptend_cp,state%psetcols, ptend%name, ls=ptend%ls, lu=ptend%lu, &
463 0 : lv=ptend%lv, lq=ptend%lq)
464 :
465 0 : ptend_cp%top_level = ptend%top_level
466 0 : ptend_cp%bot_level = ptend%bot_level
467 :
468 0 : ncol = subcol_get_ncol(state%lchnk)
469 0 : call subcol_get_indcol(state%lchnk, indcol)
470 :
471 : ! Copy the grid column data into each of the subcolumns - indcol contains the grid index for each subcolumn
472 0 : if (ptend_cp%ls) then
473 0 : ptend_cp%s(:ncol,:) = ptend%s(indcol(:ncol),:)
474 0 : ptend_cp%hflux_srf(:ncol) = ptend%hflux_srf(indcol(:ncol))
475 0 : ptend_cp%hflux_top(:ncol) = ptend%hflux_top(indcol(:ncol))
476 : end if
477 :
478 0 : if (ptend_cp%lu) then
479 0 : ptend_cp%u(:ncol,:) = ptend%u(indcol(:ncol),:)
480 0 : ptend_cp%taux_srf(:ncol) = ptend%taux_srf(indcol(:ncol))
481 0 : ptend_cp%taux_top(:ncol) = ptend%taux_top(indcol(:ncol))
482 : end if
483 :
484 0 : if (ptend_cp%lv) then
485 0 : ptend_cp%v(:ncol,:) = ptend%v(indcol(:ncol),:)
486 0 : ptend_cp%tauy_srf(:ncol) = ptend%tauy_srf(indcol(:ncol))
487 0 : ptend_cp%tauy_top(:ncol) = ptend%tauy_top(indcol(:ncol))
488 : end if
489 :
490 0 : if (any(ptend_cp%lq(:))) then
491 0 : ptend_cp%q(:ncol,:,:) = ptend%q(indcol(:ncol),:,:)
492 0 : ptend_cp%cflx_srf(:ncol,:) = ptend%cflx_srf(indcol(:ncol),:)
493 0 : ptend_cp%cflx_top(:ncol,:) = ptend%cflx_top(indcol(:ncol),:)
494 : end if
495 :
496 0 : end subroutine subcol_ptend_copy
497 :
498 0 : subroutine subcol_state_hdrinit(state, state_sc)
499 : !-----------------------------------------------------------------------
500 : ! Initialize the subcolumn state header variables
501 : !----------------------------------------------------------------------
502 : type(physics_state), intent(in) :: state
503 : type(physics_state), intent(inout) :: state_sc ! subcolumn state
504 :
505 : integer :: isize
506 : integer :: nsubcol(pcols)
507 :
508 0 : call subcol_get_nsubcol(state%lchnk, nsubcol)
509 0 : isize = state%ngrdcol
510 0 : if (size(nsubcol) < isize) then
511 0 : call endrun('subcol_state_hdrinit error: input nsubcol array must be dimensioned at least as large as state%ngrdcol')
512 : end if
513 :
514 0 : state_sc%ngrdcol = state%ngrdcol
515 0 : state_sc%lchnk = state%lchnk
516 :
517 0 : state_sc%psetcols = pcols*psubcols
518 : ! Set count to a too-large value. It should be correctly initialized in check_energy_timestep_init
519 0 : state_sc%count = pcols*psubcols*2
520 :
521 : ! Set the number of set subcolumns to the total count
522 0 : state_sc%ncol = subcol_get_ncol(state%lchnk)
523 :
524 0 : end subroutine subcol_state_hdrinit
525 :
526 : ! TYPE int,double,real
527 : ! DIMS 1,2
528 0 : subroutine subcol_field_avg_shr_{DIMS}d{TYPE}(field_sc, ngrdcol, lchnk, field, usefilter, useweight)
529 : !-----------------------------------------------------------------------
530 : ! This is the high level subcol field averaging routine which
531 : ! averages a field dimensioned pcols*psubcols to a grid one dimensioned pcols
532 : !
533 : ! Its purpose is to check filter and weight settings and to subset the
534 : ! appropriate subsection of the field array to pass on the averaging routine
535 : !-----------------------------------------------------------------------
536 :
537 : {VTYPE},intent(in) :: field_sc{DIMSTR} ! intent in
538 : integer, intent(in) :: ngrdcol ! # grid cols
539 : integer, intent(in) :: lchnk ! chunk index
540 : {VTYPE}, intent(out) :: field{DIMSTR}
541 : logical, intent(in) :: usefilter
542 : logical, intent(in) :: useweight
543 :
544 : !
545 : ! Local variables
546 : !
547 : integer :: indx, nsubcol, i, j
548 : integer :: nsubcols(pcols)
549 :
550 0 : field = ret_nan_{TYPE}
551 :
552 0 : if (usefilter .and. .not. filter_set) then
553 0 : call endrun('subcol_field_avg_shr error: Cannot specify using filters when none set')
554 : end if
555 :
556 0 : if (useweight .and. .not. weight_set) then
557 0 : call endrun('subcol_field_avg_shr error: Cannot specify using weights when none set')
558 : end if
559 :
560 0 : call subcol_get_nsubcol(lchnk, nsubcols)
561 0 : indx = 1
562 0 : do i = 1, ngrdcol
563 0 : if (nsubcols(i) .gt. 0) then
564 0 : nsubcol = nsubcols(i)
565 : #if ({DIMS} >=2)
566 0 : do j=1,size(field_sc,dim=2)
567 : #endif
568 : #if ({DIMS} == 1)
569 0 : field(i)=subcol_avg_inter(field_sc(indx:indx+nsubcol-1), lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
570 : #endif
571 : #if ({DIMS} == 2)
572 0 : field(i,j)=subcol_avg_inter(field_sc(indx:indx+nsubcol-1,j), lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
573 : #endif
574 : #if ({DIMS} >=2)
575 : end do
576 : #endif
577 0 : indx = indx + nsubcol
578 : end if
579 : end do
580 0 : end subroutine subcol_field_avg_shr_{DIMS}d{TYPE}
581 :
582 :
583 : ! TYPE int,double,real
584 : ! DIMS 1,2
585 0 : subroutine subcol_field_get_firstsubcol_{DIMS}d{TYPE}(field_sc, docheck, ngrdcol, lchnk, field)
586 : !-----------------------------------------------------------------------
587 : ! Retrieve the first subcolumn from a field dimensioned pcols*psubcols
588 : ! and assign to one with pcols, performing optional checking that all other subcolumns are identical
589 : !-----------------------------------------------------------------------
590 :
591 : {VTYPE}, intent(in) :: field_sc{DIMSTR} ! intent in
592 : logical, intent(in) :: docheck ! true=check, false=no check
593 : integer, intent(in) :: ngrdcol ! # grid cols
594 : integer, intent(in) :: lchnk ! chunk index
595 : {VTYPE}, intent(out) :: field{DIMSTR}
596 :
597 : !
598 : ! Local variables
599 : !
600 : integer :: indx, nsubcol, i, j, l
601 : integer :: nsubcols(pcols)
602 :
603 0 : call subcol_get_nsubcol(lchnk, nsubcols)
604 0 : field = 0
605 0 : indx = 1
606 0 : do i = 1, ngrdcol
607 0 : if (nsubcols(i) .gt. 0) then
608 0 : nsubcol = nsubcols(i)
609 : #if ({DIMS}>=2)
610 0 : do j=1,size(field_sc,dim=2)
611 : #endif
612 :
613 : #if ({DIMS} == 1)
614 0 : field(i) = field_sc(indx)
615 : #elif ({DIMS} == 2)
616 0 : field(i,j) = field_sc(indx,j)
617 : #endif
618 0 : if (docheck) then
619 0 : do l=1,nsubcol-1
620 : #if ({DIMS} == 1)
621 0 : if (field_sc(indx) /= field_sc(indx+l)) then
622 : #elif ({DIMS} == 2)
623 0 : if (field_sc(indx,j) /= field_sc(indx+l,j)) then
624 : #endif
625 0 : call endrun('subcol_field_get_firstsubcol error: Not all subcolumn fields are identical')
626 : end if
627 : end do
628 : end if
629 : #if ({DIMS}>=2)
630 : end do
631 : #endif
632 0 : indx = indx + nsubcol
633 : end if
634 : end do
635 0 : end subroutine subcol_field_get_firstsubcol_{DIMS}d{TYPE}
636 :
637 0 : subroutine subcol_ptend_avg_shr(ptend_sc, ngrdcol, lchnk, ptend, usefilter, useweight)
638 : !-----------------------------------------------------------------------
639 : ! Average a subcolumn ptend to a grid ptend
640 : !-----------------------------------------------------------------------
641 :
642 : type(physics_ptend), intent(in) :: ptend_sc ! subcolumn ptend
643 : integer, intent(in) :: ngrdcol ! # grid cols
644 : integer, intent(in) :: lchnk ! chunk index
645 : type(physics_ptend), intent(inout) :: ptend ! grid ptend
646 : logical, intent(in) :: usefilter
647 : logical, intent(in) :: useweight
648 :
649 : !
650 : ! Local variables
651 : !
652 : integer :: indx, i, j, k, nsubcol
653 : integer :: nsubcols(pcols)
654 :
655 0 : if (usefilter .and. .not. filter_set) then
656 0 : call endrun('subcol_ptend_avg_shr error: Cannot specify using filters when none set')
657 : end if
658 :
659 0 : if (useweight .and. .not. weight_set) then
660 0 : call endrun('subcol_ptend_avg_shr error: Cannot specify using weights when none set')
661 : end if
662 :
663 0 : call subcol_get_nsubcol(lchnk, nsubcols)
664 : ! physics_ptend_init has already been called by the master interface
665 0 : if (ptend%ls) then
666 0 : ptend%s(:,:) = 0._r8
667 0 : ptend%hflux_srf(:) = 0._r8
668 0 : ptend%hflux_top(:) = 0._r8
669 0 : indx = 1
670 0 : do i = 1, ngrdcol
671 0 : if (nsubcols(i) > 0) then
672 : nsubcol = nsubcols(i)
673 0 : do j=1,pver
674 0 : ptend%s(i,j)=subcol_avg_inter(ptend_sc%s(indx:indx+nsubcol-1,j), &
675 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
676 : end do
677 0 : ptend%hflux_srf(i) = subcol_avg_inter(ptend_sc%hflux_srf(indx:indx+nsubcol-1), &
678 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
679 0 : ptend%hflux_top(i) = subcol_avg_inter(ptend_sc%hflux_top(indx:indx+nsubcol-1), &
680 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
681 0 : indx = indx + nsubcol
682 : end if
683 : end do
684 : end if
685 0 : if (ptend%lu) then
686 0 : ptend%u(:,:) = 0._r8
687 0 : ptend%taux_srf(:) = 0._r8
688 0 : ptend%taux_top(:) = 0._r8
689 0 : indx = 1
690 0 : do i = 1, ngrdcol
691 0 : if (nsubcols(i) > 0) then
692 : nsubcol = nsubcols(i)
693 0 : do j=1,pver
694 0 : ptend%u(i,j)=subcol_avg_inter(ptend_sc%u(indx:indx+nsubcol-1,j), &
695 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
696 : end do
697 0 : ptend%taux_srf(i) = subcol_avg_inter(ptend_sc%taux_srf(indx:indx+nsubcol-1), &
698 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
699 0 : ptend%taux_top(i) = subcol_avg_inter(ptend_sc%taux_top(indx:indx+nsubcol-1), &
700 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
701 0 : indx = indx + nsubcol
702 : end if
703 : end do
704 : end if
705 0 : if (ptend%lv) then
706 0 : ptend%v(:,:) = 0._r8
707 0 : ptend%tauy_srf(:) = 0._r8
708 0 : ptend%tauy_top(:) = 0._r8
709 0 : indx = 1
710 0 : do i = 1, ngrdcol
711 0 : if (nsubcols(i) > 0) then
712 : nsubcol = nsubcols(i)
713 0 : do j=1,pver
714 0 : ptend%v(i,j)=subcol_avg_inter(ptend_sc%v(indx:indx+nsubcol-1,j), &
715 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
716 : end do
717 0 : ptend%tauy_srf(i) = subcol_avg_inter(ptend_sc%tauy_srf(indx:indx+nsubcol-1), &
718 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
719 0 : ptend%tauy_top(i) = subcol_avg_inter(ptend_sc%tauy_top(indx:indx+nsubcol-1), &
720 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
721 0 : indx = indx + nsubcol
722 : end if
723 : end do
724 : end if
725 0 : if (any(ptend%lq(:))) then
726 0 : ptend%q(:,:,:) = 0._r8
727 0 : ptend%cflx_srf(:,:) = 0._r8
728 0 : ptend%cflx_top(:,:) = 0._r8
729 0 : indx = 1
730 0 : do i = 1, ngrdcol
731 0 : if (nsubcols(i) > 0) then
732 : nsubcol = nsubcols(i)
733 0 : do j=1,pver
734 0 : do k=1, pcnst
735 0 : ptend%q(i,j,k)=subcol_avg_inter(ptend_sc%q(indx:indx+nsubcol-1,j,k), &
736 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
737 : end do
738 : end do
739 0 : do k=1,pcnst
740 0 : ptend%cflx_srf(i,k) = subcol_avg_inter(ptend_sc%cflx_srf(indx:indx+nsubcol-1,k), &
741 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
742 0 : ptend%cflx_top(i,k) = subcol_avg_inter(ptend_sc%cflx_top(indx:indx+nsubcol-1,k), &
743 0 : lchnk, i, indx, indx+nsubcol-1, usefilter, useweight)
744 : end do
745 0 : indx = indx + nsubcol
746 : end if
747 : end do
748 : end if
749 :
750 0 : end subroutine subcol_ptend_avg_shr
751 :
752 0 : subroutine subcol_ptend_get_firstsubcol(ptend_sc, docheck, ngrdcol, lchnk, ptend)
753 : !-----------------------------------------------------------------------
754 : ! Retrieve the first subcolumn from a ptend field dimensioned pcols*psubcols
755 : ! and assign to one with pcols, performing optional check that all other subcolumns are identical
756 : !-----------------------------------------------------------------------
757 :
758 : type(physics_ptend), intent(in) :: ptend_sc ! subcolumn ptend
759 : logical, intent(in) :: docheck ! perform check that all subcolumns match
760 : integer, intent(in) :: ngrdcol ! # grid cols
761 : integer, intent(in) :: lchnk ! chunk index
762 : type(physics_ptend), intent(inout) :: ptend ! grid ptend
763 :
764 : !
765 : ! Local variables
766 : !
767 : integer :: indx, i, j, l
768 : integer :: nsubcols(pcols)
769 :
770 0 : call subcol_get_nsubcol(lchnk, nsubcols)
771 : ! physics_ptend_init has already been called by the master interface
772 0 : if (ptend%ls) then
773 0 : ptend%s(:,:) = 0._r8
774 0 : ptend%hflux_srf(:) = 0._r8
775 0 : ptend%hflux_top(:) = 0._r8
776 0 : indx = 1
777 0 : do i = 1, ngrdcol
778 0 : if (docheck) then
779 0 : do l=1,nsubcols(i)-1
780 0 : if (any(ptend_sc%s(indx,:) /= ptend_sc%s(indx+l,:))) &
781 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%s')
782 0 : if (ptend_sc%hflux_srf(indx) /= ptend_sc%hflux_srf(indx+l)) &
783 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%hflux_srf')
784 0 : if (ptend_sc%hflux_top(indx) /= ptend_sc%hflux_top(indx+l)) &
785 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%hflux_top')
786 : end do
787 : end if
788 0 : ptend%s(i,:) = ptend_sc%s(indx,:)
789 0 : ptend%hflux_srf(i) = ptend_sc%hflux_srf(indx)
790 0 : ptend%hflux_top(i) = ptend_sc%hflux_top(indx)
791 0 : indx = indx + nsubcols(i)
792 : end do
793 : end if
794 0 : if (ptend%lu) then
795 0 : ptend%u(:,:) = 0._r8
796 0 : ptend%taux_srf(:) = 0._r8
797 0 : ptend%taux_top(:) = 0._r8
798 0 : indx = 1
799 0 : do i = 1, ngrdcol
800 0 : do l=1,nsubcols(i)-1
801 0 : if (any(ptend_sc%u(indx,:) /= ptend_sc%u(indx+l,:))) &
802 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%u')
803 0 : if (ptend_sc%taux_srf(indx) /= ptend_sc%taux_srf(indx+l)) &
804 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%taux_srf')
805 0 : if (ptend_sc%taux_top(indx) /= ptend_sc%taux_top(indx+l)) &
806 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%taux_top')
807 : end do
808 0 : ptend%u(i,:) = ptend_sc%u(indx,:)
809 0 : ptend%taux_srf(i) = ptend_sc%taux_srf(indx)
810 0 : ptend%taux_top(i) = ptend_sc%taux_top(indx)
811 0 : indx = indx + nsubcols(i)
812 : end do
813 : end if
814 0 : if (ptend%lv) then
815 0 : ptend%v(:,:) = 0._r8
816 0 : ptend%tauy_srf(:) = 0._r8
817 0 : ptend%tauy_top(:) = 0._r8
818 0 : indx = 1
819 0 : do i = 1, ngrdcol
820 0 : do l=1,nsubcols(i)-1
821 0 : if (any(ptend_sc%v(indx,:) /= ptend_sc%v(indx+l,:))) &
822 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%v')
823 0 : if (ptend_sc%tauy_srf(indx) /= ptend_sc%tauy_srf(indx+l)) &
824 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%tauy_srf')
825 0 : if (ptend_sc%tauy_top(indx) /= ptend_sc%tauy_top(indx+l)) &
826 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%tauy_top')
827 : end do
828 0 : ptend%v(i,:) = ptend_sc%v(indx,:)
829 0 : ptend%tauy_srf(i) = ptend_sc%tauy_srf(indx)
830 0 : ptend%tauy_top(i) = ptend_sc%tauy_top(indx)
831 0 : indx = indx + nsubcols(i)
832 : end do
833 : end if
834 0 : if (any(ptend%lq(:))) then
835 0 : ptend%q(:,:,:) = 0._r8
836 0 : ptend%cflx_srf(:,:) = 0._r8
837 0 : ptend%cflx_top(:,:) = 0._r8
838 0 : indx = 1
839 0 : do i = 1, ngrdcol
840 0 : do l=1,nsubcols(i)-1
841 0 : if (any(ptend_sc%q(indx,:,:) /= ptend_sc%q(indx+l,:,:))) &
842 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%q')
843 0 : if (any(ptend_sc%cflx_srf(indx,:) /= ptend_sc%cflx_srf(indx+l,:))) &
844 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%cflx_srf')
845 0 : if (any(ptend_sc%cflx_top(indx,:) /= ptend_sc%cflx_top(indx+l,:))) &
846 0 : call endrun('subcol_ptend_get_firstsubcol error: Not all subcolumn fields are identical for ptend%cflx_top')
847 : end do
848 0 : ptend%q(i,:,:) = ptend_sc%q(indx,:,:)
849 0 : ptend%cflx_srf(i,:) = ptend_sc%cflx_srf(indx,:)
850 0 : ptend%cflx_top(i,:) = ptend_sc%cflx_top(indx,:)
851 0 : indx = indx + nsubcols(i)
852 : end do
853 : end if
854 :
855 0 : end subroutine subcol_ptend_get_firstsubcol
856 :
857 : ! TYPE int,double,real
858 0 : {VTYPE} function subcol_avg_inter_{TYPE}(vals, lchnk, icol, indx1, indx2, usefilter, useweight) result(avgs)
859 : !------------------------------------------------------------------
860 : ! This function handles the transformation of the usefilter and useweight logicals to passing the
861 : ! actual fields which subcol_avg requires based on the values of the logicals
862 : !------------------------------------------------------------------
863 :
864 : {VTYPE},intent(in) :: vals(:)
865 : integer, intent(in) :: lchnk
866 : integer, intent(in) :: icol
867 : integer, intent(in) :: indx1
868 : integer, intent(in) :: indx2
869 : logical, intent(in) :: usefilter
870 : logical, intent(in) :: useweight
871 :
872 : integer :: nsubcol(pcols)
873 :
874 0 : call subcol_get_nsubcol(lchnk, nsubcol)
875 0 : if (usefilter .and. useweight) then
876 0 : avgs = subcol_avg(vals,nsubcol(icol),filter=filter2d(indx1:indx2,lchnk),weight=weight2d(indx1:indx2,lchnk))
877 0 : else if (useweight) then
878 0 : avgs = subcol_avg(vals,nsubcol(icol),weight=weight2d(indx1:indx2,lchnk))
879 0 : else if (usefilter) then
880 0 : avgs = subcol_avg(vals, nsubcol(icol),filter=filter2d(indx1:indx2,lchnk))
881 : else
882 0 : avgs = subcol_avg(vals,nsubcol(icol))
883 : end if
884 :
885 0 : end function subcol_avg_inter_{TYPE}
886 :
887 : ! TYPE int,double,real
888 0 : {VTYPE} function subcol_avg_{TYPE}(vals, nsubcol, filter, weight) result(avgs)
889 : !------------------------------------------------------------------
890 : ! This function performs the averaging of subcolumn fields, using the optional &
891 : ! filters and weights appropriately
892 : !------------------------------------------------------------------
893 :
894 : {VTYPE}, intent(in) :: vals(:)
895 : integer, intent(in) :: nsubcol
896 : integer, intent(in), optional :: filter(:)
897 : real(r8), intent(in), optional :: weight(:)
898 :
899 : integer :: icnt
900 : {VTYPE} :: fillval
901 :
902 0 : fillval = fillval_{TYPE}
903 :
904 0 : if (present(filter) .and. present(weight)) then
905 0 : if (any(filter==1).and. sum(weight,mask=(filter==1)) /=0 ) then
906 0 : avgs = sum(vals*weight, mask=(filter==1)) / sum(weight, mask=(filter==1))
907 : else
908 : avgs = fillval
909 : end if
910 0 : else if (present(weight)) then
911 0 : if (sum(weight) /=0 ) then
912 0 : avgs = sum(vals*weight) / sum(weight)
913 : else
914 : avgs = fillval
915 : end if
916 0 : else if (present(filter)) then
917 0 : if (any(filter==1)) then
918 0 : icnt = count(filter==1)
919 0 : avgs = sum(vals, mask=(filter==1)) / icnt
920 : else
921 : avgs = fillval
922 : end if
923 0 : else if (nsubcol /= 0) then
924 0 : avgs = sum(vals) / nsubcol
925 : else
926 : avgs = fillval
927 : end if
928 :
929 0 : end function subcol_avg_{TYPE}
930 :
931 0 : subroutine subcol_utils_init_restart(File, hdimids)
932 :
933 : use pio, only: file_desc_t, pio_int, pio_double
934 : use cam_pio_utils, only: cam_pio_def_dim, cam_pio_def_var
935 : use subcol_pack_mod, only: subcol_pack_init_restart
936 :
937 : ! Dummy arguments
938 : type(file_desc_t), intent(inout) :: File
939 : integer, intent(in) :: hdimids(:)
940 :
941 : ! Local variable
942 0 : integer, allocatable :: adimids(:)
943 :
944 0 : call subcol_pack_init_restart(File, hdimids)
945 :
946 : ! Storing filter and weight data even if not being filled by the current
947 : ! subcolumn generator. While these are 2-D arrays, they don't match
948 : ! the grid so we will have to recast them as 3-D
949 :
950 : ! We will need a dimid for the subcol dimension
951 : call cam_pio_def_dim(File, 'psubcols', psubcols, subcol_dimid, &
952 0 : existOK=.true.)
953 0 : allocate(adimids(size(hdimids,1) + 1))
954 0 : adimids(1:size(hdimids)) = hdimids(:)
955 0 : adimids(size(hdimids) + 1) = subcol_dimid
956 :
957 0 : call cam_pio_def_var(File, 'FILTER2D', pio_int, adimids, filter2d_desc)
958 0 : call cam_pio_def_var(File, 'WEIGHT2D', pio_double, adimids, weight2d_desc)
959 :
960 0 : end subroutine subcol_utils_init_restart
961 :
962 0 : subroutine subcol_utils_write_restart(File)
963 : use cam_grid_support, only: cam_grid_write_dist_array
964 0 : use cam_grid_support, only: cam_grid_id, cam_grid_dimensions
965 : use ppgrid, only: begchunk, endchunk
966 : use subcol_pack_mod, only: subcol_unpack
967 : use subcol_pack_mod, only: subcol_pack_write_restart
968 : use pio, only: file_desc_t
969 :
970 : ! Dummy argument
971 : type(file_desc_t), intent(inout) :: File
972 :
973 : ! Local variables
974 : integer :: c
975 : integer :: adimlens(3)
976 : integer :: fdimlens(3)
977 : integer :: frank
978 : integer :: grid_id
979 0 : integer, allocatable :: unpacked_i(:,:,:)
980 0 : real(r8), allocatable :: unpacked_r(:,:,:)
981 : character(len=*), parameter :: subname = 'SUBCOL_UTILS_WRITE_RESTART'
982 :
983 : ! File dimensions
984 0 : grid_id = cam_grid_id('physgrid')
985 0 : call cam_grid_dimensions(grid_id, fdimlens(1:2), frank)
986 :
987 : ! Write nsubcol2d
988 0 : call subcol_pack_write_restart(File, grid_id, fdimlens(1:frank))
989 :
990 : ! filter2d and weight2d are 3-D variables
991 0 : fdimlens(frank + 1) = psubcols
992 0 : frank = frank + 1
993 :
994 : ! Write filter2d
995 0 : adimlens(1) = pcols
996 0 : adimlens(2) = psubcols
997 0 : adimlens(3) = endchunk - begchunk + 1
998 0 : if ((pcols * psubcols) /= size(filter2d, 1)) then
999 0 : call endrun(trim(subname)//": Unsupported size for FILTER2D")
1000 : end if
1001 : ! Unpack filter2d for proper output
1002 0 : allocate(unpacked_i(pcols, psubcols, begchunk:endchunk))
1003 0 : do c = begchunk, endchunk
1004 0 : call subcol_unpack(c, filter2d(:,c), unpacked_i(:,:,c), 0)
1005 : end do
1006 : call cam_grid_write_dist_array(File, grid_id, adimlens, &
1007 0 : fdimlens(1:frank), unpacked_i, filter2d_desc)
1008 0 : deallocate(unpacked_i)
1009 :
1010 : ! Write weight2d
1011 0 : adimlens(1) = pcols
1012 0 : adimlens(2) = psubcols
1013 0 : adimlens(3) = endchunk - begchunk + 1
1014 0 : if ((pcols * psubcols) /= size(weight2d, 1)) then
1015 0 : call endrun(trim(subname)//": Unsupported size for WEIGHT2D")
1016 : end if
1017 : ! Unpack weight2d for proper output
1018 0 : allocate(unpacked_r(pcols, psubcols, begchunk:endchunk))
1019 0 : do c = begchunk, endchunk
1020 0 : call subcol_unpack(c, weight2d(:,c), unpacked_r(:,:,c), 0.0_r8)
1021 : end do
1022 : call cam_grid_write_dist_array(File, grid_id, adimlens, &
1023 0 : fdimlens(1:frank), unpacked_r, weight2d_desc)
1024 0 : deallocate(unpacked_r)
1025 :
1026 0 : end subroutine subcol_utils_write_restart
1027 :
1028 0 : subroutine subcol_utils_read_restart(File)
1029 0 : use pio, only: file_desc_t, pio_inq_varid
1030 : use cam_pio_utils, only: cam_pio_handle_error
1031 : use cam_grid_support, only: cam_grid_id, cam_grid_read_dist_array
1032 : use cam_grid_support, only: cam_grid_dimensions
1033 : use ppgrid, only: begchunk, endchunk
1034 : use subcol_pack_mod, only: subcol_pack
1035 : use subcol_pack_mod, only: subcol_pack_read_restart
1036 :
1037 : ! Dummy argument
1038 : type(file_desc_t), intent(inout) :: File
1039 :
1040 : integer :: ierr, c
1041 : integer :: adimlens(3)
1042 : integer :: fdimlens(3)
1043 : integer :: grid_id
1044 : integer :: frank
1045 0 : integer, allocatable :: unpacked_i(:,:,:)
1046 0 : real(r8), allocatable :: unpacked_r(:,:,:)
1047 : character(len=*), parameter :: subname = 'SUBCOL_UTILS_READ_RESTART'
1048 : integer, allocatable :: nsubcol2d(:,:)
1049 :
1050 0 : call subcol_allocate_internal()
1051 :
1052 : ! File dimensions
1053 0 : grid_id = cam_grid_id('physgrid')
1054 0 : call cam_grid_dimensions(grid_id, fdimlens(1:2), frank)
1055 0 : call subcol_pack_read_restart(File, grid_id, fdimlens(1:frank))
1056 :
1057 0 : ierr = pio_inq_varid(File, 'FILTER2D', filter2d_desc)
1058 0 : call cam_pio_handle_error(ierr, trim(subname)//': FILTER2D not found')
1059 : ! Array dimensions
1060 0 : adimlens(1) = pcols
1061 0 : adimlens(2) = psubcols
1062 0 : adimlens(3) = endchunk - begchunk + 1
1063 0 : if ((pcols * psubcols) /= size(filter2d, 1)) then
1064 0 : call endrun(trim(subname)//": Unsupported size for FILTER2D")
1065 : end if
1066 0 : allocate(unpacked_i(pcols, psubcols, begchunk:endchunk))
1067 : ! File dimensions (good for both filter2d and weight2d)
1068 0 : frank = frank + 1
1069 0 : fdimlens(frank) = psubcols
1070 : call cam_grid_read_dist_array(File, grid_id, adimlens, &
1071 0 : fdimlens(1:frank), unpacked_i, filter2d_desc)
1072 : ! Pack filter2d for proper output
1073 0 : do c = begchunk, endchunk
1074 0 : call subcol_pack(c, unpacked_i(:,:,c), filter2d(:,c))
1075 : end do
1076 0 : deallocate(unpacked_i)
1077 :
1078 0 : ierr = pio_inq_varid(File, 'WEIGHT2D', weight2d_desc)
1079 0 : adimlens(1) = pcols
1080 0 : adimlens(2) = psubcols
1081 0 : adimlens(3) = endchunk - begchunk + 1
1082 0 : if ((pcols * psubcols) /= size(weight2d, 1)) then
1083 0 : call endrun(trim(subname)//": Unsupported size for WEIGHT2D")
1084 : end if
1085 0 : allocate(unpacked_r(pcols, psubcols, begchunk:endchunk))
1086 0 : call cam_pio_handle_error(ierr, trim(subname)//': WEIGHT2D not found')
1087 : call cam_grid_read_dist_array(File, grid_id, adimlens, &
1088 0 : fdimlens(1:frank), unpacked_r, weight2d_desc)
1089 : ! Pack weight2d for proper output
1090 0 : do c = begchunk, endchunk
1091 0 : call subcol_pack(c, unpacked_r(:,:,c), weight2d(:,c))
1092 : end do
1093 0 : deallocate(unpacked_r)
1094 :
1095 0 : end subroutine subcol_utils_read_restart
1096 :
1097 : end module subcol_utils
|