Line data Source code
1 : module stepon
2 :
3 : use shr_kind_mod, only: r8 => shr_kind_r8
4 : use spmd_utils, only: iam, mpicom, masterproc
5 : use ppgrid, only: begchunk, endchunk
6 :
7 : use physics_types, only: physics_state, physics_tend
8 : use dyn_comp, only: dyn_import_t, dyn_export_t
9 :
10 : use perf_mod, only: t_startf, t_stopf, t_barrierf
11 : use cam_abortutils, only: endrun
12 :
13 : use parallel_mod, only: par
14 : use dimensions_mod, only: np, npsq, nlev, nelemd
15 :
16 : use aerosol_properties_mod, only: aerosol_properties
17 : use aerosol_state_mod, only: aerosol_state
18 : use microp_aero, only: aerosol_state_object, aerosol_properties_object
19 : use scamMod, only: use_iop, doiopupdate, single_column, &
20 : setiopupdate, readiopdata
21 : use se_single_column_mod, only: scm_setfield, iop_broadcast
22 : use dyn_grid, only: hvcoord
23 : use time_manager, only: get_step_size, is_first_restart_step
24 : use cam_history, only: outfld, write_camiop, addfld, add_default, horiz_only
25 : use cam_history, only: write_inithist, hist_fld_active, fieldname_len
26 :
27 : implicit none
28 : private
29 : save
30 :
31 : public stepon_init
32 : public stepon_run1
33 : public stepon_run2
34 : public stepon_run3
35 : public stepon_final
36 :
37 : class(aerosol_properties), pointer :: aero_props_obj => null()
38 : logical :: aerosols_transported = .false.
39 : logical :: iop_update_phase1
40 :
41 : !=========================================================================================
42 : contains
43 : !=========================================================================================
44 :
45 372480 : subroutine stepon_init(dyn_in, dyn_out )
46 :
47 : use constituents, only: pcnst, cnst_name, cnst_longname
48 : use dimensions_mod, only: fv_nphys, cnst_name_gll, cnst_longname_gll, qsize
49 :
50 : ! arguments
51 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
52 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
53 :
54 : ! local variables
55 : integer :: m, m_cnst
56 : !----------------------------------------------------------------------------
57 : ! These fields on dynamics grid are output before the call to d_p_coupling.
58 10752 : do m_cnst = 1, qsize
59 0 : call addfld(trim(cnst_name_gll(m_cnst))//'_gll', (/ 'lev' /), 'I', 'kg/kg', &
60 18432 : trim(cnst_longname_gll(m_cnst)), gridname='GLL')
61 0 : call addfld(trim(cnst_name_gll(m_cnst))//'dp_gll', (/ 'lev' /), 'I', 'kg/kg', &
62 19968 : trim(cnst_longname_gll(m_cnst))//'*dp', gridname='GLL')
63 : end do
64 3072 : call addfld('U_gll' ,(/ 'lev' /), 'I', 'm/s ','U wind on gll grid',gridname='GLL')
65 3072 : call addfld('V_gll' ,(/ 'lev' /), 'I', 'm/s ','V wind on gll grid',gridname='GLL')
66 3072 : call addfld('T_gll' ,(/ 'lev' /), 'I', 'K ' ,'T on gll grid' ,gridname='GLL')
67 3072 : call addfld('dp_ref_gll' ,(/ 'lev' /), 'I', ' ' ,'dp dry / dp_ref on gll grid' ,gridname='GLL')
68 1536 : call addfld('PSDRY_gll' ,horiz_only , 'I', 'Pa ' ,'psdry on gll grid' ,gridname='GLL')
69 1536 : call addfld('PS_gll' ,horiz_only , 'I', 'Pa ' ,'ps on gll grid' ,gridname='GLL')
70 1536 : call addfld('PHIS_gll' ,horiz_only , 'I', 'Pa ' ,'PHIS on gll grid' ,gridname='GLL')
71 :
72 : ! Fields for initial condition files
73 3072 : call addfld('U&IC', (/ 'lev' /), 'I', 'm/s', 'Zonal wind', gridname='GLL' )
74 3072 : call addfld('V&IC', (/ 'lev' /), 'I', 'm/s', 'Meridional wind',gridname='GLL' )
75 : ! Don't need to register U&IC V&IC as vector components since we don't interpolate IC files
76 1536 : call add_default('U&IC',0, 'I')
77 1536 : call add_default('V&IC',0, 'I')
78 :
79 1536 : call addfld('PS&IC', horiz_only, 'I', 'Pa', 'Surface pressure', gridname='GLL')
80 3072 : call addfld('T&IC', (/ 'lev' /), 'I', 'K', 'Temperature', gridname='GLL')
81 1536 : call add_default('PS&IC', 0, 'I')
82 1536 : call add_default('T&IC', 0, 'I')
83 :
84 64512 : do m_cnst = 1,pcnst
85 62976 : call addfld(trim(cnst_name(m_cnst))//'&IC', (/ 'lev' /), 'I', 'kg/kg', &
86 188928 : trim(cnst_longname(m_cnst)), gridname='GLL')
87 64512 : call add_default(trim(cnst_name(m_cnst))//'&IC', 0, 'I')
88 : end do
89 :
90 : ! get aerosol properties
91 1536 : aero_props_obj => aerosol_properties_object()
92 :
93 1536 : if (associated(aero_props_obj)) then
94 : ! determine if there are transported aerosol contistuents
95 1536 : aerosols_transported = aero_props_obj%number_transported()>0
96 : end if
97 :
98 1536 : end subroutine stepon_init
99 :
100 : !=========================================================================================
101 :
102 741888 : subroutine stepon_run1( dtime_out, phys_state, phys_tend, &
103 : pbuf2d, dyn_in, dyn_out )
104 :
105 : use dp_coupling, only: d_p_coupling
106 : use physics_buffer, only: physics_buffer_desc
107 :
108 : use se_dyn_time_mod,only: tstep ! dynamics timestep
109 :
110 : real(r8), intent(out) :: dtime_out ! Time-step
111 : type(physics_state), intent(inout) :: phys_state(begchunk:endchunk)
112 : type(physics_tend), intent(inout) :: phys_tend(begchunk:endchunk)
113 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
114 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
115 : type (physics_buffer_desc), pointer :: pbuf2d(:,:)
116 : !----------------------------------------------------------------------------
117 :
118 : integer :: c
119 : class(aerosol_state), pointer :: aero_state_obj
120 370944 : nullify(aero_state_obj)
121 :
122 370944 : dtime_out = get_step_size()
123 :
124 370944 : if (iam < par%nprocs) then
125 370944 : if (tstep <= 0) call endrun('stepon_run1: bad tstep')
126 370944 : if (dtime_out <= 0) call endrun('stepon_run1: bad dtime')
127 :
128 : ! write diagnostic fields on gll grid and initial file
129 370944 : call diag_dynvar_ic(dyn_out%elem, dyn_out%fvm)
130 : end if
131 :
132 : ! Determine whether it is time for an IOP update;
133 : ! doiopupdate set to true if model time step > next available IOP
134 :
135 :
136 370944 : if (use_iop .and. masterproc) then
137 0 : call setiopupdate
138 : end if
139 :
140 370944 : if (single_column) then
141 :
142 : ! If first restart step then ensure that IOP data is read
143 0 : if (is_first_restart_step()) then
144 0 : if (masterproc) call readiopdata( hvcoord%hyam, hvcoord%hybm, hvcoord%hyai, hvcoord%hybi, hvcoord%ps0 )
145 0 : call iop_broadcast()
146 : endif
147 :
148 0 : iop_update_phase1 = .true.
149 0 : if ((is_first_restart_step() .or. doiopupdate) .and. masterproc) then
150 0 : call readiopdata( hvcoord%hyam, hvcoord%hybm, hvcoord%hyai, hvcoord%hybi, hvcoord%ps0 )
151 : endif
152 0 : call iop_broadcast()
153 :
154 0 : call scm_setfield(dyn_out%elem,iop_update_phase1)
155 : endif
156 :
157 370944 : call t_barrierf('sync_d_p_coupling', mpicom)
158 370944 : call t_startf('d_p_coupling')
159 : ! Move data into phys_state structure.
160 370944 : call d_p_coupling(phys_state, phys_tend, pbuf2d, dyn_out )
161 370944 : call t_stopf('d_p_coupling')
162 :
163 : !----------------------------------------------------------
164 : ! update aerosol state object from CAM physics state constituents
165 : !----------------------------------------------------------
166 370944 : if (aerosols_transported) then
167 :
168 0 : do c = begchunk,endchunk
169 0 : aero_state_obj => aerosol_state_object(c)
170 : ! pass number mass or number mixing ratios of aerosol constituents
171 : ! to aerosol state object
172 0 : call aero_state_obj%set_transported(phys_state(c)%q)
173 : end do
174 :
175 : end if
176 :
177 370944 : end subroutine stepon_run1
178 :
179 : !=========================================================================================
180 :
181 738816 : subroutine stepon_run2(phys_state, phys_tend, dyn_in, dyn_out)
182 :
183 370944 : use dp_coupling, only: p_d_coupling
184 : use dyn_grid, only: TimeLevel
185 :
186 : use se_dyn_time_mod, only: TimeLevel_Qdp
187 : use control_mod, only: qsplit
188 : use prim_advance_mod, only: tot_energy_dyn
189 :
190 :
191 : ! arguments
192 : type(physics_state), intent(inout) :: phys_state(begchunk:endchunk)
193 : type(physics_tend), intent(inout) :: phys_tend(begchunk:endchunk)
194 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
195 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
196 :
197 : ! local variables
198 : integer :: tl_f, tl_fQdp
199 :
200 : integer :: c
201 : class(aerosol_state), pointer :: aero_state_obj
202 :
203 : !----------------------------------------------------------------------------
204 :
205 369408 : tl_f = TimeLevel%n0 ! timelevel which was adjusted by physics
206 369408 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_fQdp)
207 :
208 : !----------------------------------------------------------
209 : ! update physics state with aerosol constituents
210 : !----------------------------------------------------------
211 369408 : nullify(aero_state_obj)
212 :
213 369408 : if (aerosols_transported) then
214 0 : do c = begchunk,endchunk
215 0 : aero_state_obj => aerosol_state_object(c)
216 : ! get mass or number mixing ratios of aerosol constituents
217 0 : call aero_state_obj%get_transported(phys_state(c)%q)
218 : end do
219 : end if
220 :
221 369408 : call t_barrierf('sync_p_d_coupling', mpicom)
222 369408 : call t_startf('p_d_coupling')
223 : ! copy from phys structures -> dynamics structures
224 369408 : call p_d_coupling(phys_state, phys_tend, dyn_in, tl_f, tl_fQdp)
225 369408 : call t_stopf('p_d_coupling')
226 :
227 369408 : if (iam < par%nprocs) then
228 369408 : call tot_energy_dyn(dyn_in%elem,dyn_in%fvm, 1, nelemd, tl_f, tl_fQdp,'dED')
229 : end if
230 :
231 369408 : end subroutine stepon_run2
232 :
233 : !=========================================================================================
234 :
235 369408 : subroutine stepon_run3(dtime, cam_out, phys_state, dyn_in, dyn_out)
236 :
237 369408 : use camsrfexch, only: cam_out_t
238 : use dyn_comp, only: dyn_run
239 : use advect_tend, only: compute_adv_tends_xyz, compute_write_iop_fields
240 : use dyn_grid, only: TimeLevel
241 : use se_dyn_time_mod,only: TimeLevel_Qdp
242 : use control_mod, only: qsplit
243 : use constituents, only: pcnst, cnst_name
244 :
245 : ! arguments
246 : real(r8), intent(in) :: dtime ! Time-step
247 : type(cam_out_t), intent(inout) :: cam_out(:) ! Output from CAM to surface
248 : type(physics_state), intent(inout) :: phys_state(begchunk:endchunk)
249 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
250 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
251 :
252 : integer :: tl_f, tl_fQdp
253 : !--------------------------------------------------------------------------------------
254 :
255 369408 : if (single_column) then
256 : ! Update IOP properties e.g. omega, divT, divQ
257 0 : iop_update_phase1 = .false.
258 0 : if (doiopupdate) then
259 0 : if (masterproc) call readiopdata( hvcoord%hyam, hvcoord%hybm, hvcoord%hyai, hvcoord%hybi, hvcoord%ps0 )
260 0 : call iop_broadcast()
261 0 : call scm_setfield(dyn_out%elem,iop_update_phase1)
262 : endif
263 : endif
264 :
265 369408 : call t_startf('comp_adv_tends1')
266 369408 : tl_f = TimeLevel%n0
267 369408 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_fQdp)
268 369408 : call compute_adv_tends_xyz(dyn_in%elem,dyn_in%fvm,1,nelemd,tl_fQdp,tl_f)
269 369408 : if (write_camiop) call compute_write_iop_fields(dyn_in%elem,dyn_in%fvm,1,nelemd,tl_fQdp,tl_f)
270 369408 : call t_stopf('comp_adv_tends1')
271 :
272 369408 : call t_barrierf('sync_dyn_run', mpicom)
273 369408 : call t_startf('dyn_run')
274 369408 : call dyn_run(dyn_out)
275 369408 : call t_stopf('dyn_run')
276 :
277 369408 : call t_startf('comp_adv_tends2')
278 369408 : tl_f = TimeLevel%n0
279 369408 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_fQdp)
280 369408 : call compute_adv_tends_xyz(dyn_in%elem,dyn_in%fvm,1,nelemd,tl_fQdp,tl_f)
281 369408 : if (write_camiop) call compute_write_iop_fields(dyn_in%elem,dyn_in%fvm,1,nelemd,tl_fQdp,tl_f)
282 369408 : call t_stopf('comp_adv_tends2')
283 :
284 369408 : end subroutine stepon_run3
285 :
286 : !=========================================================================================
287 :
288 1536 : subroutine stepon_final(dyn_in, dyn_out)
289 :
290 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
291 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
292 :
293 369408 : end subroutine stepon_final
294 :
295 : !=========================================================================================
296 :
297 370944 : subroutine diag_dynvar_ic(elem, fvm)
298 : use constituents, only: cnst_type
299 : use dyn_grid, only: TimeLevel
300 :
301 : use se_dyn_time_mod, only: TimeLevel_Qdp ! dynamics typestep
302 : use control_mod, only: qsplit
303 : use hybrid_mod, only: config_thread_region, get_loop_ranges
304 : use hybrid_mod, only: hybrid_t
305 : use dimensions_mod, only: np, npsq, nc, nhc, fv_nphys, qsize, ntrac, nlev
306 : use dimensions_mod, only: cnst_name_gll
307 : use constituents, only: cnst_name
308 : use element_mod, only: element_t
309 : use fvm_control_volume_mod, only: fvm_struct
310 : use fvm_mapping, only: fvm2dyn
311 : use cam_thermo, only: get_sum_species, get_dp_ref, get_ps
312 : use air_composition, only: thermodynamic_active_species_idx
313 : use air_composition, only: thermodynamic_active_species_idx_dycore
314 : use hycoef, only: hyai, hybi, ps0
315 : ! arguments
316 : type(element_t) , intent(in) :: elem(1:nelemd)
317 : type(fvm_struct), intent(inout) :: fvm(:)
318 :
319 : ! local variables
320 : integer :: ie, i, j, k, m, m_cnst, nq
321 : integer :: tl_f, tl_qdp
322 : character(len=fieldname_len) :: tfname
323 :
324 : type(hybrid_t) :: hybrid
325 : integer :: nets, nete
326 370944 : real(r8), allocatable :: ftmp(:,:,:)
327 370944 : real(r8), allocatable :: fld_fvm(:,:,:,:,:), fld_gll(:,:,:,:,:)
328 370944 : real(r8), allocatable :: fld_2d(:,:)
329 : logical :: llimiter(1)
330 : real(r8) :: qtmp(np,np,nlev), dp_ref(np,np,nlev), ps_ref(np,np)
331 370944 : real(r8), allocatable :: factor_array(:,:,:)
332 : integer :: astat
333 : character(len=*), parameter :: prefix = 'diag_dynvar_ic: '
334 : !----------------------------------------------------------------------------
335 :
336 370944 : tl_f = timelevel%n0
337 370944 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_Qdp)
338 :
339 370944 : allocate(ftmp(npsq,nlev,2))
340 :
341 : ! Output tracer fields for analysis of advection schemes
342 2596608 : do m_cnst = 1, qsize
343 2225664 : tfname = trim(cnst_name_gll(m_cnst))//'_gll'
344 2596608 : if (hist_fld_active(tfname)) then
345 0 : do ie = 1, nelemd
346 0 : qtmp(:,:,:) = elem(ie)%state%Qdp(:,:,:,m_cnst,tl_qdp)/&
347 0 : elem(ie)%state%dp3d(:,:,:,tl_f)
348 0 : do j = 1, np
349 0 : do i = 1, np
350 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%Qdp(i,j,:,m_cnst,tl_qdp)/&
351 0 : elem(ie)%state%dp3d(i,j,:,tl_f)
352 : end do
353 : end do
354 0 : call outfld(tfname, ftmp(:,:,1), npsq, ie)
355 : end do
356 : end if
357 : end do
358 :
359 2596608 : do m_cnst = 1, qsize
360 2225664 : tfname = trim(cnst_name_gll(m_cnst))//'dp_gll'
361 2596608 : if (hist_fld_active(tfname)) then
362 0 : do ie = 1, nelemd
363 0 : do j = 1, np
364 0 : do i = 1, np
365 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%Qdp(i,j,:,m_cnst,tl_qdp)
366 : end do
367 : end do
368 0 : call outfld(tfname, ftmp(:,:,1), npsq, ie)
369 : end do
370 : end if
371 : end do
372 :
373 370944 : if (hist_fld_active('U_gll') .or. hist_fld_active('V_gll')) then
374 0 : do ie = 1, nelemd
375 0 : do j = 1, np
376 0 : do i = 1, np
377 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%v(i,j,1,:,tl_f)
378 0 : ftmp(i+(j-1)*np,:,2) = elem(ie)%state%v(i,j,2,:,tl_f)
379 : end do
380 : end do
381 0 : call outfld('U_gll', ftmp(:,:,1), npsq, ie)
382 0 : call outfld('V_gll', ftmp(:,:,2), npsq, ie)
383 : end do
384 : end if
385 :
386 370944 : if (hist_fld_active('T_gll')) then
387 0 : do ie = 1, nelemd
388 0 : do j = 1, np
389 0 : do i = 1, np
390 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%T(i,j,:,tl_f)
391 : end do
392 : end do
393 0 : call outfld('T_gll', ftmp(:,:,1), npsq, ie)
394 : end do
395 : end if
396 :
397 370944 : if (hist_fld_active('dp_ref_gll')) then
398 0 : do ie = 1, nelemd
399 0 : call get_dp_ref(hyai, hybi, ps0, elem(ie)%state%phis(:,:), dp_ref(:,:,:), ps_ref(:,:))
400 0 : do j = 1, np
401 0 : do i = 1, np
402 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%dp3d(i,j,:,tl_f)/dp_ref(i,j,:)
403 : end do
404 : end do
405 0 : call outfld('dp_ref_gll', ftmp(:,:,1), npsq, ie)
406 : end do
407 : end if
408 :
409 370944 : if (hist_fld_active('PSDRY_gll')) then
410 0 : do ie = 1, nelemd
411 0 : do j = 1, np
412 0 : do i = 1, np
413 0 : ftmp(i+(j-1)*np,1,1) = elem(ie)%state%psdry(i,j)
414 : end do
415 : end do
416 0 : call outfld('PSDRY_gll', ftmp(:,1,1), npsq, ie)
417 : end do
418 : end if
419 :
420 370944 : if (hist_fld_active('PS_gll')) then
421 0 : allocate(fld_2d(np,np))
422 0 : do ie = 1, nelemd
423 0 : call get_ps(elem(ie)%state%Qdp(:,:,:,:,tl_Qdp), thermodynamic_active_species_idx_dycore,&
424 0 : elem(ie)%state%dp3d(:,:,:,tl_f),fld_2d,hyai(1)*ps0)
425 0 : do j = 1, np
426 0 : do i = 1, np
427 0 : ftmp(i+(j-1)*np,1,1) = fld_2d(i,j)
428 : end do
429 : end do
430 0 : call outfld('PS_gll', ftmp(:,1,1), npsq, ie)
431 : end do
432 0 : deallocate(fld_2d)
433 : end if
434 :
435 370944 : if (hist_fld_active('PHIS_gll')) then
436 0 : do ie = 1, nelemd
437 0 : call outfld('PHIS_gll', RESHAPE(elem(ie)%state%phis, (/np*np/)), np*np, ie)
438 : end do
439 : end if
440 :
441 370944 : if (write_inithist()) then
442 0 : allocate(fld_2d(np,np))
443 0 : do ie = 1, nelemd
444 0 : call get_ps(elem(ie)%state%Qdp(:,:,:,:,tl_Qdp), thermodynamic_active_species_idx_dycore,&
445 0 : elem(ie)%state%dp3d(:,:,:,tl_f),fld_2d,hyai(1)*ps0)
446 0 : do j = 1, np
447 0 : do i = 1, np
448 0 : ftmp(i+(j-1)*np,1,1) = fld_2d(i,j)
449 : end do
450 : end do
451 0 : call outfld('PS&IC', ftmp(:,1,1), npsq, ie)
452 : end do
453 0 : deallocate(fld_2d)
454 : endif
455 :
456 370944 : deallocate(ftmp)
457 :
458 370944 : if (write_inithist()) then
459 :
460 0 : if (fv_nphys < 1) then
461 0 : allocate(factor_array(np,np,nlev),stat=astat)
462 0 : if (astat /= 0) call endrun(prefix//"Allocate factor_array failed")
463 : endif
464 :
465 0 : do ie = 1, nelemd
466 0 : call outfld('T&IC', RESHAPE(elem(ie)%state%T(:,:,:,tl_f), (/npsq,nlev/)), npsq, ie)
467 0 : call outfld('U&IC', RESHAPE(elem(ie)%state%v(:,:,1,:,tl_f), (/npsq,nlev/)), npsq, ie)
468 0 : call outfld('V&IC', RESHAPE(elem(ie)%state%v(:,:,2,:,tl_f), (/npsq,nlev/)), npsq, ie)
469 :
470 0 : if (fv_nphys < 1) then
471 0 : call get_sum_species(elem(ie)%state%Qdp(:,:,:,:,tl_qdp), &
472 0 : thermodynamic_active_species_idx_dycore, factor_array,dp_dry=elem(ie)%state%dp3d(:,:,:,tl_f))
473 0 : factor_array(:,:,:) = 1.0_r8/factor_array(:,:,:)
474 0 : do m_cnst = 1, qsize
475 0 : if (cnst_type(m_cnst) == 'wet') then
476 : call outfld(trim(cnst_name(m_cnst))//'&IC', &
477 0 : RESHAPE(factor_array(:,:,:)*elem(ie)%state%Qdp(:,:,:,m_cnst,tl_qdp)/&
478 0 : elem(ie)%state%dp3d(:,:,:,tl_f), (/npsq,nlev/)), npsq, ie)
479 : else
480 : call outfld(trim(cnst_name(m_cnst))//'&IC', &
481 0 : RESHAPE(elem(ie)%state%Qdp(:,:,:,m_cnst,tl_qdp)/&
482 0 : elem(ie)%state%dp3d(:,:,:,tl_f), (/npsq,nlev/)), npsq, ie)
483 : end if
484 : end do
485 : end if
486 : end do
487 :
488 0 : if (fv_nphys > 0) then
489 : !JMD $OMP PARALLEL NUM_THREADS(horz_num_threads), DEFAULT(SHARED), PRIVATE(hybrid,nets,nete,n)
490 : !JMD hybrid = config_thread_region(par,'horizontal')
491 0 : hybrid = config_thread_region(par,'serial')
492 0 : call get_loop_ranges(hybrid, ibeg=nets, iend=nete)
493 :
494 0 : allocate(fld_fvm(1-nhc:nc+nhc,1-nhc:nc+nhc,nlev,1,nets:nete),stat=astat)
495 0 : if (astat /= 0) call endrun(prefix//"Allocate fld_fvm failed")
496 0 : allocate(fld_gll(np,np,nlev,1,nets:nete),stat=astat)
497 0 : if (astat /= 0) call endrun(prefix//"Allocate fld_gll failed")
498 0 : allocate(factor_array(nc,nc,nlev),stat=astat)
499 0 : if (astat /= 0) call endrun(prefix//"Allocate factor_array failed")
500 :
501 0 : llimiter = .true.
502 :
503 0 : do m_cnst = 1, ntrac
504 0 : do ie = nets, nete
505 :
506 0 : call get_sum_species(fvm(ie)%c(1:nc,1:nc,:,:),thermodynamic_active_species_idx,factor_array)
507 0 : factor_array(:,:,:) = 1.0_r8/factor_array(:,:,:)
508 :
509 0 : if (cnst_type(m_cnst) == 'wet') then
510 0 : fld_fvm(1:nc,1:nc,:,1,ie) = fvm(ie)%c(1:nc,1:nc,:,m_cnst)*factor_array(:,:,:)
511 : else
512 0 : fld_fvm(1:nc,1:nc,:,1,ie) = fvm(ie)%c(1:nc,1:nc,:,m_cnst)
513 : end if
514 : end do
515 :
516 0 : call fvm2dyn(fld_fvm, fld_gll, hybrid, nets, nete, nlev, fvm(nets:nete), llimiter)
517 :
518 0 : do ie = nets, nete
519 0 : call outfld(trim(cnst_name(m_cnst))//'&IC', &
520 0 : RESHAPE(fld_gll(:,:,:,:,ie), (/npsq,nlev/)), npsq, ie)
521 : end do
522 : end do
523 :
524 0 : deallocate(fld_fvm)
525 0 : deallocate(fld_gll)
526 : end if
527 :
528 0 : deallocate(factor_array)
529 :
530 : end if ! if (write_inithist)
531 :
532 370944 : end subroutine diag_dynvar_ic
533 :
534 : !=========================================================================================
535 :
536 : end module stepon
|