Line data Source code
1 : module prim_state_mod
2 : use shr_kind_mod, only: r8=>shr_kind_r8
3 : use cam_logfile, only: iulog
4 : use dimensions_mod, only: nlev, np, nc, qsize_d, ntrac_d
5 : use parallel_mod, only: ordered
6 : use hybrid_mod, only: hybrid_t
7 : use se_dyn_time_mod, only: timelevel_t, TimeLevel_Qdp, time_at
8 : use control_mod, only: qsplit, statediag_numtrac
9 : use global_norms_mod, only: global_integrals_general
10 : use element_mod, only: element_t
11 : use reduction_mod, only: parallelmax,parallelmin
12 : use fvm_control_volume_mod, only: fvm_struct
13 :
14 : implicit none
15 : private
16 :
17 : public :: prim_printstate, adjust_nsplit
18 :
19 : CONTAINS
20 :
21 1536 : subroutine prim_printstate(elem, tl,hybrid,nets,nete, fvm, omega_cn)
22 : use dimensions_mod, only: use_cslam
23 : use constituents, only: cnst_name
24 : use air_composition, only: thermodynamic_active_species_idx_dycore, dry_air_species_num
25 : use air_composition, only: thermodynamic_active_species_num,thermodynamic_active_species_idx
26 : use cam_control_mod, only: initial_run
27 : use se_dyn_time_mod, only: tstep
28 : use control_mod, only: rsplit, qsplit
29 : use perf_mod, only: t_startf, t_stopf
30 : type (element_t), intent(inout) :: elem(:)
31 : type (TimeLevel_t), target, intent(in) :: tl
32 : type (hybrid_t), intent(in) :: hybrid
33 : integer, intent(in) :: nets,nete
34 : type(fvm_struct), intent(inout) :: fvm(:)
35 : real (kind=r8), optional, intent(in) :: omega_cn(2,nets:nete)
36 : ! Local variables...
37 : integer :: k,ie,m_cnst
38 : integer, parameter :: type=ORDERED
39 :
40 : integer, parameter :: vmax=11+2*MAX(qsize_d,ntrac_d)
41 :
42 : character(len=10) :: varname(vmax)
43 :
44 3072 : real (kind=r8), dimension(nets:nete,vmax) :: min_local,max_local
45 : real (kind=r8), dimension(vmax) :: min_p,max_p,mass,mass_chg
46 3072 : real (kind=r8), dimension(np,np,nets:nete):: moist_ps
47 3072 : real (kind=r8), dimension(nc,nc,nets:nete):: moist_ps_fvm
48 :
49 3072 : real (kind=r8) :: tmp_gll(np,np,vmax,nets:nete),tmp_mass(vmax)!
50 3072 : real (kind=r8) :: tmp_fvm(nc,nc,vmax,nets:nete)
51 : real (kind=r8) :: tmp_q(np,np,nlev)
52 : integer :: n0, n0_qdp, q, nm, nm2
53 3072 : real(kind=r8) :: da_gll(np,np,nets:nete),da_fvm(nc,nc,nets:nete)
54 :
55 : !dynamics variables in n0 are at time = 'time': time=tl%nstep*tstep
56 1536 : if (hybrid%masterthread) then
57 2 : write(iulog,*) "nstep=",tl%nstep," time=",Time_at(tl%nstep)/(24*3600)," [day]"
58 : end if
59 : ! dynamics timelevels
60 1536 : n0=tl%n0
61 1536 : call TimeLevel_Qdp( tl, qsplit, n0_qdp)
62 : ! moist surface pressure
63 1536 : if (use_cslam) then
64 12336 : do ie=nets,nete
65 9180000 : moist_ps_fvm(:,:,ie)=SUM(fvm(ie)%dp_fvm(1:nc,1:nc,:),DIM=3)
66 77136 : do q=dry_air_species_num+1,thermodynamic_active_species_num
67 64800 : m_cnst = thermodynamic_active_species_idx(q)
68 6102000 : do k=1,nlev
69 : moist_ps_fvm(:,:,ie) = moist_ps_fvm(:,:,ie)+&
70 78408000 : fvm(ie)%dp_fvm(1:nc,1:nc,k)*fvm(ie)%c(1:nc,1:nc,k,m_cnst)
71 : end do
72 : end do
73 : enddo
74 : end if
75 12336 : do ie=nets,nete
76 226800 : moist_ps(:,:,ie)=elem(ie)%state%psdry(:,:)
77 77136 : do q=dry_air_species_num+1,thermodynamic_active_species_num
78 64800 : m_cnst = thermodynamic_active_species_idx_dycore(q)
79 6102000 : do k=1,nlev
80 : moist_ps(:,:,ie) = moist_ps(:,:,ie)+&
81 126619200 : elem(ie)%state%Qdp(:,:,k,m_cnst,n0_qdp)
82 : end do
83 : end do
84 : enddo
85 : ! weights/areas for global integrals
86 12336 : do ie=nets,nete
87 228336 : da_gll(:,:,ie) = elem(ie)%mp(:,:)*elem(ie)%metdet(:,:)
88 : enddo
89 1536 : if (use_cslam) then
90 12336 : do ie=nets,nete
91 141936 : da_fvm(:,:,ie) = fvm(ie)%area_sphere(:,:)
92 : enddo
93 : end if
94 : !
95 : !*********************************************
96 : !
97 : ! min/max of u,v,T,PS,OMEGA
98 : !
99 : !*********************************************
100 : !
101 1536 : varname(1) = 'U '
102 1536 : varname(2) = 'V '
103 1536 : varname(3) = 'T '
104 1536 : varname(4) = 'OMEGA '
105 1536 : varname(5) = 'OMEGA CN '
106 1536 : if (use_cslam) then
107 1536 : varname(6) = 'PSDRY(fvm)'
108 1536 : varname(7) = 'PS(fvm) '
109 1536 : varname(8) = 'PSDRY(gll)'
110 1536 : varname(9) = 'PS(gll) '
111 1536 : nm = 9 !number of vars before tracers
112 1536 : nm2 = nm+statediag_numtrac!number of vars after tracers
113 : else
114 0 : varname(6) = 'PSDRY '
115 0 : varname(7) = 'PS '
116 0 : nm = 7 !number of vars before tracers
117 0 : nm2 = nm+statediag_numtrac!number of vars after tracers
118 : end if
119 :
120 12336 : do ie=nets,nete
121 21103200 : min_local(ie,1) = MINVAL(elem(ie)%state%v(:,:,1,:,n0))
122 21103200 : max_local(ie,1) = MAXVAL(elem(ie)%state%v(:,:,1,:,n0))
123 21103200 : min_local(ie,2) = MINVAL(elem(ie)%state%v(:,:,2,:,n0))
124 21103200 : max_local(ie,2) = MAXVAL(elem(ie)%state%v(:,:,2,:,n0))
125 21103200 : min_local(ie,3) = MINVAL(elem(ie)%state%T(:,:,:,n0))
126 21103200 : max_local(ie,3) = MAXVAL(elem(ie)%state%T(:,:,:,n0))
127 21103200 : min_local(ie,4) = MINVAL(elem(ie)%derived%Omega(:,:,:))
128 21103200 : max_local(ie,4) = MAXVAL(elem(ie)%derived%Omega(:,:,:))
129 10800 : if (present(omega_cn)) then
130 0 : min_local(ie,5) = omega_cn(1,ie)
131 0 : max_local(ie,5) = omega_cn(2,ie)
132 : else
133 10800 : min_local(ie,5) = 0.0_r8
134 10800 : max_local(ie,5) = 0.0_r8
135 : end if
136 10800 : if (use_cslam) then
137 9180000 : min_local(ie,6) = MINVAL(SUM(fvm(ie)%dp_fvm(1:nc,1:nc,:),DIM=3))
138 9180000 : max_local(ie,6) = MAXVAL(SUM(fvm(ie)%dp_fvm(1:nc,1:nc,:),DIM=3))
139 140400 : min_local(ie,7) = MINVAL(moist_ps_fvm(:,:,ie))
140 140400 : max_local(ie,7) = MINVAL(moist_ps_fvm(:,:,ie))
141 226800 : min_local(ie,8) = MINVAL(elem(ie)%state%psdry(:,:))
142 226800 : max_local(ie,8) = MAXVAL(elem(ie)%state%psdry(:,:))
143 226800 : min_local(ie,9) = MINVAL(moist_ps(:,:,ie))
144 226800 : max_local(ie,9) = MAXVAL(moist_ps(:,:,ie))
145 43200 : do q=1,statediag_numtrac
146 32400 : varname(nm+q) = TRIM(cnst_name(q))
147 39204000 : min_local(ie,nm+q) = MINVAL(fvm(ie)%c(1:nc,1:nc,:,q))
148 39214800 : max_local(ie,nm+q) = MAXVAL(fvm(ie)%c(1:nc,1:nc,:,q))
149 : end do
150 : else
151 0 : min_local(ie,6) = MINVAL(elem(ie)%state%psdry(:,:))
152 0 : max_local(ie,6) = MAXVAL(elem(ie)%state%psdry(:,:))
153 0 : min_local(ie,7) = MINVAL(moist_ps(:,:,ie))
154 0 : max_local(ie,7) = MAXVAL(moist_ps(:,:,ie))
155 0 : do q=1,statediag_numtrac
156 0 : varname(nm+q) = TRIM(cnst_name(q))
157 0 : tmp_q = elem(ie)%state%Qdp(:,:,:,q,n0_qdp)/elem(ie)%state%dp3d(:,:,:,n0)
158 0 : min_local(ie,nm+q) = MINVAL(tmp_q)
159 0 : max_local(ie,nm+q) = MAXVAL(tmp_q)
160 : end do
161 : end if
162 : !
163 : ! forcing diagnostics
164 : !
165 10800 : varname(nm2+1) = 'FT '
166 10800 : varname(nm2+2) = 'FM '
167 21103200 : min_local(ie,nm2+1) = MINVAL(elem(ie)%derived%FT(:,:,:))
168 21103200 : max_local(ie,nm2+1) = MAXVAL(elem(ie)%derived%FT(:,:,:))
169 43200000 : min_local(ie,nm2+2) = MINVAL(elem(ie)%derived%FM(:,:,:,:))
170 43200000 : max_local(ie,nm2+2) = MAXVAL(elem(ie)%derived%FM(:,:,:,:))
171 12336 : if (use_cslam) then
172 43200 : do q=1,statediag_numtrac
173 32400 : varname(nm2+2+q) = TRIM('F'//TRIM(cnst_name(q)))
174 39204000 : min_local(ie,nm2+2+q) = MINVAL(fvm(ie)%fc(1:nc,1:nc,:,q))
175 39214800 : max_local(ie,nm2+2+q) = MAXVAL(fvm(ie)%fc(1:nc,1:nc,:,q))
176 : end do
177 : else
178 0 : do q=1,statediag_numtrac
179 0 : varname(nm2+2+q) = TRIM('F'//TRIM(cnst_name(q)))
180 0 : tmp_q = elem(ie)%derived%FQ(:,:,:,q)
181 0 : min_local(ie,nm2+2+q) = MINVAL(tmp_q)
182 0 : max_local(ie,nm2+2+q) = MAXVAL(tmp_q)
183 : end do
184 : end if
185 :
186 : end do
187 : !JMD This is a Thread Safe Reduction
188 27648 : do k = 1, nm2+2+statediag_numtrac
189 26112 : if (k==1) call t_startf('parallelMin')
190 26112 : min_p(k) = ParallelMin(min_local(:,k),hybrid)
191 26112 : if (k==1) call t_stopf('parallelMin')
192 27648 : max_p(k) = ParallelMax(max_local(:,k),hybrid)
193 : end do
194 : !
195 : !*********************************************
196 : !
197 : ! Mass diagnostics
198 : !
199 : !*********************************************
200 : !
201 : ! tracers
202 : !
203 52224 : mass = -1.0_r8
204 1536 : if (use_cslam) then
205 12336 : do ie=nets,nete
206 43200 : do q=1,statediag_numtrac
207 27550800 : tmp_fvm(:,:,q,ie) = SUM(fvm(ie)%c(1:nc,1:nc,:,q)*fvm(ie)%dp_fvm(1:nc,1:nc,:),DIM=3)
208 : end do
209 10800 : q=statediag_numtrac+1
210 9180000 : tmp_fvm(:,:,q,ie) = SUM(fvm(ie)%dp_fvm(1:nc,1:nc,:),DIM=3)
211 10800 : q=statediag_numtrac+2
212 141936 : tmp_fvm(:,:,q,ie) = moist_ps_fvm(:,:,ie)
213 : end do
214 0 : call global_integrals_general(tmp_fvm(:,:,1:statediag_numtrac+2,nets:nete),hybrid,nc,da_fvm,statediag_numtrac+2,&
215 714336 : nets,nete,tmp_mass(1:statediag_numtrac+2))
216 6144 : mass(nm+1:nm+statediag_numtrac)=tmp_mass(1:statediag_numtrac)*0.01_r8
217 4608 : mass(6:7)=tmp_mass(statediag_numtrac+1:statediag_numtrac+2)*0.01_r8
218 12336 : do ie=nets,nete
219 226800 : tmp_gll(:,:,1,ie)=elem(ie)%state%psdry(:,:)
220 228336 : tmp_gll(:,:,2,ie)=moist_ps(:,:,ie)
221 : end do
222 : call global_integrals_general(tmp_gll(:,:,1:2,nets:nete),hybrid,np,da_gll,2,&
223 465936 : nets,nete,tmp_mass(1:2))
224 4608 : mass(8:9)=tmp_mass(1:2)*0.01_r8
225 : else
226 0 : do ie=nets,nete
227 0 : do q=1,statediag_numtrac
228 0 : tmp_gll(:,:,q,ie)=sum(elem(ie)%state%Qdp(:,:,:,q,n0_qdp),DIM=3)
229 : end do
230 0 : q=statediag_numtrac+1
231 0 : tmp_gll(:,:,q,ie)=elem(ie)%state%psdry(:,:)
232 0 : q=statediag_numtrac+2
233 0 : tmp_gll(:,:,q,ie)=moist_ps(:,:,ie)
234 : end do
235 0 : call global_integrals_general(tmp_gll(:,:,1:statediag_numtrac+2,nets:nete),hybrid,np,da_gll,statediag_numtrac+2,&
236 0 : nets,nete,tmp_mass(1:statediag_numtrac+2))
237 0 : mass(nm+1:nm+statediag_numtrac)=tmp_mass(1:statediag_numtrac)*0.01_r8
238 0 : mass(6:7)=tmp_mass(statediag_numtrac+1:statediag_numtrac+2)*0.01_r8
239 : end if
240 : !
241 : ! compute relative mass change
242 : !
243 1536 : if (tl%nstep==0.or..not. initial_run) then
244 1536 : mass_chg(:) = 0.0_R8
245 6144 : elem(nets)%derived%mass(nm+1:nm+statediag_numtrac) = mass(nm+1:nm+statediag_numtrac)
246 1536 : if (use_cslam) then
247 7680 : elem(nets)%derived%mass(6:9) = mass(6:9)
248 : else
249 0 : elem(nets)%derived%mass(6:7) = mass(6:7)
250 : end if
251 : else
252 0 : mass_chg(:) = 0.0_r8
253 0 : do q=1,nm2!statediag_numtrac
254 0 : if (mass(q).ne.-1.0_r8) then
255 0 : if (ABS(elem(nets)%derived%mass(q))<1.0e-12_r8) then
256 0 : mass_chg(q) =mass(q) - elem(nets)%derived%mass(q)
257 : else
258 0 : mass_chg(q) =(mass(q) - elem(nets)%derived%mass(q))/elem(nets)%derived%mass(q)
259 : end if
260 : end if
261 : end do
262 : end if
263 : !
264 : ! write diagnostics to log file
265 : !
266 1536 : if(hybrid%masterthread) then
267 2 : write(iulog,*) ' '
268 2 : write(iulog,*) 'STATE DIAGNOSTICS'
269 2 : write(iulog,*) ' '
270 2 : write(iulog,101) ' ','MIN','MAX','AVE (hPa)','REL. MASS. CHANGE'
271 26 : do k=1,nm+statediag_numtrac
272 26 : if (mass(k)==-1.0_r8) then
273 10 : write(iulog,100) varname(k),min_p(k),max_p(k)
274 : else
275 14 : write(iulog,100) varname(k),min_p(k),max_p(k),mass(k),mass_chg(k)
276 : end if
277 : end do
278 : !
279 : ! forcing diagnostics
280 : !
281 2 : write(iulog,*) ' '
282 2 : write(iulog,*) 'FORCING DIAGNOSTICS'
283 2 : write(iulog,*) ' '
284 2 : write(iulog,101) ' ','MIN','MAX'
285 12 : do k=nm2+1,nm2+2+statediag_numtrac
286 12 : write(iulog,100) varname(k),min_p(k),max_p(k)
287 : end do
288 : end if
289 :
290 : 100 format (A12,4(E23.15))
291 : 101 format (A12,A23,A23,A23,A23)
292 :
293 : #ifdef waccm_debug
294 : call prim_printstate_cslam_gamma(elem, tl,hybrid,nets,nete, fvm)
295 : #endif
296 1536 : call prim_printstate_U(elem, tl,hybrid,nets,nete, fvm)
297 3072 : end subroutine prim_printstate
298 :
299 :
300 : #ifdef waccm_debug
301 : subroutine prim_printstate_cslam_gamma(elem, tl,hybrid,nets,nete, fvm)
302 : type (element_t), intent(inout) :: elem(:)
303 : type(fvm_struct), intent(inout) :: fvm(:)
304 : type (TimeLevel_t), target, intent(in) :: tl
305 : type (hybrid_t), intent(in) :: hybrid
306 : integer, intent(in) :: nets,nete
307 :
308 : ! Local variables...
309 : integer :: k,ie
310 :
311 : real (kind=r8), dimension(nets:nete,nlev) :: max_local
312 : real (kind=r8), dimension(nlev) :: max_p
313 : integer :: n0, n0_qdp, q, nm, nm2
314 :
315 : !dt=tstep*qsplit
316 : !dt = tstep*qsplit*rsplit ! vertical REMAP timestep
317 : !dynamics variables in n0 are at time = 'time': time=tl%nstep*tstep
318 : if (hybrid%masterthread) then
319 : write(iulog,*) "nstep=",tl%nstep," time=",Time_at(tl%nstep)/(24*3600)," [day]"
320 : end if
321 : ! dynamics timelevels
322 : n0=tl%n0
323 : call TimeLevel_Qdp( tl, qsplit, n0_qdp)
324 :
325 : do ie=nets,nete
326 : do k=1,nlev
327 : max_local(ie,k) = MAXVAL(fvm(ie)%CSLAM_gamma(:,:,k,1))
328 : end do
329 : end do
330 : !JMD This is a Thread Safe Reduction
331 : do k = 1, nlev
332 : max_p(k) = Parallelmax(max_local(:,k),hybrid)
333 : end do
334 : if (hybrid%masterthread) then
335 : write(iulog,*) ' '
336 : write(iulog,*) 'Gamma max'
337 : write(iulog,*) ' '
338 : do k=1,nlev
339 : write(iulog,*) 'k,gamma= ',k,max_p(k)
340 : end do
341 : end if
342 : end subroutine prim_printstate_cslam_gamma
343 : #endif
344 :
345 0 : subroutine adjust_nsplit(elem, tl,hybrid,nets,nete, fvm, omega_cn)
346 : use dimensions_mod, only: ksponge_end
347 1536 : use dimensions_mod, only: fvm_supercycling, fvm_supercycling_jet
348 : use se_dyn_time_mod, only: tstep
349 : use control_mod, only: rsplit, qsplit
350 : use perf_mod, only: t_startf, t_stopf
351 : use se_dyn_time_mod, only: nsplit, nsplit_baseline,rsplit_baseline
352 : use control_mod, only: qsplit, rsplit
353 : use time_manager, only: get_step_size
354 : use cam_abortutils, only: endrun
355 : use control_mod, only: nu_top
356 : !
357 : type (element_t), intent(inout) :: elem(:)
358 : type (TimeLevel_t), target, intent(in) :: tl
359 : type (hybrid_t), intent(in) :: hybrid
360 : integer, intent(in) :: nets,nete
361 : type(fvm_struct), intent(inout) :: fvm(:)
362 : real (kind=r8), intent(in) :: omega_cn(2,nets:nete)
363 : ! Local variables...
364 : integer :: k,ie
365 : real (kind=r8), dimension(1) :: min_o
366 : real (kind=r8), dimension(1) :: max_o
367 : real (kind=r8) :: dtime
368 : character(len=128) :: errmsg
369 : real (kind=r8) :: threshold=0.90_r8
370 0 : real (kind=r8) :: max_abs_omega_cn(nets:nete)
371 : real (kind=r8) :: min_abs_omega_cn(nets:nete)
372 : !
373 : ! The threshold values for when to double nsplit are empirical.
374 : ! In FW2000climo runs the Courant numbers are large in the sponge
375 : !
376 : ! The model was found to be stable if regular del4 is increased
377 : ! in the sponge and nu_top is increased (when nsplit doubles)
378 : !
379 : !
380 0 : do ie=nets,nete
381 0 : max_abs_omega_cn(ie) = MAXVAL(ABS(omega_cn(:,ie)))
382 : end do
383 :
384 : !JMD This is a Thread Safe Reduction
385 0 : do k = 1,1
386 0 : max_o(k) = ParallelMax(max_abs_omega_cn(:),hybrid)
387 : ! min_o(k) = ParallelMin(min_abs_omega_cn(:),hybrid)
388 : end do
389 0 : if (max_o(1)>threshold.and.nsplit==nsplit_baseline) then
390 : !
391 : ! change vertical remap time-step
392 : !
393 0 : nsplit=2*nsplit_baseline
394 0 : fvm_supercycling = rsplit
395 0 : fvm_supercycling_jet = rsplit
396 0 : nu_top=2.0_r8*nu_top
397 : !
398 : ! write diagnostics to log file
399 : !
400 0 : if(hybrid%masterthread) then
401 : !dynamics variables in n0 are at time = 'time': time=tl%nstep*tstep
402 : !dt=tstep*qsplit
403 : ! dt_remap = tstep*qsplit*rsplit ! vertical REMAP timestep
404 : !
405 0 : write(iulog,*) 'adj. nsplit: doubling nsplit; t=',Time_at(tl%nstep)/(24*3600)," [day]; max OMEGA",max_o(1)
406 : end if
407 0 : dtime = get_step_size()
408 0 : tstep = dtime / real(nsplit*qsplit*rsplit, r8)
409 :
410 0 : else if (nsplit.ne.nsplit_baseline.and.max_o(1)<0.4_r8*threshold) then
411 : !
412 : ! should nsplit be reduced again?
413 : !
414 0 : nsplit=nsplit_baseline
415 0 : rsplit=rsplit_baseline
416 0 : fvm_supercycling = rsplit
417 0 : fvm_supercycling_jet = rsplit
418 0 : nu_top=nu_top/2.0_r8
419 :
420 : ! nu_div_scale_top(:) = 1.0_r8
421 :
422 0 : dtime = get_step_size()
423 0 : tstep = dtime / real(nsplit*qsplit*rsplit, r8)
424 0 : if(hybrid%masterthread) then
425 0 : write(iulog,*) 'adj. nsplit: reset nsplit ; t=',Time_at(tl%nstep)/(24*3600)," [day]; max OMEGA",max_o(1)
426 : end if
427 : end if
428 0 : end subroutine adjust_nsplit
429 :
430 1536 : subroutine prim_printstate_U(elem, tl,hybrid,nets,nete, fvm)
431 : type (element_t), intent(inout) :: elem(:)
432 : type(fvm_struct), intent(inout) :: fvm(:)
433 : type (TimeLevel_t), target, intent(in) :: tl
434 : type (hybrid_t), intent(in) :: hybrid
435 : integer, intent(in) :: nets,nete
436 :
437 : ! Local variables...
438 : integer :: k,ie
439 :
440 3072 : real (kind=r8), dimension(nets:nete,nlev) :: max_local
441 3072 : real (kind=r8), dimension(nets:nete,nlev) :: min_local
442 : real (kind=r8), dimension(nlev) :: max_p
443 : real (kind=r8), dimension(nlev) :: min_p
444 : integer :: n0, n0_qdp, q, nm, nm2
445 :
446 : !dt=tstep*qsplit
447 : !dt = tstep*qsplit*rsplit ! vertical REMAP timestep
448 : !dynamics variables in n0 are at time = 'time': time=tl%nstep*tstep
449 1536 : if (hybrid%masterthread) then
450 2 : write(iulog,*) "nstep=",tl%nstep," time=",Time_at(tl%nstep)/(24*3600)," [day]"
451 : end if
452 : ! dynamics timelevels
453 1536 : n0=tl%n0
454 1536 : call TimeLevel_Qdp( tl, qsplit, n0_qdp)
455 :
456 12336 : do ie=nets,nete
457 1016736 : do k=1,nlev
458 43189200 : max_local(ie,k) = MAXVAL(elem(ie)%state%v(:,:,:,k,n0))
459 43200000 : min_local(ie,k) = MINVAL(elem(ie)%state%v(:,:,:,k,n0))
460 : end do
461 : end do
462 : !JMD This is a Thread Safe Reduction
463 144384 : do k = 1, nlev
464 142848 : max_p(k) = Parallelmax(max_local(:,k),hybrid)
465 144384 : min_p(k) = Parallelmin(min_local(:,k),hybrid)
466 : end do
467 1536 : if (hybrid%masterthread) then
468 2 : write(iulog,*) ' '
469 2 : write(iulog,*) 'min/max of wind components in each layer'
470 2 : write(iulog,*) ' '
471 188 : do k=1,nlev
472 188 : write(iulog,*) 'k,V (min max)= ',k,min_p(k),max_p(k)
473 : end do
474 : end if
475 0 : end subroutine prim_printstate_U
476 : end module prim_state_mod
|