Line data Source code
1 : !#define _DBG_ print *,"file: ",__FILE__," line: ",__LINE__," ithr: ",hybrid%ithr
2 : #define _DBG_
3 : module prim_driver_mod
4 : use shr_kind_mod, only: r8=>shr_kind_r8
5 : use cam_logfile, only: iulog
6 : use cam_abortutils, only: endrun
7 : use dimensions_mod, only: np, nlev, nelem, nelemd, GlobalUniqueCols, qsize, nc,nhc
8 : use hybrid_mod, only: hybrid_t, config_thread_region, PrintHybrid
9 : use derivative_mod, only: derivative_t
10 : use fvm_control_volume_mod, only: fvm_struct
11 :
12 : use element_mod, only: element_t, timelevels, allocate_element_desc
13 : use thread_mod , only: horz_num_threads, vert_num_threads, tracer_num_threads
14 : use thread_mod , only: omp_set_nested
15 : use perf_mod, only: t_startf, t_stopf
16 : use prim_init, only: gp, fvm_corners, fvm_points
17 :
18 : implicit none
19 : private
20 : public :: prim_init2, prim_run_subcycle, prim_finalize
21 : public :: prim_set_dry_mass
22 : contains
23 :
24 : !=============================================================================!
25 :
26 1536 : subroutine prim_init2(elem, fvm, hybrid, nets, nete, tl, hvcoord)
27 : use dimensions_mod, only: irecons_tracer, fvm_supercycling
28 : use dimensions_mod, only: fv_nphys, nc
29 : use parallel_mod, only: syncmp
30 : use se_dyn_time_mod, only: timelevel_t, tstep, phys_tscale, nsplit, TimeLevel_Qdp
31 : use se_dyn_time_mod, only: nsplit_baseline,rsplit_baseline
32 : use prim_state_mod, only: prim_printstate
33 : use control_mod, only: runtype, topology, rsplit, qsplit, rk_stage_user, &
34 : nu, nu_q, nu_div, hypervis_subcycle, hypervis_subcycle_q, &
35 : hypervis_subcycle_sponge, variable_nsplit
36 : use fvm_mod, only: fill_halo_fvm,ghostBufQnhc_h
37 : use thread_mod, only: omp_get_thread_num
38 : use global_norms_mod, only: print_cfl
39 : use hybvcoord_mod, only: hvcoord_t
40 : use prim_advection_mod, only: prim_advec_init2,deriv
41 : use prim_advance_mod, only: compute_omega
42 : use physconst, only: rga, cappa, cpair, tref, lapse_rate
43 : use cam_thermo, only: get_dp_ref
44 : use physconst, only: pstd
45 :
46 : type (element_t), intent(inout) :: elem(:)
47 : type (fvm_struct), intent(inout) :: fvm(:)
48 : type (hybrid_t), intent(in) :: hybrid
49 :
50 : type (TimeLevel_t), intent(inout) :: tl ! time level struct
51 : type (hvcoord_t), intent(inout) :: hvcoord ! hybrid vertical coordinate struct
52 :
53 : integer, intent(in) :: nets ! starting thread element number (private)
54 : integer, intent(in) :: nete ! ending thread element number (private)
55 :
56 :
57 : ! ==================================
58 : ! Local variables
59 : ! ==================================
60 :
61 : ! variables used to calculate CFL
62 : real (kind=r8) :: dtnu ! timestep*viscosity parameter
63 : real (kind=r8) :: dt_dyn_del2_sponge
64 : real (kind=r8) :: dt_tracer_vis ! viscosity timestep used in tracers
65 : real (kind=r8) :: dt_dyn_vis ! viscosity timestep
66 : real (kind=r8) :: dt_remap ! remapping timestep
67 :
68 : real (kind=r8) :: dp,dp0,T1,T0,pmid_ref(np,np)
69 3072 : real (kind=r8) :: ps_ref(np,np,nets:nete)
70 :
71 : integer :: i,j,k,ie,t,q
72 : integer :: n0,n0_qdp
73 :
74 :
75 12336 : do ie=nets,nete
76 43200000 : elem(ie)%derived%FM=0.0_r8
77 21103200 : elem(ie)%derived%FT=0.0_r8
78 126631536 : elem(ie)%derived%FQ=0.0_r8
79 : end do
80 :
81 : ! ==========================
82 : ! begin executable code
83 : ! ==========================
84 : !call prim_advance_init(hybrid%par,elem)
85 :
86 : ! compute most restrictive dt*nu for use by variable res viscosity:
87 : ! compute timestep seen by viscosity operator:
88 1536 : dt_dyn_vis = tstep
89 1536 : dt_dyn_del2_sponge = tstep
90 1536 : dt_tracer_vis=tstep*qsplit
91 1536 : dt_remap=dt_tracer_vis*rsplit
92 : ! compute most restrictive condition:
93 : ! note: dtnu ignores subcycling
94 1536 : dtnu=max(dt_dyn_vis*max(nu,nu_div), dt_tracer_vis*nu_q)
95 : ! compute actual viscosity timesteps with subcycling
96 1536 : dt_tracer_vis = dt_tracer_vis/hypervis_subcycle_q
97 1536 : dt_dyn_vis = dt_dyn_vis/hypervis_subcycle
98 1536 : dt_dyn_del2_sponge = dt_dyn_del2_sponge/hypervis_subcycle_sponge
99 1536 : if (variable_nsplit) then
100 0 : nsplit_baseline=nsplit
101 0 : rsplit_baseline=rsplit
102 : end if
103 : ! ==================================
104 : ! Initialize derivative structure
105 : ! ==================================
106 1536 : call Prim_Advec_Init2(fvm_corners, fvm_points)
107 1536 : if (fv_nphys>0.and.nc.ne.fv_nphys) then
108 : !
109 : ! need to fill halo for dp_coupling for fvm2phys mapping
110 : !
111 0 : call fill_halo_fvm(ghostBufQnhc_h,elem,fvm,hybrid,nets,nete,nhc,1,nlev,nlev)
112 : end if
113 : ! !$OMP BARRIER
114 : ! if (hybrid%ithr==0) then
115 : ! call syncmp(hybrid%par)
116 : ! end if
117 : ! !$OMP BARRIER
118 :
119 1536 : if (topology /= "cube") then
120 0 : call endrun('Error: only cube topology supported for primaitve equations')
121 : endif
122 :
123 : ! CAM has set tstep based on dtime before calling prim_init2(),
124 : ! so only now does HOMME learn the timstep. print them out:
125 : call print_cfl(elem,hybrid,nets,nete,dtnu,&
126 : !p top and p mid levels
127 : hvcoord%hyai(1)*hvcoord%ps0,hvcoord%hyam(:)*hvcoord%ps0+hvcoord%hybm(:)*pstd,&
128 : !dt_remap,dt_tracer_fvm,dt_tracer_se
129 : tstep*qsplit*rsplit,tstep*qsplit*fvm_supercycling,tstep*qsplit,&
130 : !dt_dyn,dt_dyn_visco,dt_tracer_visco, dt_phys
131 144384 : tstep,dt_dyn_vis,dt_dyn_del2_sponge,dt_tracer_vis,tstep*nsplit*qsplit*rsplit)
132 :
133 1536 : if (hybrid%masterthread) then
134 2 : if (phys_tscale/=0) then
135 0 : write(iulog,'(a,2f9.2)') "CAM physics timescale: ",phys_tscale
136 : endif
137 2 : write(iulog,'(a,2f9.2)') "CAM dtime (dt_phys): ",tstep*nsplit*qsplit*rsplit
138 :
139 2 : write(iulog,*) "CAM-SE uses dry-mass vertical coordinates"
140 : end if
141 :
142 1536 : n0=tl%n0
143 1536 : call TimeLevel_Qdp( tl, qsplit, n0_qdp)
144 1536 : call compute_omega(hybrid,n0,n0_qdp,elem,deriv,nets,nete,dt_remap,hvcoord)
145 : !
146 : ! pre-compute pressure-level thickness reference profile
147 : !
148 12336 : do ie=nets,nete
149 10800 : call get_dp_ref(hvcoord%hyai, hvcoord%hybi, hvcoord%ps0, elem(ie)%state%phis(:,:), &
150 23136 : elem(ie)%derived%dp_ref(:,:,:), ps_ref(:,:,ie))
151 : end do
152 : !
153 : ! pre-compute reference temperature profile (Simmons and Jiabin, 1991, QJRMS, Section 2a
154 : ! doi: https://doi.org/10.1002/qj.49711749703c)
155 : !
156 : ! Tref = T0+T1*Exner
157 : ! T1 = .0065*Tref*Cp/g ! = ~191
158 : ! T0 = Tref-T1 ! = ~97
159 : !
160 1536 : T1 = lapse_rate*Tref*cpair*rga
161 1536 : T0 = Tref-T1
162 12336 : do ie=nets,nete
163 1016736 : do k=1,nlev
164 21092400 : pmid_ref =hvcoord%hyam(k)*hvcoord%ps0 + hvcoord%hybm(k)*ps_ref(:,:,ie)
165 1004400 : dp0 = ( hvcoord%hyai(k+1) - hvcoord%hyai(k) )*hvcoord%ps0 + &
166 1004400 : ( hvcoord%hybi(k+1) - hvcoord%hybi(k) )*hvcoord%ps0
167 1015200 : if (hvcoord%hybm(k)>0) then
168 11340000 : elem(ie)%derived%T_ref(:,:,k) = T0+T1*(pmid_ref/hvcoord%ps0)**cappa
169 : !
170 : ! pel@ucar.edu: resolved noise issue over Antartica
171 : !
172 11340000 : elem(ie)%derived%dp_ref(:,:,k) = elem(ie)%derived%dp_ref(:,:,k)-dp0
173 : else
174 9752400 : elem(ie)%derived%T_ref(:,:,k) = 0.0_r8
175 : end if
176 : end do
177 : end do
178 :
179 1536 : if (hybrid%masterthread) write(iulog,*) "initial state:"
180 1536 : call prim_printstate(elem, tl, hybrid,nets,nete, fvm)
181 :
182 1536 : end subroutine prim_init2
183 :
184 : !=======================================================================================================!
185 :
186 :
187 29184 : subroutine prim_run_subcycle(elem, fvm, hybrid,nets,nete, dt, tl, hvcoord,nsubstep, omega_cn)
188 : !
189 : ! advance all variables (u,v,T,ps,Q,C) from time t to t + dt_q
190 : !
191 : ! for the RK schemes:
192 : ! input:
193 : ! tl%nm1 not used
194 : ! tl%n0 data at time t
195 : ! tl%np1 new values at t+dt_q
196 : !
197 : ! then we update timelevel pointers:
198 : ! tl%nm1 = tl%n0
199 : ! tl%n0 = tl%np1
200 : ! so that:
201 : ! tl%nm1 tracers: t dynamics: t+(qsplit-1)*dt
202 : ! tl%n0 time t + dt_q
203 : !
204 : ! for the implicit schemes:
205 : !
206 : ! input:
207 : ! tl%nm1 variables at t-1 level are stored fro BDF2 scheme
208 : ! tl%n0 data at time t
209 : ! tl%np1 new values at t+dt_q
210 : ! generally dt_q = t for BDF2, so its t+1
211 : !
212 : ! then we update timelevel pointers:
213 : ! tl%nm1 = tl%n0
214 : ! tl%n0 = tl%np1
215 : ! so that:
216 : ! tl%nm1 tracers: t dynamics: t+(qsplit-1)*dt
217 : ! tl%n0 time t + dt_q
218 : !
219 : !
220 1536 : use hybvcoord_mod, only : hvcoord_t
221 : use se_dyn_time_mod, only: TimeLevel_t, timelevel_update, timelevel_qdp, nsplit
222 : use control_mod, only: statefreq,qsplit, rsplit, variable_nsplit, dribble_in_rsplit_loop
223 : use prim_advance_mod, only: applycamforcing
224 : use prim_advance_mod, only: tot_energy_dyn,compute_omega
225 : use prim_state_mod, only: prim_printstate, adjust_nsplit
226 : use prim_advection_mod, only: vertical_remap, deriv
227 : use thread_mod, only: omp_get_thread_num
228 : use perf_mod , only: t_startf, t_stopf
229 : use fvm_mod , only: fill_halo_fvm, ghostBufQnhc_h
230 : use dimensions_mod, only: use_cslam,fv_nphys
231 : use fvm_mapping, only: cslam2gll
232 : type (element_t) , intent(inout) :: elem(:)
233 : type(fvm_struct), intent(inout) :: fvm(:)
234 : type (hybrid_t), intent(in) :: hybrid ! distributed parallel structure (shared)
235 : type (hvcoord_t), intent(in) :: hvcoord ! hybrid vertical coordinate struct
236 : integer, intent(in) :: nets ! starting thread element number (private)
237 : integer, intent(in) :: nete ! ending thread element number (private)
238 : real(kind=r8), intent(in) :: dt ! "timestep dependent" timestep
239 : type (TimeLevel_t), intent(inout):: tl
240 : integer, intent(in) :: nsubstep ! nsubstep = 1 .. nsplit
241 : real (kind=r8) , intent(inout):: omega_cn(2,nets:nete) !min and max of vertical Courant number
242 :
243 : real(kind=r8) :: dt_q, dt_remap, dt_phys
244 : integer :: ie, q,k,n0_qdp,np1_qdp,r, nstep_end,region_num_threads,i,j
245 : real (kind=r8) :: dp_np1(np,np)
246 58368 : real (kind=r8) :: dp_start(np,np,nlev+1,nets:nete),dp_end(np,np,nlev,nets:nete)
247 : logical :: compute_diagnostics
248 : ! ===================================
249 : ! Main timestepping loop
250 : ! ===================================
251 29184 : dt_q = dt*qsplit
252 29184 : nstep_end = tl%nstep + qsplit
253 29184 : dt_remap=dt_q*rsplit
254 29184 : nstep_end = tl%nstep + qsplit*rsplit ! nstep at end of this routine
255 29184 : dt_phys = nsplit*dt_remap
256 :
257 : ! compute diagnostics for STDOUT
258 29184 : compute_diagnostics=.false.
259 :
260 29184 : if (statefreq>0) then
261 0 : if (MODULO(nstep_end,statefreq)==0 .or. nstep_end==tl%nstep0) then
262 0 : compute_diagnostics=.true.
263 : endif
264 : end if
265 : !
266 : ! initialize variables for computing vertical Courant number
267 : !
268 29184 : if (variable_nsplit.or.compute_diagnostics) then
269 0 : if (nsubstep==1) then
270 0 : do ie=nets,nete
271 0 : omega_cn(1,ie) = 0.0_r8
272 0 : omega_cn(2,ie) = 0.0_r8
273 : end do
274 : end if
275 0 : do ie=nets,nete
276 0 : dp_start(:,:,1:nlev,ie) = elem(ie)%state%dp3d(:,:,:,tl%n0)
277 0 : dp_start(:,:,nlev+1,ie) = elem(ie)%state%dp3d(:,:,nlev,tl%n0)
278 : end do
279 : endif
280 :
281 :
282 29184 : call TimeLevel_Qdp( tl, qsplit, n0_qdp)
283 :
284 29184 : if (dribble_in_rsplit_loop==0) then
285 29184 : call tot_energy_dyn(elem,fvm,nets,nete,tl%n0,n0_qdp,'dAF')
286 29184 : call ApplyCAMForcing(elem,fvm,tl%n0,n0_qdp,dt_remap,dt_phys,nets,nete,nsubstep)
287 29184 : call tot_energy_dyn(elem,fvm,nets,nete,tl%n0,n0_qdp,'dBD')
288 : end if
289 116736 : do r=1,rsplit
290 87552 : if (r.ne.1) call TimeLevel_update(tl,"leapfrog")
291 : !
292 : ! if nsplit==1 and physics time-step is long then there will be noise in the
293 : ! pressure field; hence "dripple" in tendencies
294 : !
295 87552 : if (dribble_in_rsplit_loop==1) then
296 0 : call tot_energy_dyn(elem,fvm,nets,nete,tl%n0,n0_qdp,'dAF')
297 0 : call ApplyCAMForcing(elem,fvm,tl%n0,n0_qdp,dt,dt_phys,nets,nete,MAX(nsubstep,r))
298 0 : call tot_energy_dyn(elem,fvm,nets,nete,tl%n0,n0_qdp,'dBD')
299 : end if
300 : !
301 : ! right after physics overwrite Qdp with CSLAM values
302 : !
303 87552 : if (use_cslam.and.nsubstep==1.and.r==1) then
304 14592 : call tot_energy_dyn(elem,fvm,nets,nete,tl%n0,n0_qdp,'dAF')
305 14592 : call cslam2gll(elem, fvm, hybrid,nets,nete, tl%n0, n0_qdp)
306 14592 : call tot_energy_dyn(elem,fvm,nets,nete,tl%n0,n0_qdp,'dBD')
307 : end if
308 87552 : call tot_energy_dyn(elem,fvm,nets,nete,tl%n0,n0_qdp,'dBL')
309 87552 : call prim_step(elem, fvm, hybrid,nets,nete, dt, tl, hvcoord,r,nsubstep==nsplit,dt_remap)
310 116736 : call tot_energy_dyn(elem,fvm,nets,nete,tl%np1,n0_qdp,'dAL')
311 : enddo
312 :
313 :
314 : ! defer final timelevel update until after remap and diagnostics
315 29184 : call TimeLevel_Qdp( tl, qsplit, n0_qdp, np1_qdp)
316 :
317 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
318 : !
319 : ! apply vertical remap
320 : ! always for tracers
321 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
322 :
323 29184 : call tot_energy_dyn(elem,fvm,nets,nete,tl%np1,np1_qdp,'dAD')
324 :
325 29184 : if (variable_nsplit.or.compute_diagnostics) then
326 : !
327 : ! initialize variables for computing vertical Courant number
328 : !
329 0 : do ie=nets,nete
330 0 : dp_end(:,:,:,ie) = elem(ie)%state%dp3d(:,:,:,tl%np1)
331 : end do
332 : end if
333 29184 : call t_startf('vertical_remap')
334 29184 : call vertical_remap(hybrid,elem,fvm,hvcoord,tl%np1,np1_qdp,nets,nete)
335 29184 : call t_stopf('vertical_remap')
336 :
337 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
338 : ! time step is complete.
339 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
340 29184 : call tot_energy_dyn(elem,fvm,nets,nete,tl%np1,np1_qdp,'dAR')
341 :
342 29184 : if (nsubstep==nsplit) then
343 14592 : call compute_omega(hybrid,tl%np1,np1_qdp,elem,deriv,nets,nete,dt_remap,hvcoord)
344 : end if
345 :
346 : ! now we have:
347 : ! u(nm1) dynamics at t+dt_remap - 2*dt
348 : ! u(n0) dynamics at t+dt_remap - dt
349 : ! u(np1) dynamics at t+dt_remap
350 : !
351 : ! Q(1) Q at t+dt_remap
352 :
353 :
354 :
355 : ! =================================
356 : ! update dynamics time level pointers
357 : ! =================================
358 29184 : call TimeLevel_update(tl,"leapfrog")
359 : ! note: time level update for fvm tracers takes place in fvm_mod
360 :
361 : ! now we have:
362 : ! u(nm1) dynamics at t+dt_remap - dt (Robert-filtered)
363 : ! u(n0) dynamics at t+dt_remap
364 : ! u(np1) undefined
365 :
366 :
367 : !
368 : ! Compute vertical Courant numbers
369 : !
370 29184 : if (variable_nsplit.or.compute_diagnostics) then
371 0 : do ie=nets,nete
372 0 : do k=1,nlev
373 0 : do j=1,np
374 0 : do i=1,np
375 0 : if (dp_end(i,j,k,ie)<dp_start(i,j,k,ie)) then
376 0 : omega_cn(1,ie) = MIN((dp_end(i,j,k,ie)-dp_start(i,j,k,ie))/dp_start(i,j,k,ie),omega_cn(1,ie))
377 0 : omega_cn(2,ie) = MAX((dp_end(i,j,k,ie)-dp_start(i,j,k,ie))/dp_start(i,j,k,ie),omega_cn(2,ie))
378 : else
379 0 : omega_cn(1,ie) = MIN((dp_end(i,j,k,ie)-dp_start(i,j,k,ie))/dp_start(i,j,k+1,ie),omega_cn(1,ie))
380 0 : omega_cn(2,ie) = MAX((dp_end(i,j,k,ie)-dp_start(i,j,k,ie))/dp_start(i,j,k+1,ie),omega_cn(2,ie))
381 : end if
382 : end do
383 : end do
384 : end do
385 : end do
386 0 : if (nsubstep==nsplit.and.variable_nsplit) then
387 0 : call t_startf('adjust_nsplit')
388 0 : call adjust_nsplit(elem, tl, hybrid,nets,nete, fvm, omega_cn)
389 0 : call t_stopf('adjust_nsplit')
390 : end if
391 : end if
392 :
393 : ! ============================================================
394 : ! Print some diagnostic information
395 : ! ============================================================
396 0 : if (compute_diagnostics) then
397 0 : call prim_printstate(elem, tl, hybrid,nets,nete, fvm, omega_cn)
398 : end if
399 :
400 29184 : if (use_cslam.and.nsubstep==nsplit.and.nc.ne.fv_nphys) then
401 : !
402 : ! fill the fvm halo for mapping in d_p_coupling if
403 : ! physics grid resolution is different than fvm resolution
404 : !
405 0 : call fill_halo_fvm(ghostBufQnhc_h, elem,fvm,hybrid,nets,nete,nhc,1,nlev,nlev)
406 : end if
407 :
408 29184 : end subroutine prim_run_subcycle
409 :
410 :
411 87552 : subroutine prim_step(elem, fvm, hybrid,nets,nete, dt, tl, hvcoord, rstep, last_step,dt_remap)
412 : !
413 : ! Take qsplit dynamics steps and one tracer step
414 : ! for vertically lagrangian option, this subroutine does only the horizontal step
415 : !
416 : ! input:
417 : ! tl%nm1 not used
418 : ! tl%n0 data at time t
419 : ! tl%np1 new values at t+dt_q
420 : !
421 : ! then we update timelevel pointers:
422 : ! tl%nm1 = tl%n0
423 : ! tl%n0 = tl%np1
424 : ! so that:
425 : ! tl%nm1 tracers: t dynamics: t+(qsplit-1)*dt
426 : ! tl%n0 time t + dt_q
427 : !
428 29184 : use hybvcoord_mod, only: hvcoord_t
429 : use se_dyn_time_mod, only: TimeLevel_t, timelevel_update
430 : use control_mod, only: statefreq, qsplit, nu_p
431 : use thread_mod, only: omp_get_thread_num
432 : use prim_advance_mod, only: prim_advance_exp
433 : use prim_advection_mod, only: prim_advec_tracers_remap, prim_advec_tracers_fvm, deriv
434 : use derivative_mod, only: subcell_integration
435 : use hybrid_mod, only: set_region_num_threads, config_thread_region, get_loop_ranges
436 : use dimensions_mod, only: use_cslam,fvm_supercycling,fvm_supercycling_jet
437 : use dimensions_mod, only: kmin_jet, kmax_jet
438 : use fvm_mod, only: ghostBufQnhc_vh,ghostBufQ1_vh, ghostBufFlux_vh
439 : use fvm_mod, only: ghostBufQ1_h,ghostBufQnhcJet_h, ghostBufFluxJet_h
440 : use se_dyn_time_mod, only: timelevel_qdp
441 : use fvm_mapping, only: cslam2gll
442 : #ifdef waccm_debug
443 : use cam_history, only: outfld
444 : #endif
445 :
446 :
447 : type (element_t) , intent(inout) :: elem(:)
448 : type(fvm_struct), intent(inout) :: fvm(:)
449 : type (hybrid_t), intent(in) :: hybrid ! distributed parallel structure (shared)
450 : type (hvcoord_t), intent(in) :: hvcoord ! hybrid vertical coordinate struct
451 : integer, intent(in) :: nets ! starting thread element number (private)
452 : integer, intent(in) :: nete ! ending thread element number (private)
453 : real(kind=r8), intent(in) :: dt ! "timestep dependent" timestep
454 : type (TimeLevel_t), intent(inout) :: tl
455 : integer, intent(in) :: rstep ! vertical remap subcycling step
456 : logical, intent(in) :: last_step! last step before d_p_coupling
457 : real(kind=r8), intent(in) :: dt_remap
458 :
459 : type (hybrid_t):: hybridnew,hybridnew2
460 : real(kind=r8) :: st, st1, dp, dt_q
461 : integer :: ie,t,q,k,i,j,n, n_Q
462 : integer :: ithr
463 : integer :: region_num_threads
464 : integer :: kbeg,kend
465 : integer :: n0_qdp, np1_qdp
466 :
467 : real (kind=r8) :: tempdp3d(np,np), x
468 : real (kind=r8) :: tempmass(nc,nc)
469 : real (kind=r8) :: tempflux(nc,nc,4)
470 :
471 : real (kind=r8) :: dp_np1(np,np)
472 :
473 :
474 87552 : dt_q = dt*qsplit
475 : ! ===============
476 : ! initialize mean flux accumulation variables and save some variables at n0
477 : ! for use by advection
478 : ! ===============
479 703152 : do ie=nets,nete
480 2462400000 : elem(ie)%derived%vn0=0 ! mean horizontal mass flux
481 1202882400 : elem(ie)%derived%omega=0
482 615600 : if (nu_p>0) then
483 1202882400 : elem(ie)%derived%dpdiss_ave=0
484 1202882400 : elem(ie)%derived%dpdiss_biharmonic=0
485 : endif
486 :
487 : ! dp at time t: use floating lagrangian levels:
488 1202969952 : elem(ie)%derived%dp(:,:,:)=elem(ie)%state%dp3d(:,:,:,tl%n0)
489 : enddo
490 :
491 : ! ===============
492 : ! Dynamical Step
493 : ! ===============
494 87552 : n_Q = tl%n0 ! n_Q = timelevel of FV tracers at time t. need to save this
495 : ! SE tracers only carry 2 timelevels
496 :
497 87552 : call t_startf('prim_advance_exp')
498 : ! ithr = 0 ! omp_get_thread_num()
499 : ! vybrid = hybrid_create(hybrid%par,ithr)
500 :
501 : call prim_advance_exp(elem, fvm, deriv, hvcoord, &
502 87552 : hybrid, dt, tl, nets, nete)
503 :
504 87552 : call t_stopf('prim_advance_exp')
505 :
506 87552 : do n=2,qsplit
507 0 : call TimeLevel_update(tl,"leapfrog")
508 :
509 0 : call t_startf('prim_advance_exp')
510 :
511 : call prim_advance_exp(elem, fvm, deriv, hvcoord, &
512 0 : hybrid, dt, tl, nets, nete)
513 :
514 87552 : call t_stopf('prim_advance_exp')
515 :
516 : ! defer final timelevel update until after Q update.
517 : enddo
518 : #ifdef HOMME_TEST_SUB_ELEMENT_MASS_FLUX
519 : if (use_cslam.and.rstep==1) then
520 : do ie=nets,nete
521 : do k=1,nlev
522 : tempdp3d = elem(ie)%state%dp3d(:,:,k,tl%np1) - &
523 : elem(ie)%derived%dp(:,:,k)
524 : call subcell_integration(tempdp3d, np, nc, elem(ie)%metdet,tempmass)
525 : tempflux = dt_q*elem(ie)%sub_elem_mass_flux(:,:,:,k)
526 : do i=1,nc
527 : do j=1,nc
528 : x = SUM(tempflux(i,j,:))
529 : if (ABS(tempmass(i,j)).lt.1e-11_r8 .and. 1e-11_r8.lt.ABS(x)) then
530 : write(iulog,*) __FILE__,__LINE__,"**CSLAM mass-flux ERROR***",ie,k,i,j,tempmass(i,j),x
531 : call endrun('**CSLAM mass-flux ERROR***')
532 : elseif (1e-5_r8.lt.ABS((tempmass(i,j)-x)/tempmass(i,j))) then
533 : write(iulog,*) __FILE__,__LINE__,"**CSLAM mass-flux ERROR**",ie,k,i,j,tempmass(i,j),x,&
534 : ABS((tempmass(i,j)-x)/tempmass(i,j))
535 : call endrun('**CSLAM mass-flux ERROR**')
536 : endif
537 : end do
538 : end do
539 : end do
540 : end do
541 : end if
542 : #endif
543 : ! current dynamics state variables:
544 : ! derived%dp = dp at start of timestep
545 : ! derived%vn0 = mean horiz. flux: U*dp
546 : ! rsplit>0
547 : ! state%v(:,:,:,np1) = velocity on lagrangian levels
548 : ! state%dp3d(:,:,:,np1) = dp3d
549 : !
550 :
551 :
552 : ! ===============
553 : ! Tracer Advection.
554 : ! in addition, this routine will apply the DSS to:
555 : ! derived%omega =
556 : ! Tracers are always vertically lagrangian.
557 : ! ===============
558 : ! Advect tracers if their count is > 0.
559 : ! special case in CAM: if CSLAM tracers are turned on , qsize=1 but this tracer should
560 : ! not be advected. This will be cleaned up when the physgrid is merged into CAM trunk
561 : ! Currently advecting all species
562 87552 : if (.not.use_cslam) then
563 0 : call t_startf('prim_advec_tracers_remap')
564 0 : region_num_threads=tracer_num_threads
565 0 : call omp_set_nested(.true.)
566 : !$OMP PARALLEL NUM_THREADS(region_num_threads), DEFAULT(SHARED), PRIVATE(hybridnew)
567 0 : hybridnew = config_thread_region(hybrid,'tracer')
568 0 : call Prim_Advec_Tracers_remap(elem, deriv,hvcoord,hybridnew,dt_q,tl,nets,nete)
569 : !$OMP END PARALLEL
570 0 : call omp_set_nested(.false.)
571 0 : call t_stopf('prim_advec_tracers_remap')
572 : else
573 : !
574 : ! only run fvm transport every fvm_supercycling rstep
575 : !
576 : ! FVM transport
577 : !
578 87552 : if ((mod(rstep,fvm_supercycling) == 0).and.(mod(rstep,fvm_supercycling_jet) == 0)) then
579 :
580 : ! call omp_set_nested(.true.)
581 : ! !$OMP PARALLEL NUM_THREADS(vert_num_threads), DEFAULT(SHARED), PRIVATE(hybridnew2,kbeg,kend)
582 : ! hybridnew2 = config_thread_region(hybrid,'vertical')
583 : ! call get_loop_ranges(hybridnew2,kbeg=kbeg,kend=kend)
584 : call Prim_Advec_Tracers_fvm(elem,fvm,hvcoord,hybrid,&
585 29184 : dt_q,tl,nets,nete,ghostBufQnhc_vh,ghostBufQ1_vh, ghostBufFlux_vh,1,nlev)
586 : ! !$OMP END PARALLEL
587 : ! call omp_set_nested(.false.)
588 : !
589 : ! to avoid accumulation of truncation error overwrite CSLAM surface pressure with SE
590 : ! surface pressure
591 : !
592 234384 : do ie=nets,nete
593 : !
594 : ! overwrite PSDRY on CSLAM grid with SE PSDRY integrated over CSLAM control volume
595 : !
596 : ! call subcell_integration(elem(ie)%state%psdry(:,:), np, nc, elem(ie)%metdet,fvm(ie)%psc)
597 : ! fvm(ie)%psc = fvm(ie)%psc*fvm(ie)%inv_se_area_sphere
598 : !
599 : ! Update CSLAM surface pressure
600 : !
601 849984 : do j=1,nc
602 2667600 : do i=1,nc
603 174214800 : fvm(ie)%psc(i,j) = sum(fvm(ie)%dp_fvm(i,j,:)) + hvcoord%hyai(1)*hvcoord%ps0
604 : end do
605 : end do
606 : end do
607 29184 : call TimeLevel_Qdp( tl, qsplit, n0_qdp, np1_qdp)
608 29184 : if (.not.last_step) call cslam2gll(elem, fvm, hybrid,nets,nete, tl%np1, np1_qdp)
609 58368 : else if ((mod(rstep,fvm_supercycling_jet) == 0)) then
610 : !
611 : ! shorter fvm time-step in jet region
612 : !
613 : call Prim_Advec_Tracers_fvm(elem,fvm,hvcoord,hybrid,&
614 0 : dt_q,tl,nets,nete,ghostBufQnhcJet_h,ghostBufQ1_h, ghostBufFluxJet_h,kmin_jet,kmax_jet)
615 : end if
616 :
617 : #ifdef waccm_debug
618 : do ie=nets,nete
619 : call outfld('CSLAM_gamma', RESHAPE(fvm(ie)%CSLAM_gamma(:,:,:,1), &
620 : (/nc*nc,nlev/)), nc*nc, ie)
621 : end do
622 : #endif
623 : endif
624 :
625 87552 : end subroutine prim_step
626 :
627 :
628 : !=======================================================================================================!
629 :
630 :
631 0 : subroutine prim_finalize(hybrid)
632 : type (hybrid_t), intent(in) :: hybrid ! distributed parallel structure (shared)
633 :
634 : ! ==========================
635 : ! end of the hybrid program
636 : ! ==========================
637 87552 : end subroutine prim_finalize
638 :
639 : !=========================================================================================
640 :
641 768 : subroutine prim_set_dry_mass(elem, hvcoord,initial_global_ave_dry_ps,q)
642 : use element_mod, only: element_t
643 : use hybvcoord_mod , only: hvcoord_t
644 : use dimensions_mod, only: nelemd, nlev, np
645 : use constituents, only: cnst_type, qmin, pcnst
646 : use cam_logfile, only: iulog
647 : use spmd_utils, only: masterproc
648 :
649 : type (element_t) , intent(inout):: elem(:)
650 : type (hvcoord_t) , intent(in) :: hvcoord
651 : real (kind=r8), intent(in) :: initial_global_ave_dry_ps
652 : real (kind=r8), intent(inout):: q(np,np,nlev,nelemd,pcnst)
653 :
654 : ! local
655 : real (kind=r8) :: global_ave_ps_inic,dp_tmp, factor(np,np,nlev)
656 : integer :: ie, i, j ,k, m_cnst
657 :
658 768 : if (initial_global_ave_dry_ps == 0) return;
659 :
660 768 : call get_global_ave_surface_pressure(elem, global_ave_ps_inic)
661 :
662 6168 : do ie=1,nelemd
663 113400 : elem(ie)%state%psdry(:,:)=elem(ie)%state%psdry(:,:)*(initial_global_ave_dry_ps/global_ave_ps_inic)
664 507600 : do k=1,nlev
665 2516400 : do j = 1,np
666 10546200 : do i = 1,np
667 16070400 : dp_tmp = ((hvcoord%hyai(k+1) - hvcoord%hyai(k))*hvcoord%ps0)+&
668 24105600 : ((hvcoord%hybi(k+1) - hvcoord%hybi(k))*elem(ie)%state%psdry(i,j))
669 8035200 : factor(i,j,k) = elem(ie)%state%dp3d(i,j,k,1)/dp_tmp
670 34149600 : elem(ie)%state%dp3d(i,j,k,:) = dp_tmp
671 : end do
672 : end do
673 : end do
674 : !
675 : ! conserve initial condition mass of 'wet' tracers (following dryairm.F90 for FV dycore)
676 : ! and conserve mixing ratio (not mass) of 'dry' tracers
677 : !
678 227568 : do m_cnst=1,pcnst
679 226800 : if (cnst_type(m_cnst).ne.'dry') then
680 5583600 : do k=1,nlev
681 27680400 : do j = 1,np
682 116008200 : do i = 1,np
683 88387200 : q(i,j,k,ie,m_cnst) = q(i,j,k,ie,m_cnst)*factor(i,j,k)
684 110484000 : q(i,j,k,ie,m_cnst) = max(qmin(m_cnst),q(i,j,k,ie,m_cnst))
685 : end do
686 : end do
687 : end do
688 : end if
689 : end do
690 : end do
691 768 : if (masterproc) then
692 1 : write (iulog,*) "------ info from prim_set_dry_mass -----------------------------------------------------------"
693 1 : write (iulog,*) "Scaling dry surface pressure to global average of = ",&
694 2 : initial_global_ave_dry_ps/100.0_r8,"hPa"
695 1 : write (iulog,*) "Average dry surface pressure in initial condition = ",global_ave_ps_inic/100.0_r8,"hPa"
696 1 : write (iulog,*) "Average dry surface pressure change = ",&
697 2 : initial_global_ave_dry_ps-global_ave_ps_inic,"Pa"
698 1 : write (iulog,*) "Mixing ratios that are wet have been scaled so that total mass of tracer is conserved"
699 1 : write (iulog,*) "Mixing ratios that are dry have not been changed (mass not conserved in scaling process)"
700 1 : write (iulog,*) "------ end info from prim_set_dry_mass -------------------------------------------------------"
701 : endif
702 : end subroutine prim_set_dry_mass
703 :
704 768 : subroutine get_global_ave_surface_pressure(elem, global_ave_ps_inic)
705 : use element_mod , only : element_t
706 : use dimensions_mod , only : np
707 : use global_norms_mod , only : global_integral
708 : use hybrid_mod , only : config_thread_region, get_loop_ranges, hybrid_t
709 : use parallel_mod , only : par
710 :
711 : type (element_t) , intent(in) :: elem(:)
712 : real (kind=r8), intent(out) :: global_ave_ps_inic
713 :
714 : ! local
715 768 : real (kind=r8), allocatable :: tmp(:,:,:)
716 : type (hybrid_t) :: hybrid
717 : integer :: ie, nets, nete
718 :
719 : !JMD $OMP PARALLEL NUM_THREADS(horz_num_threads), DEFAULT(SHARED), PRIVATE(hybrid,nets,nete,n)
720 : !JMD hybrid = config_thread_region(par,'horizontal')
721 768 : hybrid = config_thread_region(par,'serial')
722 768 : call get_loop_ranges(hybrid,ibeg=nets,iend=nete)
723 2304 : allocate(tmp(np,np,nets:nete))
724 :
725 6168 : do ie=nets,nete
726 114168 : tmp(:,:,ie)=elem(ie)%state%psdry(:,:)
727 : enddo
728 768 : global_ave_ps_inic = global_integral(elem, tmp(:,:,nets:nete),hybrid,np,nets,nete)
729 768 : deallocate(tmp)
730 768 : end subroutine get_global_ave_surface_pressure
731 :
732 : end module prim_driver_mod
|