Line data Source code
1 : module global_norms_mod
2 :
3 : use shr_kind_mod, only: r8=>shr_kind_r8
4 : use cam_logfile, only: iulog
5 : use edgetype_mod, only: EdgeBuffer_t
6 :
7 : implicit none
8 : private
9 : save
10 :
11 : public :: l1_snorm
12 : public :: l2_snorm
13 : public :: linf_snorm
14 :
15 : public :: l1_vnorm
16 : public :: l2_vnorm
17 : public :: linf_vnorm
18 :
19 : public :: print_cfl
20 : public :: global_integral
21 : public :: global_integrals_general
22 : public :: wrap_repro_sum
23 :
24 : private :: global_maximum
25 : type (EdgeBuffer_t), private :: edgebuf
26 :
27 : interface global_integral
28 : module procedure global_integral_elem
29 : module procedure global_integral_fvm
30 : end interface global_integral
31 :
32 : contains
33 :
34 :
35 : subroutine global_integrals(elem,fld,hybrid,npts,num_flds,nets,nete,I_sphere)
36 : use hybrid_mod, only: hybrid_t
37 : use element_mod, only: element_t
38 : use dimensions_mod, only: np
39 : use physconst, only: pi
40 : use parallel_mod, only: global_shared_buf, global_shared_sum
41 :
42 : type(element_t) , intent(in) :: elem(:)
43 : integer , intent(in) :: npts,nets,nete,num_flds
44 : real (kind=r8), intent(in) :: fld(npts,npts,num_flds,nets:nete)
45 : type (hybrid_t) , intent(in) :: hybrid
46 :
47 : real (kind=r8) :: I_sphere(num_flds)
48 : !
49 : ! Local variables
50 : !
51 : integer :: ie,j,i,q
52 :
53 : real (kind=r8) :: da
54 : real (kind=r8) :: J_tmp(nets:nete,num_flds)
55 : !
56 : ! This algorithm is independent of thread count and task count.
57 : ! This is a requirement of consistancy checking in cam.
58 : !
59 : J_tmp = 0.0_r8
60 :
61 : do ie=nets,nete
62 : do q=1,num_flds
63 : do j=1,np
64 : do i=1,np
65 : da = elem(ie)%mp(i,j)*elem(ie)%metdet(i,j)
66 : J_tmp(ie,q) = J_tmp(ie,q) + da*fld(i,j,q,ie)
67 : end do
68 : end do
69 : end do
70 : end do
71 : do ie=nets,nete
72 : global_shared_buf(ie,1:num_flds) = J_tmp(ie,:)
73 : enddo
74 : call wrap_repro_sum(nvars=num_flds, comm=hybrid%par%comm)
75 : I_sphere(:) =global_shared_sum(1:num_flds) /(4.0_r8*PI)
76 : end subroutine global_integrals
77 :
78 3072 : subroutine global_integrals_general(fld,hybrid,npts,da,num_flds,nets,nete,I_sphere)
79 : use hybrid_mod, only: hybrid_t
80 : use physconst, only: pi
81 : use parallel_mod, only: global_shared_buf, global_shared_sum
82 :
83 : integer, intent(in) :: npts,nets,nete,num_flds
84 : real (kind=r8), intent(in) :: fld(npts,npts,num_flds,nets:nete)
85 : type (hybrid_t), intent(in) :: hybrid
86 : real (kind=r8), intent(in) :: da(npts,npts,nets:nete)
87 :
88 : real (kind=r8) :: I_sphere(num_flds)
89 : !
90 : ! Local variables
91 : !
92 : integer :: ie,j,i,q
93 :
94 6144 : real (kind=r8) :: J_tmp(nets:nete,num_flds)
95 : !
96 : ! This algorithm is independent of thread count and task count.
97 : ! This is a requirement of consistancy checking in cam.
98 : !
99 89424 : J_tmp = 0.0_r8
100 :
101 24672 : do ie=nets,nete
102 100272 : do q=1,num_flds
103 345600 : do j=1,npts
104 1155600 : do i=1,npts
105 1080000 : J_tmp(ie,q) = J_tmp(ie,q) + da(i,j,ie)*fld(i,j,q,ie)
106 : end do
107 : end do
108 : end do
109 : end do
110 24672 : do ie=nets,nete
111 100272 : global_shared_buf(ie,1:num_flds) = J_tmp(ie,:)
112 : enddo
113 3072 : call wrap_repro_sum(nvars=num_flds, comm=hybrid%par%comm)
114 16896 : I_sphere(:) =global_shared_sum(1:num_flds) /(4.0_r8*PI)
115 3072 : end subroutine global_integrals_general
116 :
117 :
118 : ! ================================
119 : ! global_integral:
120 : !
121 : ! eq 81 in Williamson, et. al. p 218
122 : ! for spectral elements
123 : !
124 : ! ================================
125 : ! --------------------------
126 2304 : function global_integral_elem(elem,fld,hybrid,npts,nets,nete) result(I_sphere)
127 : use hybrid_mod, only: hybrid_t
128 : use element_mod, only: element_t
129 : use dimensions_mod, only: np
130 : use physconst, only: pi
131 : use parallel_mod, only: global_shared_buf, global_shared_sum
132 :
133 : type(element_t) , intent(in) :: elem(:)
134 : integer , intent(in) :: npts,nets,nete
135 : real (kind=r8), intent(in) :: fld(npts,npts,nets:nete)
136 : type (hybrid_t) , intent(in) :: hybrid
137 :
138 : real (kind=r8) :: I_sphere
139 :
140 : ! Local variables
141 :
142 : integer :: ie,j,i
143 : real(kind=r8) :: I_tmp(1)
144 :
145 : real (kind=r8) :: da
146 2304 : real (kind=r8) :: J_tmp(nets:nete)
147 : !
148 : ! This algorithm is independent of thread count and task count.
149 : ! This is a requirement of consistancy checking in cam.
150 : !
151 18504 : J_tmp = 0.0_r8
152 :
153 18504 : do ie=nets,nete
154 83304 : do j=1,np
155 340200 : do i=1,np
156 259200 : da = elem(ie)%mp(i,j)*elem(ie)%metdet(i,j)
157 324000 : J_tmp(ie) = J_tmp(ie) + da*fld(i,j,ie)
158 : end do
159 : end do
160 : end do
161 18504 : do ie=nets,nete
162 18504 : global_shared_buf(ie,1) = J_tmp(ie)
163 : enddo
164 2304 : call wrap_repro_sum(nvars=1, comm=hybrid%par%comm)
165 4608 : I_tmp = global_shared_sum(1)
166 2304 : I_sphere = I_tmp(1)/(4.0_r8*PI)
167 :
168 2304 : end function global_integral_elem
169 :
170 0 : function global_integral_fvm(fvm,fld,hybrid,npts,nets,nete) result(I_sphere)
171 : use hybrid_mod, only: hybrid_t
172 : use fvm_control_volume_mod, only: fvm_struct
173 : use physconst, only: pi
174 : use parallel_mod, only: global_shared_buf, global_shared_sum
175 :
176 : type (fvm_struct) , intent(in) :: fvm(:)
177 : integer , intent(in) :: npts,nets,nete
178 : real (kind=r8), intent(in) :: fld(npts,npts,nets:nete)
179 : type (hybrid_t) , intent(in) :: hybrid
180 :
181 : real (kind=r8) :: I_sphere
182 :
183 : ! Local variables
184 :
185 : integer :: ie,j,i
186 : real(kind=r8) :: I_tmp(1)
187 :
188 : real (kind=r8) :: da
189 0 : real (kind=r8) :: J_tmp(nets:nete)
190 : !
191 : ! This algorithm is independent of thread count and task count.
192 : ! This is a requirement of consistancy checking in cam.
193 : !
194 0 : J_tmp = 0.0_r8
195 0 : do ie=nets,nete
196 0 : do j=1,npts
197 0 : do i=1,npts
198 0 : da = fvm(ie)%area_sphere(i,j)
199 0 : J_tmp(ie) = J_tmp(ie) + da*fld(i,j,ie)
200 : end do
201 : end do
202 : end do
203 0 : do ie=nets,nete
204 0 : global_shared_buf(ie,1) = J_tmp(ie)
205 : enddo
206 0 : call wrap_repro_sum(nvars=1, comm=hybrid%par%comm)
207 0 : I_tmp = global_shared_sum(1)
208 0 : I_sphere = I_tmp(1)/(4.0_r8*PI)
209 :
210 0 : end function global_integral_fvm
211 :
212 : !------------------------------------------------------------------------------------
213 :
214 : ! ================================
215 : ! print_cfl:
216 : !
217 : ! Calculate / output CFL info
218 : ! (both advective and based on
219 : ! viscosity or hyperviscosity)
220 : !
221 : ! ================================
222 :
223 1536 : subroutine print_cfl(elem,hybrid,nets,nete,dtnu,ptop,pmid,&
224 : dt_remap_actual,dt_tracer_fvm_actual,dt_tracer_se_actual,&
225 : dt_dyn_actual,dt_dyn_visco_actual,dt_dyn_del2_actual,dt_tracer_visco_actual,dt_phys)
226 : !
227 : ! estimate various CFL limits
228 : ! also, for variable resolution viscosity coefficient, make sure
229 : ! worse viscosity CFL (given by dtnu) is not violated by reducing
230 : ! viscosity coefficient in regions where CFL is violated
231 : !
232 : use hybrid_mod, only: hybrid_t
233 : use element_mod, only: element_t
234 : use dimensions_mod, only: np,ne,nelem,nc,nhe,use_cslam,nlev,large_Courant_incr
235 : use dimensions_mod, only: nu_scale_top,nu_div_lev,nu_lev,nu_t_lev
236 :
237 : use quadrature_mod, only: gausslobatto, quadrature_t
238 :
239 : use reduction_mod, only: ParallelMin,ParallelMax
240 : use physconst, only: ra, rearth, pi
241 : use control_mod, only: nu, nu_div, nu_q, nu_p, nu_t, nu_top, fine_ne, max_hypervis_courant
242 : use control_mod, only: tstep_type, hypervis_power, hypervis_scaling
243 : use control_mod, only: sponge_del4_nu_div_fac, sponge_del4_nu_fac, sponge_del4_lev
244 : use cam_abortutils, only: endrun
245 : use parallel_mod, only: global_shared_buf, global_shared_sum
246 : use edge_mod, only: initedgebuffer, FreeEdgeBuffer, edgeVpack, edgeVunpack
247 : use bndry_mod, only: bndry_exchange
248 : use mesh_mod, only: MeshUseMeshFile
249 : use dimensions_mod, only: ksponge_end, kmvis_ref, kmcnd_ref,rho_ref
250 : use physconst, only: cpair
251 : use std_atm_profile,only: std_atm_height
252 : type(element_t) , intent(inout) :: elem(:)
253 : integer , intent(in) :: nets,nete
254 : type (hybrid_t) , intent(in) :: hybrid
255 : real (kind=r8), intent(in) :: dtnu, ptop, pmid(nlev)
256 : !
257 : ! actual time-steps
258 : !
259 : real (kind=r8), intent(in) :: dt_remap_actual,dt_tracer_fvm_actual,dt_tracer_se_actual,&
260 : dt_dyn_actual,dt_dyn_visco_actual,dt_dyn_del2_actual, &
261 : dt_tracer_visco_actual, dt_phys
262 :
263 : ! Element statisics
264 : real (kind=r8) :: max_min_dx,min_min_dx,min_max_dx,max_unif_dx ! used for normalizing scalar HV
265 : real (kind=r8) :: max_normDinv, min_normDinv ! used for CFL
266 : real (kind=r8) :: min_area, max_area,max_ratio !min/max element area
267 : real (kind=r8) :: avg_area, avg_min_dx,tot_area,tot_area_rad
268 : real (kind=r8) :: min_hypervis, max_hypervis, avg_hypervis, stable_hv
269 : real (kind=r8) :: normDinv_hypervis
270 : real (kind=r8) :: x, y, noreast, nw, se, sw
271 3072 : real (kind=r8), dimension(np,np,nets:nete) :: zeta
272 : real (kind=r8) :: lambda_max, lambda_vis, min_gw, lambda,umax, ugw
273 : real (kind=r8) :: scale1, max_laplace,z(nlev)
274 : integer :: ie, i, j, rowind, colind, k
275 : type (quadrature_t) :: gp
276 : character(LEN=256) :: rk_str
277 :
278 : real (kind=r8) :: s_laplacian, s_hypervis, s_rk, s_rk_tracer !Stability region
279 : real (kind=r8) :: dt_max_adv, dt_max_gw, dt_max_tracer_se, dt_max_tracer_fvm
280 : real (kind=r8) :: dt_max_hypervis, dt_max_hypervis_tracer, dt_max_laplacian_top
281 :
282 : real(kind=r8) :: I_sphere, nu_max, nu_div_max
283 3072 : real(kind=r8) :: fld(np,np,nets:nete)
284 :
285 : logical :: top_000_032km, top_032_042km, top_042_090km, top_090_140km, top_140_600km ! model top location ranges
286 : logical :: nu_set,div_set,lev_set
287 :
288 : ! Eigenvalues calculated by folks at UMich (Paul U & Jared W)
289 : select case (np)
290 : case (2)
291 : lambda_max = 0.5_r8
292 : lambda_vis = 0.0_r8 ! need to compute this
293 : case (3)
294 : lambda_max = 1.5_r8
295 : lambda_vis = 12.0_r8
296 : case (4)
297 1536 : lambda_max = 2.74_r8
298 1536 : lambda_vis = 30.0_r8
299 : case (5)
300 : lambda_max = 4.18_r8
301 : lambda_vis = 91.6742_r8
302 : case (6)
303 : lambda_max = 5.86_r8
304 : lambda_vis = 190.1176_r8
305 : case (7)
306 : lambda_max = 7.79_r8
307 : lambda_vis = 374.7788_r8
308 : case (8)
309 : lambda_max = 10.0_r8
310 : lambda_vis = 652.3015_r8
311 : case DEFAULT
312 : lambda_max = 0.0_r8
313 : lambda_vis = 0.0_r8
314 : end select
315 :
316 1536 : if ((lambda_max.eq.0_r8).and.(hybrid%masterthread)) then
317 : print*, "lambda_max not calculated for NP = ",np
318 : print*, "Estimate of gravity wave timestep will be incorrect"
319 : end if
320 : if ((lambda_vis.eq.0_r8).and.(hybrid%masterthread)) then
321 : print*, "lambda_vis not calculated for NP = ",np
322 : print*, "Estimate of viscous CFLs will be incorrect"
323 : end if
324 :
325 12336 : do ie=nets,nete
326 228336 : elem(ie)%variable_hyperviscosity = 1.0_r8
327 : end do
328 :
329 1536 : gp=gausslobatto(np)
330 9216 : min_gw = minval(gp%weights)
331 : !
332 : !******************************************************************************************
333 : !
334 : ! compute some local and global grid metrics
335 : !
336 : !******************************************************************************************
337 : !
338 228336 : fld(:,:,nets:nete)=1.0_r8
339 : ! Calculate surface area by integrating 1.0_r8 over sphere and dividing by 4*PI (Should be 1)
340 1536 : I_sphere = global_integral(elem, fld(:,:,nets:nete),hybrid,np,nets,nete)
341 :
342 1536 : min_normDinv = 1E99_r8
343 1536 : max_normDinv = 0
344 1536 : min_max_dx = 1E99_r8
345 1536 : min_min_dx = 1E99_r8
346 1536 : max_min_dx = 0
347 1536 : min_area = 1E99_r8
348 1536 : max_area = 0
349 1536 : max_ratio = 0
350 12336 : do ie=nets,nete
351 10800 : max_normDinv = max(max_normDinv,elem(ie)%normDinv)
352 10800 : min_normDinv = min(min_normDinv,elem(ie)%normDinv)
353 10800 : min_min_dx = min(min_min_dx,elem(ie)%dx_short)
354 10800 : max_min_dx = max(max_min_dx,elem(ie)%dx_short)
355 10800 : min_max_dx = min(min_max_dx,elem(ie)%dx_long)
356 :
357 226800 : elem(ie)%area = sum(elem(ie)%spheremp(:,:))
358 10800 : min_area = min(min_area,elem(ie)%area)
359 10800 : max_area = max(max_area,elem(ie)%area)
360 10800 : max_ratio = max(max_ratio,elem(ie)%dx_long/elem(ie)%dx_short)
361 :
362 10800 : global_shared_buf(ie,1) = elem(ie)%area
363 12336 : global_shared_buf(ie,2) = elem(ie)%dx_short
364 : enddo
365 1536 : call wrap_repro_sum(nvars=2, comm=hybrid%par%comm)
366 1536 : avg_area = global_shared_sum(1)/dble(nelem)
367 1536 : tot_area_rad = global_shared_sum(1)
368 1536 : avg_min_dx = global_shared_sum(2)/dble(nelem)
369 :
370 1536 : min_area = ParallelMin(min_area,hybrid)
371 1536 : max_area = ParallelMax(max_area,hybrid)
372 1536 : min_normDinv = ParallelMin(min_normDinv,hybrid)
373 1536 : max_normDinv = ParallelMax(max_normDinv,hybrid)
374 1536 : min_min_dx = ParallelMin(min_min_dx,hybrid)
375 1536 : max_min_dx = ParallelMax(max_min_dx,hybrid)
376 1536 : min_max_dx = ParallelMin(min_max_dx,hybrid)
377 1536 : max_ratio = ParallelMax(max_ratio,hybrid)
378 : ! Physical units for area (unit sphere to Earth sphere)
379 1536 : min_area = min_area*rearth*rearth/1000000._r8 !m2 (rearth is in units of km)
380 1536 : max_area = max_area*rearth*rearth/1000000._r8 !m2 (rearth is in units of km)
381 1536 : avg_area = avg_area*rearth*rearth/1000000._r8 !m2 (rearth is in units of km)
382 1536 : tot_area = tot_area_rad*rearth*rearth/1000000._r8!m2 (rearth is in units of km)
383 1536 : if (hybrid%masterthread) then
384 2 : write(iulog,* )""
385 2 : write(iulog,* )"Running Global Integral Diagnostic..."
386 2 : write(iulog,*)"Area of unit sphere is",I_sphere
387 2 : write(iulog,*)"Should be 1.0 to round off..."
388 2 : write(iulog,'(a,f9.3)') 'Element area: max/min',(max_area/min_area)
389 2 : write(iulog,'(a,E23.15)') 'Total Grid area: ',(tot_area)
390 2 : write(iulog,'(a,E23.15)') 'Total Grid area rad^2: ',(tot_area_rad)
391 2 : if (.not.MeshUseMeshFile) then
392 2 : write(iulog,'(a,f6.3,f8.2)') "Average equatorial node spacing (deg, km) = ", &
393 4 : dble(90)/dble(ne*(np-1)), PI*rearth/(2000.0_r8*dble(ne*(np-1)))
394 : end if
395 2 : write(iulog,'(a,2f9.3)') 'norm of Dinv (min, max): ', min_normDinv, max_normDinv
396 2 : write(iulog,'(a,1f8.2)') 'Max Dinv-based element distortion: ', max_ratio
397 2 : write(iulog,'(a,3f8.2)') 'dx based on Dinv svd: ave,min,max = ', avg_min_dx, min_min_dx, max_min_dx
398 2 : write(iulog,'(a,3f8.2)') "dx based on sqrt element area: ave,min,max = ", &
399 4 : sqrt(avg_area)/(np-1),sqrt(min_area)/(np-1),sqrt(max_area)/(np-1)
400 : end if
401 :
402 :
403 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
404 : ! SCALAR, RESOLUTION-AWARE HYPERVISCOSITY
405 : ! this block of code initializes the variable_hyperviscsoity() array
406 : ! based on largest length scale in each element and user specified scaling
407 : ! it then limits the coefficient if the user specifed a max CFL
408 : ! this limiting is based on the smallest length scale of each element
409 : ! since that controls the CFL.
410 : ! Mike Levy
411 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
412 1536 : if (hypervis_power /= 0) then
413 :
414 0 : min_hypervis = 1d99
415 0 : max_hypervis = 0
416 0 : avg_hypervis = 0
417 :
418 :
419 0 : max_unif_dx = min_max_dx ! use this for average resolution, unless:
420 : ! viscosity in namelist specified for smallest element:
421 0 : if (fine_ne>0) then
422 : ! viscosity in namelist specified for regions with a resolution
423 : ! equivilant to a uniform grid with ne=fine_ne
424 : if (np /= 4 ) call endrun('ERROR: setting fine_ne only supported with NP=4')
425 0 : max_unif_dx = (111.28_r8*30)/dble(fine_ne) ! in km
426 : endif
427 :
428 : !
429 : ! note: if L = eigenvalue of metinv, then associated length scale (km) is
430 : ! dx = 1.0_r8/( sqrt(L)*0.5_r8*dble(np-1)*ra*1000.0_r8)
431 : !
432 : ! for viscosity *tensor*, we take at each point:
433 : ! nu1 = nu*(dx1/max_unif_dx)**3.2 dx1 associated with eigenvalue 1
434 : ! nu2 = nu*(dx2/max_unif_dx)**3.2 dx2 associated with eigenvalue 2
435 : ! with this approach:
436 : ! - with this formula, no need to adjust for CFL violations
437 : ! - if nu comes from a 3.2 scaling that is stable for coarse and fine resolutions,
438 : ! this formulat will be stable.
439 : ! - gives the correct answer in long skinny rectangles:
440 : ! large viscosity in the long direction, small viscosity in the short direction
441 : !
442 0 : normDinv_hypervis = 0
443 0 : do ie=nets,nete
444 : ! variable viscosity based on map from ulatlon -> ucontra
445 :
446 : ! dx_long
447 0 : elem(ie)%variable_hyperviscosity = sqrt((elem(ie)%dx_long/max_unif_dx) ** hypervis_power)
448 : elem(ie)%hv_courant = dtnu*(elem(ie)%variable_hyperviscosity(1,1)**2) * &
449 0 : (lambda_vis**2) * ((ra*elem(ie)%normDinv)**4)
450 :
451 : ! Check to see if this is stable
452 0 : if (elem(ie)%hv_courant.gt.max_hypervis_courant) then
453 : stable_hv = sqrt( max_hypervis_courant / &
454 0 : ( dtnu * (lambda_vis)**2 * (ra*elem(ie)%normDinv)**4 ) )
455 :
456 : #if 0
457 : ! Useful print statements for debugging the adjustments to hypervis
458 : print*, "Adjusting hypervis on elem ", elem(ie)%GlobalId
459 : print*, "From ", nu*elem(ie)%variable_hyperviscosity(1,1)**2, " to ", nu*stable_hv
460 : print*, "Difference = ", nu*(/elem(ie)%variable_hyperviscosity(1,1)**2-stable_hv/)
461 : print*, "Factor of ", elem(ie)%variable_hyperviscosity(1,1)**2/stable_hv
462 : print*, " "
463 : #endif
464 : ! make sure that: elem(ie)%hv_courant <= max_hypervis_courant
465 0 : elem(ie)%variable_hyperviscosity = stable_hv
466 0 : elem(ie)%hv_courant = dtnu*(stable_hv**2) * (lambda_vis)**2 * (ra*elem(ie)%normDinv)**4
467 : end if
468 0 : normDinv_hypervis = max(normDinv_hypervis, elem(ie)%hv_courant/dtnu)
469 :
470 0 : min_hypervis = min(min_hypervis, elem(ie)%variable_hyperviscosity(1,1))
471 0 : max_hypervis = max(max_hypervis, elem(ie)%variable_hyperviscosity(1,1))
472 0 : global_shared_buf(ie,1) = elem(ie)%variable_hyperviscosity(1,1)
473 : end do
474 :
475 0 : min_hypervis = ParallelMin(min_hypervis, hybrid)
476 0 : max_hypervis = ParallelMax(max_hypervis, hybrid)
477 0 : call wrap_repro_sum(nvars=1, comm=hybrid%par%comm)
478 0 : avg_hypervis = global_shared_sum(1)/dble(nelem)
479 :
480 0 : normDinv_hypervis = ParallelMax(normDinv_hypervis, hybrid)
481 :
482 : ! apply DSS (aka assembly procedure) to variable_hyperviscosity (makes continuous)
483 0 : call initEdgeBuffer(hybrid%par,edgebuf,elem,1)
484 0 : do ie=nets,nete
485 0 : zeta(:,:,ie) = elem(ie)%variable_hyperviscosity(:,:)*elem(ie)%spheremp(:,:)
486 0 : call edgeVpack(edgebuf,zeta(1,1,ie),1,0,ie)
487 : end do
488 0 : call bndry_exchange(hybrid,edgebuf,location='print_cfl #1')
489 0 : do ie=nets,nete
490 0 : call edgeVunpack(edgebuf,zeta(1,1,ie),1,0,ie)
491 0 : elem(ie)%variable_hyperviscosity(:,:) = zeta(:,:,ie)*elem(ie)%rspheremp(:,:)
492 : end do
493 0 : call FreeEdgeBuffer(edgebuf)
494 :
495 : ! replace hypervis w/ bilinear based on continuous corner values
496 0 : do ie=nets,nete
497 0 : noreast = elem(ie)%variable_hyperviscosity(np,np)
498 0 : nw = elem(ie)%variable_hyperviscosity(1,np)
499 0 : se = elem(ie)%variable_hyperviscosity(np,1)
500 0 : sw = elem(ie)%variable_hyperviscosity(1,1)
501 0 : do i=1,np
502 0 : x = gp%points(i)
503 0 : do j=1,np
504 0 : y = gp%points(j)
505 0 : elem(ie)%variable_hyperviscosity(i,j) = 0.25_r8*( &
506 : (1.0_r8-x)*(1.0_r8-y)*sw + &
507 : (1.0_r8-x)*(y+1.0_r8)*nw + &
508 : (x+1.0_r8)*(1.0_r8-y)*se + &
509 0 : (x+1.0_r8)*(y+1.0_r8)*noreast)
510 : end do
511 : end do
512 : end do
513 1536 : else if (hypervis_scaling/=0) then
514 : ! tensorHV. New eigenvalues are the eigenvalues of the tensor V
515 : ! formulas here must match what is in cube_mod.F90
516 : ! for tensorHV, we scale out the rearth dependency
517 0 : lambda = max_normDinv**2
518 : normDinv_hypervis = (lambda_vis**2) * (max_normDinv**4) * &
519 0 : (lambda**(-hypervis_scaling/2) )
520 : else
521 : ! constant coefficient formula:
522 1536 : normDinv_hypervis = (lambda_vis**2) * (ra*max_normDinv)**4
523 : endif
524 :
525 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
526 : ! TENSOR, RESOLUTION-AWARE HYPERVISCOSITY
527 : ! The tensorVisc() array is computed in cube_mod.F90
528 : ! this block of code will DSS it so the tensor if C0
529 : ! and also make it bilinear in each element.
530 : ! Oksana Guba
531 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
532 1536 : if (hypervis_scaling /= 0) then
533 :
534 0 : call initEdgeBuffer(hybrid%par,edgebuf,elem,1)
535 0 : do rowind=1,2
536 0 : do colind=1,2
537 0 : do ie=nets,nete
538 0 : zeta(:,:,ie) = elem(ie)%tensorVisc(:,:,rowind,colind)*elem(ie)%spheremp(:,:)
539 0 : call edgeVpack(edgebuf,zeta(1,1,ie),1,0,ie)
540 : end do
541 :
542 0 : call bndry_exchange(hybrid,edgebuf)
543 0 : do ie=nets,nete
544 0 : call edgeVunpack(edgebuf,zeta(1,1,ie),1,0,ie)
545 0 : elem(ie)%tensorVisc(:,:,rowind,colind) = zeta(:,:,ie)*elem(ie)%rspheremp(:,:)
546 : end do
547 : enddo !rowind
548 : enddo !colind
549 0 : call FreeEdgeBuffer(edgebuf)
550 :
551 : !IF BILINEAR MAP OF V NEEDED
552 :
553 0 : do rowind=1,2
554 0 : do colind=1,2
555 : ! replace hypervis w/ bilinear based on continuous corner values
556 0 : do ie=nets,nete
557 0 : noreast = elem(ie)%tensorVisc(np,np,rowind,colind)
558 0 : nw = elem(ie)%tensorVisc(1,np,rowind,colind)
559 0 : se = elem(ie)%tensorVisc(np,1,rowind,colind)
560 0 : sw = elem(ie)%tensorVisc(1,1,rowind,colind)
561 0 : do i=1,np
562 0 : x = gp%points(i)
563 0 : do j=1,np
564 0 : y = gp%points(j)
565 0 : elem(ie)%tensorVisc(i,j,rowind,colind) = 0.25_r8*( &
566 : (1.0_r8-x)*(1.0_r8-y)*sw + &
567 : (1.0_r8-x)*(y+1.0_r8)*nw + &
568 : (x+1.0_r8)*(1.0_r8-y)*se + &
569 0 : (x+1.0_r8)*(y+1.0_r8)*noreast)
570 : end do
571 : end do
572 : end do
573 : enddo !rowind
574 : enddo !colind
575 : endif
576 1536 : deallocate(gp%points)
577 1536 : deallocate(gp%weights)
578 :
579 1536 : call automatically_set_viscosity_coefficients(hybrid,ne,max_min_dx,min_min_dx,nu_p ,1.0_r8 ,'_p ')
580 1536 : call automatically_set_viscosity_coefficients(hybrid,ne,max_min_dx,min_min_dx,nu ,1.0_r8,' ')
581 1536 : call automatically_set_viscosity_coefficients(hybrid,ne,max_min_dx,min_min_dx,nu_div,2.5_r8 ,'_div')
582 :
583 1536 : if (nu_q<0) nu_q = nu_p ! necessary for consistency
584 1536 : if (nu_t<0) nu_t = nu_p ! temperature damping is always equal to nu_p
585 :
586 144384 : nu_div_lev(:) = nu_div
587 144384 : nu_lev(:) = nu
588 144384 : nu_t_lev(:) = nu_p
589 :
590 : !
591 : ! sponge layer strength needed for stability depends on model top location
592 : !
593 1536 : top_000_032km = .false.
594 1536 : top_032_042km = .false.
595 1536 : top_042_090km = .false.
596 1536 : top_090_140km = .false.
597 1536 : top_140_600km = .false.
598 1536 : nu_set = sponge_del4_nu_fac < 0
599 1536 : div_set = sponge_del4_nu_div_fac < 0
600 1536 : lev_set = sponge_del4_lev < 0
601 1536 : if (ptop>1000.0_r8) then
602 : !
603 : ! low top; usually idealized test cases
604 : !
605 0 : top_000_032km = .true.
606 0 : if (hybrid%masterthread) write(iulog,* )"Model top damping configuration: top_000_032km"
607 1536 : else if (ptop>100.0_r8) then
608 : !
609 : ! CAM6 top (~225 Pa) or CAM7 low top
610 : !
611 0 : top_032_042km = .true.
612 0 : if (hybrid%masterthread) write(iulog,* )"Model top damping configuration: top_032_042km"
613 1536 : else if (ptop>1e-1_r8) then
614 : !
615 : ! CAM7 top (~4.35e-1 Pa)
616 : !
617 1536 : top_042_090km = .true.
618 1536 : if (hybrid%masterthread) write(iulog,* )"Model top damping configuration: top_042_090km"
619 0 : else if (ptop>1E-4_r8) then
620 : !
621 : ! WACCM top (~4.5e-4 Pa)
622 : !
623 0 : top_090_140km = .true.
624 0 : if (hybrid%masterthread) write(iulog,* )"Model top damping configuration: top_090_140km"
625 : else
626 : !
627 : ! WACCM-x - geospace (~4e-7 Pa)
628 : !
629 0 : top_140_600km = .true.
630 0 : if (hybrid%masterthread) write(iulog,* )"Model top damping configuration: top_140_600km"
631 : end if
632 : !
633 : ! Logging text for sponge layer configuration
634 : !
635 1536 : if (hybrid%masterthread .and. (nu_set .or. div_set .or. lev_set)) then
636 2 : write(iulog,* )""
637 2 : write(iulog,* )"Sponge layer del4 coefficient defaults based on model top location:"
638 : end if
639 : !
640 : ! if user or namelist is not specifying sponge del4 settings here are best guesses (empirically determined)
641 : !
642 1536 : if (top_090_140km.or.top_140_600km) then ! defaults for waccm(x)
643 0 : if (sponge_del4_lev <0) sponge_del4_lev = 20
644 0 : if (sponge_del4_nu_fac <0) sponge_del4_nu_fac = 5.0_r8
645 0 : if (sponge_del4_nu_div_fac<0) sponge_del4_nu_div_fac = 10.0_r8
646 : else
647 1536 : if (sponge_del4_lev <0) sponge_del4_lev = 1
648 1536 : if (sponge_del4_nu_fac <0) sponge_del4_nu_fac = 1.0_r8
649 1536 : if (sponge_del4_nu_div_fac<0) sponge_del4_nu_div_fac = 1.0_r8
650 : end if
651 :
652 : ! set max wind speed for diagnostics
653 1536 : umax = 120.0_r8
654 1536 : if (top_042_090km) then
655 1536 : umax = 240._r8
656 0 : else if (top_090_140km) then
657 0 : umax = 300._r8
658 0 : else if (top_140_600km) then
659 0 : umax = 800._r8
660 : end if
661 : !
662 : ! Log sponge layer configuration
663 : !
664 1536 : if (hybrid%masterthread) then
665 2 : if (nu_set) then
666 2 : write(iulog, '(a,e9.2)') ' sponge_del4_nu_fac = ',sponge_del4_nu_fac
667 : end if
668 :
669 2 : if (div_set) then
670 2 : write(iulog, '(a,e9.2)') ' sponge_del4_nu_div_fac = ',sponge_del4_nu_div_fac
671 : end if
672 :
673 2 : if (lev_set) then
674 2 : write(iulog, '(a,i0)') ' sponge_del4_lev = ',sponge_del4_lev
675 : end if
676 2 : write(iulog,* )""
677 : end if
678 :
679 1536 : nu_max = sponge_del4_nu_fac*nu_p
680 1536 : nu_div_max = sponge_del4_nu_div_fac*nu_p
681 144384 : do k=1,nlev
682 : ! Vertical profile from FV dycore (see Lauritzen et al. 2012 DOI:10.1177/1094342011410088)
683 142848 : scale1 = 0.5_r8*(1.0_r8+tanh(2.0_r8*log(pmid(sponge_del4_lev)/pmid(k))))
684 142848 : if (sponge_del4_nu_div_fac /= 1.0_r8) then
685 0 : nu_div_lev(k) = (1.0_r8-scale1)*nu_div+scale1*nu_div_max
686 : end if
687 144384 : if (sponge_del4_nu_fac /= 1.0_r8) then
688 0 : nu_lev(k) = (1.0_r8-scale1)*nu +scale1*nu_max
689 0 : nu_t_lev(k) = (1.0_r8-scale1)*nu_p +scale1*nu_max
690 : end if
691 : end do
692 :
693 1536 : if (hybrid%masterthread)then
694 2 : write(iulog,*) "z computed from barometric formula (using US std atmosphere)"
695 2 : call std_atm_height(pmid(:),z(:))
696 2 : write(iulog,*) "k,pmid_ref,z,nu_lev,nu_t_lev,nu_div_lev"
697 188 : do k=1,nlev
698 188 : write(iulog,'(i3,5e11.4)') k,pmid(k),z(k),nu_lev(k),nu_t_lev(k),nu_div_lev(k)
699 : end do
700 2 : if (nu_top>0) then
701 2 : write(iulog,*) ": ksponge_end = ",ksponge_end
702 2 : write(iulog,*) ": sponge layer Laplacian damping"
703 2 : write(iulog,*) "k, p, z, nu_scale_top, nu (actual Laplacian damping coefficient)"
704 :
705 10 : do k=1,ksponge_end
706 10 : write(iulog,'(i3,4e11.4)') k,pmid(k),z(k),nu_scale_top(k),nu_scale_top(k)*nu_top
707 : end do
708 : end if
709 : end if
710 :
711 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
712 : !
713 : ! time-step information
714 : !
715 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
716 : !
717 : ! S=time-step stability region (i.e. advection w/leapfrog: S=1, viscosity w/forward Euler: S=2)
718 : !
719 1536 : if (tstep_type==1) then
720 0 : S_rk = 2.0_r8
721 0 : rk_str = ' * RK2-SSP 3 stage (same as tracers)'
722 1536 : elseif (tstep_type==2) then
723 0 : S_rk = 2.0_r8
724 0 : rk_str = ' * classic RK3'
725 1536 : elseif (tstep_type==3) then
726 0 : S_rk = 2.0_r8
727 0 : rk_str = ' * Kinnmark&Gray RK4'
728 1536 : elseif (tstep_type==4) then
729 1536 : S_rk = 3.0_r8
730 1536 : rk_str = ' * Kinnmark&Gray RK3 5 stage (3rd order)'
731 : end if
732 1536 : if (hybrid%masterthread) then
733 2 : write(iulog,'(a,f12.8,a)') 'Model top is ',ptop,'Pa'
734 2 : write(iulog,'(a)') ' '
735 2 : write(iulog,'(a)') 'Timestepping methods used in dynamical core:'
736 2 : write(iulog,'(a)')
737 2 : write(iulog,*) rk_str
738 2 : write(iulog,'(a)') ' * Spectral-element advection uses SSP preservation RK3'
739 2 : write(iulog,'(a)') ' * Viscosity operators use forward Euler'
740 : end if
741 1536 : S_laplacian = 2.0_r8 !using forward Euler for sponge diffusion
742 1536 : S_hypervis = 2.0_r8 !using forward Euler for hyperviscosity
743 1536 : S_rk_tracer = 2.0_r8
744 :
745 1536 : ugw = 342.0_r8 !max gravity wave speed
746 :
747 1536 : dt_max_adv = S_rk/(umax*max_normDinv*lambda_max*ra)
748 1536 : dt_max_gw = S_rk/(ugw*max_normDinv*lambda_max*ra)
749 1536 : dt_max_tracer_se = S_rk_tracer*min_gw/(umax*max_normDinv*ra)
750 1536 : if (use_cslam) then
751 1536 : if (large_Courant_incr) then
752 1536 : dt_max_tracer_fvm = dble(nhe)*(4.0_r8*pi*Rearth/dble(4.0_r8*ne*nc))/umax
753 : else
754 0 : dt_max_tracer_fvm = dble(nhe)*(2.0_r8*pi*Rearth/dble(4.0_r8*ne*nc))/umax
755 : end if
756 : else
757 0 : dt_max_tracer_fvm = -1.0_r8
758 : end if
759 434688 : nu_max = MAX(MAXVAL(nu_div_lev(:)),MAXVAL(nu_lev(:)),MAXVAL(nu_t_lev(:)))
760 1536 : dt_max_hypervis = s_hypervis/(nu_max*normDinv_hypervis)
761 1536 : dt_max_hypervis_tracer = s_hypervis/(nu_q*normDinv_hypervis)
762 :
763 290304 : max_laplace = MAX(MAXVAL(nu_scale_top(:))*nu_top,MAXVAL(kmvis_ref(:)/rho_ref(:)))
764 145920 : max_laplace = MAX(max_laplace,MAXVAL(kmcnd_ref(:)/(cpair*rho_ref(:))))
765 1536 : dt_max_laplacian_top = 1.0_r8/(max_laplace*((ra*max_normDinv)**2)*lambda_vis)
766 :
767 1536 : if (hybrid%masterthread) then
768 2 : write(iulog,'(a,f10.2,a)') ' '
769 2 : write(iulog,'(a,f10.2,a)') 'Estimates for maximum stable and actual time-steps for different aspects of algorithm:'
770 2 : write(iulog,'(a,f12.8,a)') '(assume max wind is ',umax,'m/s)'
771 2 : write(iulog,'(a)') '(assume max gravity wave speed is 342m/s)'
772 2 : write(iulog,'(a,f10.2,a)') ' '
773 2 : write(iulog,'(a,f10.2,a,f10.2,a)') '* dt_dyn (time-stepping dycore ; u,v,T,dM) < ',&
774 4 : MIN(dt_max_adv,dt_max_gw),'s ',dt_dyn_actual,'s'
775 2 : if (dt_dyn_actual>MIN(dt_max_adv,dt_max_gw)) write(iulog,*) 'WARNING: dt_dyn theoretically unstable'
776 :
777 2 : write(iulog,'(a,f10.2,a,f10.2,a)') '* dt_dyn_vis (hyperviscosity) ; u,v,T,dM) < ',dt_max_hypervis,&
778 4 : 's ',dt_dyn_visco_actual,'s'
779 2 : if (dt_dyn_visco_actual>dt_max_hypervis) write(iulog,*) 'WARNING: dt_dyn_vis theoretically unstable'
780 2 : if (.not.use_cslam) then
781 0 : write(iulog,'(a,f10.2,a,f10.2,a)') '* dt_tracer_se (time-stepping tracers ; q ) < ',dt_max_tracer_se,'s ',&
782 0 : dt_tracer_se_actual,'s'
783 0 : if (dt_tracer_se_actual>dt_max_tracer_se) write(iulog,*) 'WARNING: dt_tracer_se theoretically unstable'
784 0 : write(iulog,'(a,f10.2,a,f10.2,a)') '* dt_tracer_vis (hyperviscosity tracers; q ) < ',dt_max_hypervis_tracer,'s',&
785 0 : dt_tracer_visco_actual,'s'
786 0 : if (dt_tracer_visco_actual>dt_max_hypervis_tracer) write(iulog,*) 'WARNING: dt_tracer_hypervis theoretically unstable'
787 : end if
788 2 : if (use_cslam) then
789 2 : write(iulog,'(a,f10.2,a,f10.2,a)') '* dt_tracer_fvm (time-stepping tracers ; q ) < ',dt_max_tracer_fvm,&
790 4 : 's ',dt_tracer_fvm_actual
791 2 : if (dt_tracer_fvm_actual>dt_max_tracer_fvm) write(iulog,*) 'WARNING: dt_tracer_fvm theortically unstable'
792 : end if
793 2 : write(iulog,'(a,f10.2)') '* dt_remap (vertical remap dt) ',dt_remap_actual
794 10 : do k=1,ksponge_end
795 8 : max_laplace = MAX(nu_scale_top(k)*nu_top,kmvis_ref(k)/rho_ref(k))
796 8 : max_laplace = MAX(max_laplace,kmcnd_ref(k)/(cpair*rho_ref(k)))
797 8 : dt_max_laplacian_top = 1.0_r8/(max_laplace*((ra*max_normDinv)**2)*lambda_vis)
798 :
799 8 : write(iulog,'(a,f10.2,a,f10.2,a)') '* dt (del2 sponge ; u,v,T,dM) < ',&
800 16 : dt_max_laplacian_top,'s',dt_dyn_del2_actual,'s'
801 10 : if (dt_dyn_del2_actual>dt_max_laplacian_top) then
802 0 : if (k==1) then
803 0 : write(iulog,*) 'WARNING: theoretically unstable in sponge; increase se_hypervis_subcycle_sponge',&
804 0 : ' (this WARNING can sometimes be ignored in level 1)'
805 : else
806 0 : write(iulog,*) 'WARNING: theoretically unstable in sponge; increase se_hypervis_subcycle_sponge'
807 : endif
808 : end if
809 : end do
810 2 : write(iulog,*) ' '
811 2 : if (hypervis_power /= 0) then
812 0 : write(iulog,'(a,3e11.4)')'Scalar hyperviscosity (dynamics): ave,min,max = ', &
813 0 : nu*(/avg_hypervis**2,min_hypervis**2,max_hypervis**2/)
814 : end if
815 2 : write(iulog,*) 'tstep_type = ',tstep_type
816 : end if
817 1536 : end subroutine print_cfl
818 :
819 : !
820 : ! ============================
821 : ! global_maximum:
822 : !
823 : ! Find global maximum on sphere
824 : !
825 : ! ================================
826 :
827 0 : function global_maximum(fld,hybrid,npts,nets,nete) result(Max_sphere)
828 :
829 1536 : use hybrid_mod, only : hybrid_t
830 : use reduction_mod, only : red_max, pmax_mt
831 :
832 : integer , intent(in) :: npts,nets,nete
833 : real (kind=r8), intent(in) :: fld(npts,npts,nets:nete)
834 : type (hybrid_t) , intent(in) :: hybrid
835 :
836 : real (kind=r8) :: Max_sphere
837 :
838 : ! Local variables
839 :
840 : real (kind=r8) :: redp(1)
841 :
842 0 : Max_sphere = MAXVAL(fld(:,:,nets:nete))
843 :
844 0 : redp(1) = Max_sphere
845 0 : call pmax_mt(red_max,redp,1,hybrid)
846 0 : Max_sphere = red_max%buf(1)
847 :
848 0 : end function global_maximum
849 :
850 : ! ==========================================================
851 : ! l1_snorm:
852 : !
853 : ! computes the l1 norm per Williamson et al, p. 218 eq(8)
854 : ! for a scalar quantity
855 : ! ===========================================================
856 :
857 0 : function l1_snorm(elem,fld,fld_exact,hybrid,npts,nets,nete) result(l1)
858 :
859 : use element_mod, only : element_t
860 : use hybrid_mod, only : hybrid_t
861 :
862 : type(element_t) , intent(in) :: elem(:)
863 : integer , intent(in) :: npts,nets,nete
864 : real (kind=r8), intent(in) :: fld(npts,npts,nets:nete) ! computed soln
865 : real (kind=r8), intent(in) :: fld_exact(npts,npts,nets:nete) ! true soln
866 : type (hybrid_t) , intent(in) :: hybrid
867 : real (kind=r8) :: l1
868 :
869 : ! Local variables
870 :
871 0 : real (kind=r8) :: dfld_abs(npts,npts,nets:nete)
872 0 : real (kind=r8) :: fld_exact_abs(npts,npts,nets:nete)
873 : real (kind=r8) :: dfld_abs_int
874 : real (kind=r8) :: fld_exact_abs_int
875 : integer i,j,ie
876 :
877 0 : do ie=nets,nete
878 0 : do j=1,npts
879 0 : do i=1,npts
880 0 : dfld_abs(i,j,ie) = ABS(fld(i,j,ie)-fld_exact(i,j,ie))
881 0 : fld_exact_abs(i,j,ie) = ABS(fld_exact(i,j,ie))
882 : end do
883 : end do
884 : end do
885 :
886 0 : dfld_abs_int = global_integral(elem, dfld_abs(:,:,nets:nete),hybrid,npts,nets,nete)
887 0 : fld_exact_abs_int = global_integral(elem, fld_exact_abs(:,:,nets:nete),hybrid,npts,nets,nete)
888 :
889 0 : l1 = dfld_abs_int/fld_exact_abs_int
890 :
891 0 : end function l1_snorm
892 :
893 : ! ===========================================================
894 : ! l1_vnorm:
895 : !
896 : ! computes the l1 norm per Williamson et al, p. 218 eq(97),
897 : ! for a contravariant vector quantity on the velocity grid.
898 : !
899 : ! ===========================================================
900 :
901 0 : function l1_vnorm(elem, v,vt,hybrid,npts,nets,nete) result(l1)
902 : use element_mod, only : element_t
903 : use hybrid_mod, only : hybrid_t
904 :
905 : type(element_t) , intent(in), target :: elem(:)
906 : integer , intent(in) :: npts,nets,nete
907 : real (kind=r8), intent(in) :: v(npts,npts,2,nets:nete) ! computed soln
908 : real (kind=r8), intent(in) :: vt(npts,npts,2,nets:nete) ! true soln
909 : type (hybrid_t) , intent(in) :: hybrid
910 : real (kind=r8) :: l1
911 :
912 : ! Local variables
913 :
914 0 : real (kind=r8), dimension(:,:,:,:), pointer :: met
915 0 : real (kind=r8) :: dvsq(npts,npts,nets:nete)
916 0 : real (kind=r8) :: vtsq(npts,npts,nets:nete)
917 0 : real (kind=r8) :: dvco(npts,npts,2) ! covariant velocity
918 0 : real (kind=r8) :: vtco(npts,npts,2) ! covariant velocity
919 : real (kind=r8) :: dv1,dv2
920 : real (kind=r8) :: vt1,vt2
921 : real (kind=r8) :: dvsq_int
922 : real (kind=r8) :: vtsq_int
923 :
924 : integer i,j,ie
925 :
926 0 : do ie=nets,nete
927 0 : met => elem(ie)%met
928 0 : do j=1,npts
929 0 : do i=1,npts
930 :
931 0 : dv1 = v(i,j,1,ie)-vt(i,j,1,ie)
932 0 : dv2 = v(i,j,2,ie)-vt(i,j,2,ie)
933 :
934 0 : vt1 = vt(i,j,1,ie)
935 0 : vt2 = vt(i,j,2,ie)
936 :
937 0 : dvco(i,j,1) = met(i,j,1,1)*dv1 + met(i,j,1,2)*dv2
938 0 : dvco(i,j,2) = met(i,j,2,1)*dv1 + met(i,j,2,2)*dv2
939 :
940 0 : vtco(i,j,1) = met(i,j,1,1)*vt1 + met(i,j,1,2)*vt2
941 0 : vtco(i,j,2) = met(i,j,2,1)*vt1 + met(i,j,2,2)*vt2
942 :
943 0 : dvsq(i,j,ie) = SQRT(dvco(i,j,1)*dv1 + dvco(i,j,2)*dv2)
944 0 : vtsq(i,j,ie) = SQRT(vtco(i,j,1)*vt1 + vtco(i,j,2)*vt2)
945 :
946 : end do
947 : end do
948 : end do
949 :
950 0 : dvsq_int = global_integral(elem, dvsq(:,:,nets:nete),hybrid,npts,nets,nete)
951 0 : vtsq_int = global_integral(elem, vtsq(:,:,nets:nete),hybrid,npts,nets,nete)
952 :
953 0 : l1 = dvsq_int/vtsq_int
954 :
955 0 : end function l1_vnorm
956 :
957 : ! ==========================================================
958 : ! l2_snorm:
959 : !
960 : ! computes the l2 norm per Williamson et al, p. 218 eq(83)
961 : ! for a scalar quantity on the pressure grid.
962 : !
963 : ! ===========================================================
964 :
965 0 : function l2_snorm(elem,fld,fld_exact,hybrid,npts,nets,nete) result(l2)
966 : use element_mod, only : element_t
967 : use hybrid_mod, only : hybrid_t
968 :
969 : type(element_t), intent(in) :: elem(:)
970 : integer , intent(in) :: npts,nets,nete
971 : real (kind=r8), intent(in) :: fld(npts,npts,nets:nete) ! computed soln
972 : real (kind=r8), intent(in) :: fld_exact(npts,npts,nets:nete) ! true soln
973 : type (hybrid_t) , intent(in) :: hybrid
974 : real (kind=r8) :: l2
975 :
976 : ! Local variables
977 :
978 0 : real (kind=r8) :: dh2(npts,npts,nets:nete)
979 0 : real (kind=r8) :: fld_exact2(npts,npts,nets:nete)
980 : real (kind=r8) :: dh2_int
981 : real (kind=r8) :: fld_exact2_int
982 : integer i,j,ie
983 :
984 0 : do ie=nets,nete
985 0 : do j=1,npts
986 0 : do i=1,npts
987 0 : dh2(i,j,ie)=(fld(i,j,ie)-fld_exact(i,j,ie))**2
988 0 : fld_exact2(i,j,ie)=fld_exact(i,j,ie)**2
989 : end do
990 : end do
991 : end do
992 :
993 0 : dh2_int = global_integral(elem,dh2(:,:,nets:nete),hybrid,npts,nets,nete)
994 0 : fld_exact2_int = global_integral(elem,fld_exact2(:,:,nets:nete),hybrid,npts,nets,nete)
995 :
996 0 : l2 = SQRT(dh2_int)/SQRT(fld_exact2_int)
997 :
998 0 : end function l2_snorm
999 :
1000 : ! ==========================================================
1001 : ! l2_vnorm:
1002 : !
1003 : ! computes the l2 norm per Williamson et al, p. 219 eq(98)
1004 : ! for a contravariant vector quantity on the velocity grid.
1005 : !
1006 : ! ===========================================================
1007 :
1008 0 : function l2_vnorm(elem, v,vt,hybrid,npts,nets,nete) result(l2)
1009 : use element_mod, only : element_t
1010 : use hybrid_mod, only : hybrid_t
1011 :
1012 : type(element_t) , intent(in), target :: elem(:)
1013 : integer , intent(in) :: npts,nets,nete
1014 : real (kind=r8), intent(in) :: v(npts,npts,2,nets:nete) ! computed soln
1015 : real (kind=r8), intent(in) :: vt(npts,npts,2,nets:nete) ! true soln
1016 : type (hybrid_t) , intent(in) :: hybrid
1017 : real (kind=r8) :: l2
1018 :
1019 : ! Local variables
1020 :
1021 0 : real (kind=r8), dimension(:,:,:,:), pointer :: met
1022 0 : real (kind=r8) :: dvsq(npts,npts,nets:nete)
1023 0 : real (kind=r8) :: vtsq(npts,npts,nets:nete)
1024 0 : real (kind=r8) :: dvco(npts,npts,2) ! covariant velocity
1025 0 : real (kind=r8) :: vtco(npts,npts,2) ! covariant velocity
1026 : real (kind=r8) :: dv1,dv2
1027 : real (kind=r8) :: vt1,vt2
1028 : real (kind=r8) :: dvsq_int
1029 : real (kind=r8) :: vtsq_int
1030 : integer i,j,ie
1031 :
1032 0 : do ie=nets,nete
1033 0 : met => elem(ie)%met
1034 0 : do j=1,npts
1035 0 : do i=1,npts
1036 :
1037 0 : dv1 = v(i,j,1,ie)-vt(i,j,1,ie)
1038 0 : dv2 = v(i,j,2,ie)-vt(i,j,2,ie)
1039 :
1040 0 : vt1 = vt(i,j,1,ie)
1041 0 : vt2 = vt(i,j,2,ie)
1042 :
1043 0 : dvco(i,j,1) = met(i,j,1,1)*dv1 + met(i,j,1,2)*dv2
1044 0 : dvco(i,j,2) = met(i,j,2,1)*dv1 + met(i,j,2,2)*dv2
1045 :
1046 0 : vtco(i,j,1) = met(i,j,1,1)*vt1 + met(i,j,1,2)*vt2
1047 0 : vtco(i,j,2) = met(i,j,2,1)*vt1 + met(i,j,2,2)*vt2
1048 :
1049 0 : dvsq(i,j,ie) = dvco(i,j,1)*dv1 + dvco(i,j,2)*dv2
1050 0 : vtsq(i,j,ie) = vtco(i,j,1)*vt1 + vtco(i,j,2)*vt2
1051 :
1052 : end do
1053 : end do
1054 : end do
1055 :
1056 0 : dvsq_int = global_integral(elem, dvsq(:,:,nets:nete),hybrid,npts,nets,nete)
1057 0 : vtsq_int = global_integral(elem, vtsq(:,:,nets:nete),hybrid,npts,nets,nete)
1058 :
1059 0 : l2 = SQRT(dvsq_int)/SQRT(vtsq_int)
1060 :
1061 0 : end function l2_vnorm
1062 :
1063 : ! ==========================================================
1064 : ! linf_snorm:
1065 : !
1066 : ! computes the l infinity norm per Williamson et al, p. 218 eq(84)
1067 : ! for a scalar quantity on the pressure grid...
1068 : !
1069 : ! ===========================================================
1070 :
1071 0 : function linf_snorm(fld,fld_exact,hybrid,npts,nets,nete) result(linf)
1072 : use hybrid_mod, only : hybrid_t
1073 : integer , intent(in) :: npts,nets,nete
1074 : real (kind=r8), intent(in) :: fld(npts,npts,nets:nete) ! computed soln
1075 : real (kind=r8), intent(in) :: fld_exact(npts,npts,nets:nete) ! true soln
1076 : type (hybrid_t) , intent(in) :: hybrid
1077 : real (kind=r8) :: linf
1078 :
1079 : ! Local variables
1080 :
1081 0 : real (kind=r8) :: dfld_abs(npts,npts,nets:nete)
1082 0 : real (kind=r8) :: fld_exact_abs(npts,npts,nets:nete)
1083 : real (kind=r8) :: dfld_abs_max
1084 : real (kind=r8) :: fld_exact_abs_max
1085 : integer i,j,ie
1086 :
1087 0 : do ie=nets,nete
1088 0 : do j=1,npts
1089 0 : do i=1,npts
1090 0 : dfld_abs(i,j,ie)=ABS(fld(i,j,ie)-fld_exact(i,j,ie))
1091 0 : fld_exact_abs(i,j,ie)=ABS(fld_exact(i,j,ie))
1092 : end do
1093 : end do
1094 : end do
1095 :
1096 0 : dfld_abs_max = global_maximum(dfld_abs(:,:,nets:nete),hybrid,npts,nets,nete)
1097 0 : fld_exact_abs_max = global_maximum(fld_exact_abs(:,:,nets:nete),hybrid,npts,nets,nete)
1098 :
1099 0 : linf = dfld_abs_max/fld_exact_abs_max
1100 :
1101 0 : end function linf_snorm
1102 :
1103 :
1104 : ! ==========================================================
1105 : ! linf_vnorm:
1106 : !
1107 : ! computes the linf norm per Williamson et al, p. 218 eq(99),
1108 : ! for a contravariant vector quantity on the velocity grid.
1109 : !
1110 : ! ===========================================================
1111 :
1112 0 : function linf_vnorm(elem,v,vt,hybrid,npts,nets,nete) result(linf)
1113 : use hybrid_mod, only : hybrid_t
1114 : use element_mod, only : element_t
1115 :
1116 : type(element_t) , intent(in), target :: elem(:)
1117 : integer , intent(in) :: npts,nets,nete
1118 : real (kind=r8), intent(in) :: v(npts,npts,2,nets:nete) ! computed soln
1119 : real (kind=r8), intent(in) :: vt(npts,npts,2,nets:nete) ! true soln
1120 : type (hybrid_t) , intent(in) :: hybrid
1121 : real (kind=r8) :: linf
1122 :
1123 : ! Local variables
1124 :
1125 0 : real (kind=r8), dimension(:,:,:,:), pointer :: met
1126 0 : real (kind=r8) :: dvsq(npts,npts,nets:nete)
1127 0 : real (kind=r8) :: vtsq(npts,npts,nets:nete)
1128 0 : real (kind=r8) :: dvco(npts,npts,2) ! covariant velocity
1129 0 : real (kind=r8) :: vtco(npts,npts,2) ! covariant velocity
1130 : real (kind=r8) :: dv1,dv2
1131 : real (kind=r8) :: vt1,vt2
1132 : real (kind=r8) :: dvsq_max
1133 : real (kind=r8) :: vtsq_max
1134 : integer i,j,ie
1135 :
1136 0 : do ie=nets,nete
1137 0 : met => elem(ie)%met
1138 :
1139 0 : do j=1,npts
1140 0 : do i=1,npts
1141 :
1142 0 : dv1 = v(i,j,1,ie)-vt(i,j,1,ie)
1143 0 : dv2 = v(i,j,2,ie)-vt(i,j,2,ie)
1144 :
1145 0 : vt1 = vt(i,j,1,ie)
1146 0 : vt2 = vt(i,j,2,ie)
1147 :
1148 0 : dvco(i,j,1) = met(i,j,1,1)*dv1 + met(i,j,1,2)*dv2
1149 0 : dvco(i,j,2) = met(i,j,2,1)*dv1 + met(i,j,2,2)*dv2
1150 :
1151 0 : vtco(i,j,1) = met(i,j,1,1)*vt1 + met(i,j,1,2)*vt2
1152 0 : vtco(i,j,2) = met(i,j,2,1)*vt1 + met(i,j,2,2)*vt2
1153 :
1154 0 : dvsq(i,j,ie) = SQRT(dvco(i,j,1)*dv1 + dvco(i,j,2)*dv2)
1155 0 : vtsq(i,j,ie) = SQRT(vtco(i,j,1)*vt1 + vtco(i,j,2)*vt2)
1156 :
1157 : end do
1158 : end do
1159 : end do
1160 :
1161 0 : dvsq_max = global_maximum(dvsq(:,:,nets:nete),hybrid,npts,nets,nete)
1162 0 : vtsq_max = global_maximum(vtsq(:,:,nets:nete),hybrid,npts,nets,nete)
1163 :
1164 0 : linf = dvsq_max/vtsq_max
1165 :
1166 0 : end function linf_vnorm
1167 :
1168 6912 : subroutine wrap_repro_sum (nvars, comm, nsize)
1169 : use dimensions_mod, only: nelemd
1170 : use shr_reprosum_mod, only: repro_sum => shr_reprosum_calc
1171 : use cam_abortutils, only: endrun
1172 : use parallel_mod, only: global_shared_buf, global_shared_sum, nrepro_vars
1173 :
1174 : integer :: nvars ! number of variables to be summed (cannot exceed nrepro_vars)
1175 : integer :: comm ! mpi communicator
1176 : integer, optional :: nsize ! local buffer size (defaults to nelemd - number of elements in mpi task)
1177 :
1178 : integer nsize_use
1179 :
1180 6912 : if (present(nsize)) then
1181 0 : nsize_use = nsize
1182 : else
1183 6912 : nsize_use = nelemd
1184 : endif
1185 6912 : if (nvars .gt. nrepro_vars) call endrun('ERROR: repro_sum_buffer_size exceeded')
1186 :
1187 : ! Repro_sum contains its own OpenMP, so only one thread should call it (AAM)
1188 :
1189 : !$OMP BARRIER
1190 : !$OMP MASTER
1191 :
1192 6912 : call repro_sum(global_shared_buf, global_shared_sum, nsize_use, nelemd, nvars, commid=comm)
1193 :
1194 :
1195 : !$OMP END MASTER
1196 : !$OMP BARRIER
1197 :
1198 6912 : end subroutine wrap_repro_sum
1199 :
1200 4608 : subroutine automatically_set_viscosity_coefficients(hybrid,ne,max_min_dx,min_min_dx,nu,factor,str)
1201 6912 : use physconst, only: rearth
1202 : use control_mod, only: hypervis_scaling,hypervis_power
1203 : use hybrid_mod, only: hybrid_t
1204 : use cam_abortutils, only: endrun
1205 :
1206 : type (hybrid_t), intent(in) :: hybrid
1207 : integer , intent(in) :: ne
1208 : real (kind=r8), intent(in) :: max_min_dx,min_min_dx,factor
1209 : real (kind=r8), intent(inout) :: nu
1210 : character(len=4), intent(in) :: str
1211 :
1212 : real(r8) :: uniform_res_hypervis_scaling,nu_fac
1213 : real(kind=r8) :: nu_min, nu_max
1214 : !
1215 : !************************************************************************************************************
1216 : !
1217 : ! automatically set viscosity coefficients
1218 : !
1219 : !
1220 : ! Use scaling from
1221 : !
1222 : ! - Boville, B. A., 1991: Sensitivity of simulated climate to
1223 : ! model resolution. J. Climate, 4, 469-485.
1224 : !
1225 : ! - TAKAHASHI ET AL., 2006: GLOBAL SIMULATION OF MESOSCALE SPECTRUM
1226 : !
1227 4608 : uniform_res_hypervis_scaling = 1.0_r8/log10(2.0_r8)
1228 : !
1229 : ! compute factor so that at ne30 resolution nu=1E15
1230 : ! scale so that scaling works for other planets
1231 : !
1232 : ! grid spacing in meters = max_min_dx*1000.0_r8
1233 : !
1234 4608 : nu_fac = (rearth/6.37122E6_r8)*1.0E15_r8/(110000.0_r8**uniform_res_hypervis_scaling)
1235 :
1236 4608 : if (nu < 0) then
1237 4608 : if (ne <= 0) then
1238 0 : if (hypervis_power/=0) then
1239 0 : call endrun('ERROR: Automatic scaling of scalar viscosity not implemented')
1240 0 : else if (hypervis_scaling/=0) then
1241 0 : nu_min = factor*nu_fac*(max_min_dx*1000.0_r8)**uniform_res_hypervis_scaling
1242 0 : nu_max = factor*nu_fac*(min_min_dx*1000.0_r8)**uniform_res_hypervis_scaling
1243 0 : nu = factor*nu_min
1244 0 : if (hybrid%masterthread) then
1245 0 : write(iulog,'(a,a)') "Automatically setting nu",TRIM(str)
1246 0 : write(iulog,'(a,2e9.2,a,2f9.2)') "Value at min/max grid spacing: ",nu_min,nu_max,&
1247 0 : " Max/min grid spacing (km) = ",max_min_dx,min_min_dx
1248 : end if
1249 0 : nu = nu_min*(2.0_r8*rearth/(3.0_r8*max_min_dx*1000.0_r8))**hypervis_scaling/(rearth**4)
1250 0 : if (hybrid%masterthread) &
1251 0 : write(iulog,'(a,a,a,e9.3)') "Nu_tensor",TRIM(str)," = ",nu
1252 : end if
1253 : else
1254 4608 : nu = factor*nu_fac*((30.0_r8/ne)*110000.0_r8)**uniform_res_hypervis_scaling
1255 4608 : if (hybrid%masterthread) then
1256 6 : write(iulog,'(a,a,a,e9.2)') "Automatically setting nu",TRIM(str)," =",nu
1257 : end if
1258 : end if
1259 : end if
1260 4608 : end subroutine automatically_set_viscosity_coefficients
1261 : end module global_norms_mod
|