Line data Source code
1 : module subcol
2 : !---------------------------------------------------------------------------
3 : ! Purpose:
4 : !
5 : ! Provides the infrastructure for handling the relationship between
6 : ! grid-box-averaged (GBA) fields and Sub-column (subcol) fields
7 : !
8 : ! Different methods for generating and averaging sub-column fields are
9 : ! called sub-column schemes and are designated with the 'subcol_scheme'
10 : ! namelist variable.
11 : !
12 : ! This module provides several public interfaces (see below) which operate
13 : ! based on the scheme. In order to implement a new scheme, you need to
14 : ! follow the following steps:
15 : ! I) Implement a sub-column scheme in a separate module (subcol_<scheme>)
16 : ! I.1) Implement a subcol_register_<scheme> function to register any scheme-
17 : ! specific fields. Any pbuf_add_field and/or pbuf_register_subcol calls
18 : ! need to go here.
19 : ! This step is optional
20 : ! I.2) Implement a subcol_init_<scheme> function to initialize any scheme-
21 : ! specific variables or fields.
22 : ! This step is optional
23 : ! I.3) Implement a subcol_gen_<scheme> function to generate the appropriate
24 : ! subcol fields based on the existing GBA fields and other scheme data.
25 : ! I.4) Implement a subcol_field_avg_<scheme> function to average subcol fields
26 : ! back into the appropriate GBA fields.
27 : ! This step is optional
28 : ! I.5) Implement a subcol_ptend_avg_<scheme> function to average the
29 : ! sub-column ptend to a grid ptend
30 : ! This step is optional
31 : ! II) Add necessary cases in the master subcol module (this file)
32 : ! II.1) Add a case for your scheme name in subcol_register if you are calling
33 : ! your own subcol_<scheme> registration function.
34 : ! II.2) Add a case for your scheme name in subcol_init if you are calling your
35 : ! own subcol_<scheme> initialization function.
36 : ! II.3) Add a case for your scheme name in subcol_gen and call your
37 : ! subcol_<scheme> subcol generator.
38 : ! II.4) Add a case for your scheme name in subcol_field_avg if you are calling your
39 : ! own subcol_<scheme> field-averaging function.
40 : ! II.5) Add a case for your scheme name in subcol_ptend_avg if you are calling your
41 : ! own subcol_<scheme> ptend-averaging function.
42 : !
43 : ! New schemes should be implemented in a separate file which is used by
44 : ! this module (and thus may not 'use' any subcol module variable, function,
45 : ! or subroutine).
46 : !
47 : !---------------------------------------------------------------------------
48 :
49 : use shr_kind_mod, only: r8=>shr_kind_r8, r4=>shr_kind_r4, i4=>shr_kind_i4
50 : use physics_types, only: physics_state, physics_tend, physics_ptend
51 : use ppgrid, only: pcols, psubcols, pver, pverp
52 : use cam_abortutils, only: endrun
53 : use subcol_utils, only: subcol_field_avg_shr, subcol_ptend_avg_shr, &
54 : subcol_field_get_firstsubcol, subcol_ptend_get_firstsubcol, &
55 : is_filter_set, is_weight_set
56 : use subcol_tstcp , only: subcol_gen_tstcp, subcol_register_tstcp, subcol_field_avg_tstcp, subcol_ptend_avg_tstcp
57 :
58 : ! CloudObj, SILHS and vamp are currently being developed
59 : ! References are being left in for convenience and demonstration purposes
60 : ! use subcol_CloudObj, only: cloudobj_scheme_name, subcol_register_CloudObj, &
61 : ! subcol_init_CloudObj, subcol_gen_CloudObj
62 : ! use subcol_CloudObj, only: subcol_ptend_avg_CloudObj
63 : use subcol_SILHS, only: subcol_register_SILHS, subcol_init_SILHS, &
64 : subcol_gen_SILHS, subcol_ptend_avg_SILHS
65 : ! use subcol_vamp, only: subcol_gen_vamp, subcol_register_vamp, subcol_init_vamp
66 :
67 : implicit none
68 :
69 : private
70 : save
71 :
72 : !! Public interface functions which implement a sub-column scheme
73 : public :: subcol_register ! Read scheme from namelist and initialize any scheme-global variables
74 : public :: subcol_init ! Initialize any variables or fields specific to the active scheme
75 : public :: subcol_gen ! Generate subcol fields from GBA fields
76 : public :: subcol_field_avg ! Average subcol fields back into GBA fields
77 : public :: subcol_ptend_avg ! Average sub-column ptend to grid ptend
78 : public :: subcol_readnl ! Namelist reader for subcolumns
79 : public :: subcol_init_restart ! Initialize restart with subcolumn specific fields
80 : public :: subcol_read_restart ! Read subcolumn specific fields from restart
81 : public :: subcol_write_restart ! Write subcolumn specific fields for restart
82 :
83 :
84 : interface subcol_field_avg
85 : module procedure subcol_field_avg_1dr
86 : module procedure subcol_field_avg_1di
87 : module procedure subcol_field_avg_2dr
88 : end interface
89 :
90 : contains
91 1536 : subroutine subcol_readnl(nlfile)
92 : use subcol_utils, only: subcol_get_scheme, subcol_utils_readnl
93 : use subcol_tstcp, only: subcol_readnl_tstcp
94 : use subcol_SILHS, only: subcol_readnl_SILHS
95 : ! use subcol_vamp, only: subcol_readnl_vamp
96 :
97 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
98 :
99 : !
100 : ! Local variables
101 : !
102 : character(len=16) :: subcol_scheme_init ! Name of subcolumn schem
103 : !-----------------------------------------------------------------------------
104 :
105 1536 : call subcol_utils_readnl(nlfile)
106 1536 : subcol_scheme_init = subcol_get_scheme()
107 :
108 3072 : select case(trim(subcol_scheme_init))
109 : case('tstcp')
110 0 : call subcol_readnl_tstcp(nlfile)
111 : case ('SILHS')
112 0 : call subcol_readnl_SILHS(nlfile)
113 : ! case (cloudobj_scheme_name)
114 : ! call subcol_readnl_CloudObj(nlfile)
115 : ! case ('vamp')
116 : ! call subcol_readnl_vamp(nlfile)
117 : case ('off')
118 : ! No namelist for off
119 : case default
120 3072 : call endrun('subcol_register error: unsupported subcol_scheme specified')
121 : end select
122 :
123 1536 : end subroutine subcol_readnl
124 :
125 1536 : subroutine subcol_register()
126 1536 : use phys_control, only: phys_getopts
127 : use physics_buffer, only: pbuf_add_field, dtype_i4
128 : use subcol_utils, only: subcol_get_scheme
129 :
130 0 : select case(subcol_get_scheme())
131 : case('tstcp')
132 0 : call subcol_register_tstcp()
133 : case ('SILHS')
134 0 : call subcol_register_SILHS()
135 : ! case (cloudobj_scheme_name)
136 : ! call subcol_register_CloudObj()
137 : ! case ('vamp')
138 : ! call subcol_register_vamp()
139 : case ('off')
140 : ! No registration called
141 : case default
142 1536 : call endrun('subcol_register error: unsupported subcol_scheme specified')
143 : end select
144 :
145 1536 : end subroutine subcol_register
146 :
147 0 : subroutine subcol_init_restart(File, hdimids)
148 1536 : use subcol_utils, only: subcol_utils_init_restart
149 : use pio, only: file_desc_t
150 :
151 : type(file_desc_t),intent(inout) :: File
152 : integer ,intent(in) :: hdimids(:)
153 :
154 0 : call subcol_utils_init_restart(File, hdimids)
155 : ! Put scheme-specific calls here (in select statement)
156 :
157 0 : end subroutine subcol_init_restart
158 :
159 0 : subroutine subcol_write_restart(File)
160 0 : use subcol_utils, only: subcol_utils_write_restart
161 : use pio, only: file_desc_t
162 :
163 : type(file_desc_t), intent(inout) :: File
164 :
165 0 : call subcol_utils_write_restart(File)
166 : ! Put scheme-specific calls here (in select statement)
167 :
168 0 : end subroutine subcol_write_restart
169 :
170 0 : subroutine subcol_read_restart(File)
171 0 : use subcol_utils, only: subcol_utils_read_restart
172 : use pio, only: file_desc_t
173 :
174 : type(file_desc_t), intent(inout) :: File
175 :
176 0 : call subcol_utils_read_restart(File)
177 : ! Put scheme-specific calls here (in select statement)
178 :
179 0 : end subroutine subcol_read_restart
180 :
181 1536 : subroutine subcol_init(pbuf2d, subcol_scheme_in)
182 0 : use physics_buffer, only: physics_buffer_desc
183 : use cam_history_support, only: add_hist_coord
184 : use subcol_utils, only: subcol_utils_init, subcol_get_scheme
185 : use time_manager, only: is_first_step, is_first_restart_step
186 :
187 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
188 : character(len=*), optional, intent(in) :: subcol_scheme_in ! Name of subcolumn generator
189 :
190 : !
191 : ! Local variables
192 : !
193 : character(len=16) :: subcol_scheme_init ! Name of subcolumn scheme
194 :
195 : ! Set the subcol_scheme_gen to the one passed in , otherwise use the module scheme read from the namelist
196 1536 : if (present(subcol_scheme_in)) then
197 0 : subcol_scheme_init = trim(subcol_scheme_in)
198 : else
199 : ! By default, use the module scheme read from the namelist
200 1536 : subcol_scheme_init = subcol_get_scheme()
201 : end if
202 :
203 1536 : if (is_first_step() .and. .not. is_first_restart_step()) then
204 : ! Initialize the subcol utility data only at the beginning of the run and not at restart
205 768 : call subcol_utils_init(subcol_scheme_init)
206 : end if
207 :
208 : ! Set the psubcols history coordinate for output
209 1536 : if (trim(subcol_scheme_init) /= 'off') then
210 0 : call add_hist_coord('psubcols', psubcols, 'Subcolumn Index')
211 : end if
212 :
213 : ! Call the appropriate subcol init method
214 3072 : select case(trim(subcol_scheme_init))
215 : case ('tstcp')
216 : ! none needed for this scheme
217 : case ('SILHS')
218 0 : call subcol_init_SILHS(pbuf2d)
219 : ! case (cloudobj_scheme_name)
220 : ! call subcol_init_CloudObj(pbuf2d)
221 : ! case ('vamp')
222 : ! call subcol_init_vamp()
223 : case ('off')
224 : ! No initialization needed for off
225 : case default
226 3072 : call endrun('subcol_init error: unsupported subcol_scheme specified')
227 : end select
228 1536 : end subroutine subcol_init
229 :
230 0 : subroutine subcol_gen(state, tend, state_sc, tend_sc, pbuf, subcol_scheme_in)
231 1536 : use physics_buffer, only: physics_buffer_desc
232 : use subcol_utils, only: subcol_get_scheme
233 :
234 :
235 : !-----------------------------------
236 : ! sub-column generator
237 : !-----------------------------------
238 : type(physics_state), intent(inout) :: state
239 : type(physics_tend), intent(inout) :: tend
240 : type(physics_state), intent(inout) :: state_sc ! sub-column state
241 : type(physics_tend), intent(inout) :: tend_sc ! sub-column tend
242 : type(physics_buffer_desc), pointer :: pbuf(:)
243 : character(len=*), optional, intent(in) :: subcol_scheme_in ! Name of subcolumn generator
244 :
245 : !
246 : ! Local variables
247 : !
248 : character(len=16) :: subcol_scheme_gen ! Name of subcolumn scheme
249 :
250 : ! Set the subcol_scheme_gen to the one passed in , otherwise use the module scheme read from the namelist
251 0 : if (present(subcol_scheme_in)) then
252 0 : subcol_scheme_gen = trim(subcol_scheme_in)
253 : else
254 0 : subcol_scheme_gen = subcol_get_scheme()
255 : end if
256 :
257 0 : if (subcol_scheme_gen /= 'off') then
258 0 : if (.not. allocated(state_sc%lat)) then
259 0 : call endrun('subcol_gen error: state_sc must be allocated before calling subcol_gen')
260 : end if
261 : end if
262 :
263 0 : if (state_sc%psetcols /= (pcols * psubcols)) then
264 0 : call endrun('subcol_gen error: state_sc%psetcols must be pcols * psubcols')
265 : end if
266 :
267 0 : select case(trim(subcol_scheme_gen))
268 : case('tstcp')
269 0 : call subcol_gen_tstcp(state, tend, state_sc, tend_sc, pbuf)
270 : case ('SILHS')
271 0 : call subcol_gen_SILHS(state, tend, state_sc, tend_sc, pbuf)
272 : ! case (cloudobj_scheme_name)
273 : ! call subcol_gen_CloudObj(state, tend, state_sc, tend_sc, pbuf)
274 : ! case ('vamp')
275 : ! call subcol_gen_vamp(state, tend, state_sc, tend_sc, pbuf)
276 : case default
277 0 : call endrun('subcol_gen error: unsupported subcol_scheme specified')
278 : end select
279 :
280 0 : end subroutine subcol_gen
281 :
282 0 : subroutine subcol_field_avg_1dr (field_sc, ngrdcol, lchnk, field, subcol_scheme_in)
283 0 : use physics_buffer, only: physics_buffer_desc
284 : use subcol_utils, only: subcol_get_scheme
285 :
286 : !-----------------------------------
287 : ! Average the subcolumns dimension (pcols*psubcols) to the grid dimension (pcols)
288 : !-----------------------------------
289 :
290 : real(r8), intent(in) :: field_sc(:) ! intent in
291 : integer, intent(in) :: ngrdcol ! # grid cols
292 : integer, intent(in) :: lchnk ! chunk index
293 : real(r8), intent(out) :: field(:)
294 : character(len=*), optional, intent(in) :: subcol_scheme_in ! Name of subcolumn generator
295 :
296 : !
297 : ! Local variables
298 : !
299 : character(len=16) :: subcol_scheme_avg ! Name of subcolumn scheme
300 :
301 0 : if (present(subcol_scheme_in)) then
302 0 : subcol_scheme_avg = trim(subcol_scheme_in)
303 : else
304 : ! By default, use the module scheme read from the namelist
305 0 : subcol_scheme_avg = subcol_get_scheme()
306 : end if
307 :
308 0 : if (size(field_sc,dim=1) .ne. pcols*psubcols) then
309 0 : call endrun('subcol_field_avg error: only fields with first dimension pcols*psubcols may use this routine')
310 : end if
311 :
312 0 : select case(trim(subcol_scheme_avg))
313 : ! Example of specialized averaging for specific subcolumn scheme
314 : case ('tstcp')
315 0 : call subcol_field_avg_tstcp(field_sc, ngrdcol, lchnk, field)
316 : ! Unless specialized averaging is needed, most subcolumn schemes will be handled by the default
317 : ! If filters and/or weights have been set, they are automatically used by this averager
318 : case default
319 0 : call subcol_field_avg_shr(field_sc, ngrdcol, lchnk, field, is_filter_set(), is_weight_set())
320 : end select
321 :
322 0 : end subroutine subcol_field_avg_1dr
323 :
324 0 : subroutine subcol_field_avg_1di (field_sc, ngrdcol, lchnk, field, subcol_scheme_in)
325 0 : use physics_buffer, only: physics_buffer_desc
326 : use subcol_utils, only: subcol_get_scheme
327 :
328 : !-----------------------------------
329 : ! Average the subcolumns dimension (pcols*psubcols) to the grid dimension (pcols)
330 : !-----------------------------------
331 :
332 : integer, intent(in) :: field_sc(:) ! intent in
333 : integer, intent(in) :: ngrdcol ! # grid cols
334 : integer, intent(in) :: lchnk ! chunk index
335 : integer, intent(out) :: field(:)
336 : character(len=*), optional, intent(in) :: subcol_scheme_in ! Name of subcolumn generator
337 :
338 : !
339 : ! Local variables
340 : !
341 : character(len=16) :: subcol_scheme_avg ! Name of subcolumn scheme
342 :
343 0 : if (present(subcol_scheme_in)) then
344 0 : subcol_scheme_avg = trim(subcol_scheme_in)
345 : else
346 : ! By default, use the module scheme read from the namelist
347 0 : subcol_scheme_avg = subcol_get_scheme()
348 : end if
349 :
350 0 : if (size(field_sc,dim=1) .ne. pcols*psubcols) then
351 0 : call endrun('subcol_field_avg error: only fields with first dimension pcols*psubcols may use this routine')
352 : end if
353 :
354 0 : select case(trim(subcol_scheme_avg))
355 : ! Example of specialized averaging for specific subcolumn scheme
356 : case ('tstcp')
357 0 : call subcol_field_avg_tstcp(field_sc, ngrdcol, lchnk, field)
358 : ! Unless specialized averaging is needed, most subcolumn schemes will be handled by the default
359 : ! If filters and/or weights have been set, they are automatically used by this averager
360 : case default
361 0 : call subcol_field_avg_shr(field_sc, ngrdcol, lchnk, field, is_filter_set(), is_weight_set())
362 : end select
363 :
364 0 : end subroutine subcol_field_avg_1di
365 :
366 0 : subroutine subcol_field_avg_2dr (field_sc, ngrdcol, lchnk, field, subcol_scheme_in)
367 0 : use physics_buffer, only: physics_buffer_desc
368 : use subcol_utils, only: subcol_get_scheme
369 :
370 : !-----------------------------------
371 : ! Average the subcolumns dimension (pcols*psubcols) to the grid dimension (pcols)
372 : !-----------------------------------
373 :
374 : real(r8), intent(in) :: field_sc(:,:) ! intent in
375 : integer, intent(in) :: ngrdcol ! # grid cols
376 : integer, intent(in) :: lchnk ! chunk index
377 : real(r8), intent(out) :: field(:,:)
378 : character(len=*), optional, intent(in) :: subcol_scheme_in ! Name of subcolumn generator
379 :
380 : !
381 : ! Local variables
382 : !
383 : character(len=16) :: subcol_scheme_avg ! Name of subcolumn scheme
384 :
385 0 : if (present(subcol_scheme_in)) then
386 0 : subcol_scheme_avg = trim(subcol_scheme_in)
387 : else
388 : ! By default, use the module scheme read from the namelist
389 0 : subcol_scheme_avg = subcol_get_scheme()
390 : end if
391 :
392 0 : if (size(field_sc,dim=1) .ne. pcols*psubcols) then
393 0 : call endrun('subcol_field_avg error: only fields with first dimension pcols*psubcols may use this routine')
394 : end if
395 :
396 0 : select case(trim(subcol_scheme_avg))
397 : ! Example of specialized averaging for specific subcolumn scheme
398 : case ('tstcp')
399 0 : call subcol_field_avg_tstcp(field_sc, ngrdcol, lchnk, field)
400 : ! Unless specialized averaging is needed, most subcolumn schemes will be handled with the default
401 : ! If filters and/or weights have been set, they are automatically used by this averager
402 : case default
403 0 : call subcol_field_avg_shr(field_sc, ngrdcol, lchnk, field, is_filter_set(), is_weight_set())
404 : end select
405 :
406 0 : end subroutine subcol_field_avg_2dr
407 :
408 0 : subroutine subcol_ptend_avg(ptend_sc, ngrdcol, lchnk, ptend, subcol_scheme_in)
409 0 : use physics_buffer, only: physics_buffer_desc
410 : use physics_types, only: physics_ptend_init
411 : use subcol_utils, only: subcol_get_scheme
412 :
413 : !-----------------------------------------------------------------------
414 : ! Average a sub-column ptend to a grid ptend
415 : !-----------------------------------------------------------------------
416 :
417 : type(physics_ptend), intent(in) :: ptend_sc ! sub-column ptend
418 : integer, intent(in) :: ngrdcol ! # grid cols
419 : integer, intent(in) :: lchnk ! chunk index
420 : type(physics_ptend), intent(inout) :: ptend ! grid ptend
421 : character(len=*), optional, intent(in) :: subcol_scheme_in ! Name of subcolumn generator
422 :
423 : !
424 : ! Local variables
425 : !
426 : character(len=16) :: subcol_scheme_avg ! Name of subcolumn scheme
427 : integer :: indx, i, j
428 :
429 0 : if (present(subcol_scheme_in)) then
430 0 : subcol_scheme_avg = trim(subcol_scheme_in)
431 : else
432 : ! By default, use the module scheme read from the namelist
433 0 : subcol_scheme_avg = subcol_get_scheme()
434 : end if
435 :
436 : !-----------------------------------------------------------------------
437 :
438 : call physics_ptend_init(ptend, pcols, name=ptend_sc%name, ls=ptend_sc%ls, lu=ptend_sc%lu, &
439 0 : lv=ptend_sc%lv, lq=ptend_sc%lq)
440 :
441 0 : select case(trim(subcol_scheme_avg))
442 : case ('tstcp')
443 0 : call subcol_ptend_avg_tstcp(ptend_sc, ngrdcol, lchnk, ptend)
444 : case ('SILHS')
445 0 : call subcol_ptend_avg_SILHS(ptend_sc, ngrdcol, lchnk, ptend)
446 : case default
447 : ! If filters and/or weights have been set, they are automatically used by this averager
448 0 : call subcol_ptend_avg_shr(ptend_sc, ngrdcol, lchnk, ptend, is_filter_set(), is_weight_set())
449 : end select
450 :
451 0 : end subroutine subcol_ptend_avg
452 :
453 : end module subcol
454 :
455 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
456 :
457 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
458 : !!
459 : !! subcol_field_avg_handler is an external routine used by outfld
460 : !! It is outside the module to prevent a dependency loop
461 : !!
462 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
463 : !! MAKE SURE TO UPDATE THE INTERFACE IN OUTFLD IF THIS INTERFACE IS CHANGED
464 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
465 :
466 0 : subroutine subcol_field_avg_handler(idim, field_in, c, field_out)
467 : use shr_kind_mod, only: r8 => shr_kind_r8
468 : use ppgrid, only: pcols, psubcols
469 : use phys_grid, only: get_ncols_p
470 : use subcol, only: subcol_field_avg
471 :
472 : implicit none
473 :
474 : !! Dummy arguments
475 : integer, intent(in) :: idim
476 : real(r8), intent(in) :: field_in(idim, *)
477 : integer, intent(in) :: c
478 : real(r8), intent(inout) :: field_out(:,:)
479 :
480 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
481 :
482 : !! Local variables
483 0 : real(r8), allocatable :: field_sc(:,:)
484 : integer :: i, j, ngrdcol, dim2
485 :
486 0 : dim2 = size(field_out, 2)
487 0 : allocate(field_sc(pcols*psubcols, dim2))
488 :
489 0 : do j = 1, dim2
490 0 : do i = 1, idim
491 0 : field_sc(i, j) = field_in(i, j)
492 : end do
493 : end do
494 0 : if (idim < (pcols * psubcols)) then
495 0 : field_sc(idim+1:pcols*psubcols,:) = 0.0_r8
496 : end if
497 :
498 0 : ngrdcol = get_ncols_p(c)
499 0 : call subcol_field_avg(field_sc, ngrdcol, c, field_out)
500 :
501 0 : deallocate(field_sc)
502 0 : end subroutine subcol_field_avg_handler
|