Line data Source code
1 : module viscosity_mod
2 : !
3 : ! This module should be renamed "global_deriv_mod.F90"
4 : !
5 : ! It is a collection of derivative operators that must be applied to the field
6 : ! over the sphere (as opposed to derivative operators that can be applied element
7 : ! by element)
8 : !
9 : !
10 : use shr_kind_mod, only: r8=>shr_kind_r8
11 : use thread_mod, only: max_num_threads, omp_get_num_threads
12 : use dimensions_mod, only: np, nc, nlev,nlevp, qsize,nelemd
13 : use hybrid_mod, only: hybrid_t, get_loop_ranges, config_thread_region
14 : use parallel_mod, only: parallel_t
15 : use element_mod, only: element_t
16 : use derivative_mod, only: derivative_t, laplace_sphere_wk, vlaplace_sphere_wk, vorticity_sphere, derivinit, divergence_sphere
17 : use edgetype_mod, only: EdgeBuffer_t, EdgeDescriptor_t
18 : use edge_mod, only: edgevpack, edgevunpack, edgeVunpackmin, edgeSunpackmin, &
19 : edgeVunpackmax, initEdgeBuffer, FreeEdgeBuffer, edgeSunpackmax, edgeSpack
20 : use bndry_mod, only: bndry_exchange, bndry_exchange_start,bndry_exchange_finish
21 : use control_mod, only: hypervis_scaling, nu, nu_div
22 : use thread_mod, only: vert_num_threads
23 :
24 : implicit none
25 : save
26 :
27 : public :: biharmonic_wk_scalar
28 : public :: biharmonic_wk_omega
29 : public :: neighbor_minmax, neighbor_minmax_start,neighbor_minmax_finish
30 :
31 : !
32 : ! compute vorticity/divergence and then project to make continious
33 : ! high-level routines uses only for I/O
34 : public :: compute_zeta_C0
35 : public :: compute_div_C0
36 :
37 : interface compute_zeta_C0
38 : module procedure compute_zeta_C0_hybrid ! hybrid version
39 : module procedure compute_zeta_C0_par ! single threaded
40 : end interface compute_zeta_C0
41 : interface compute_div_C0
42 : module procedure compute_div_C0_hybrid
43 : module procedure compute_div_C0_par
44 : end interface compute_div_C0
45 :
46 : public :: compute_zeta_C0_contra ! for older versions of sweq which carry
47 : public :: compute_div_C0_contra ! velocity around in contra-coordinates
48 :
49 : type (EdgeBuffer_t) :: edge1
50 :
51 : CONTAINS
52 :
53 262656 : subroutine biharmonic_wk_dp3d(elem,dptens,dpflux,ttens,vtens,deriv,edge3,hybrid,nt,nets,nete,kbeg,kend,hvcoord)
54 : use derivative_mod, only : subcell_Laplace_fluxes
55 : use dimensions_mod, only : use_cslam, nu_div_lev,nu_lev
56 : use hybvcoord_mod, only : hvcoord_t
57 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
58 : ! compute weak biharmonic operator
59 : ! input: h,v (stored in elem()%, in lat-lon coordinates
60 : ! output: ttens,vtens overwritten with weak biharmonic of h,v (output in lat-lon coordinates)
61 : !
62 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
63 : type (hybrid_t) , intent(in) :: hybrid
64 : type (element_t) , intent(inout), target :: elem(:)
65 : integer , intent(in) :: nt,nets,nete
66 : integer , intent(in) :: kbeg, kend
67 : real (kind=r8), intent(out), dimension(nc,nc,4,nlev,nets:nete) :: dpflux
68 : real (kind=r8), dimension(np,np,2,nlev,nets:nete) :: vtens
69 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: ttens,dptens
70 : type (EdgeBuffer_t) , intent(inout) :: edge3
71 : type (derivative_t) , intent(in) :: deriv
72 : type (hvcoord_t) , intent(in) :: hvcoord
73 : ! local
74 : integer :: i,j,k,kptr,ie,kblk
75 : ! real (kind=r8), dimension(:,:), pointer :: rspheremv
76 : real (kind=r8), dimension(np,np) :: tmp
77 : real (kind=r8), dimension(np,np) :: tmp2
78 : real (kind=r8), dimension(np,np,2) :: v
79 :
80 : real (kind=r8), dimension(np,np,nlev) :: lap_p_wk
81 : real (kind=r8), dimension(np,np,nlevp) :: T_i
82 :
83 :
84 : real (kind=r8) :: nu_ratio1, nu_ratio2, dp_thresh
85 : logical var_coef1
86 :
87 262656 : kblk = kend - kbeg + 1
88 :
89 9104986656 : if (use_cslam) dpflux = 0
90 : !if tensor hyperviscosity with tensor V is used, then biharmonic operator is (\grad\cdot V\grad) (\grad \cdot \grad)
91 : !so tensor is only used on second call to laplace_sphere_wk
92 262656 : var_coef1 = .true.
93 262656 : if(hypervis_scaling > 0) var_coef1 = .false.
94 262656 : dp_thresh=.025_r8 ! tunable coefficient
95 2109456 : do ie=nets,nete
96 : !$omp parallel do num_threads(vert_num_threads) private(k,tmp)
97 173599200 : do k=kbeg,kend
98 171752400 : nu_ratio1=1
99 171752400 : nu_ratio2=1
100 171752400 : if (nu_div_lev(k)/=nu_lev(k)) then
101 171752400 : if(hypervis_scaling /= 0) then
102 : ! we have a problem with the tensor in that we cant seperate
103 : ! div and curl components. So we do, with tensor V:
104 : ! nu * (del V del ) * ( nu_ratio * grad(div) - curl(curl))
105 0 : nu_ratio1=nu_div_lev(k)/nu_lev(k)
106 : nu_ratio2=1
107 : else
108 171752400 : nu_ratio1=sqrt(nu_div_lev(k)/nu_lev(k))
109 171752400 : nu_ratio2=sqrt(nu_div_lev(k)/nu_lev(k))
110 : endif
111 : endif
112 :
113 3606800400 : tmp=elem(ie)%state%T(:,:,k,nt)-elem(ie)%derived%T_ref(:,:,k)
114 171752400 : call laplace_sphere_wk(tmp,deriv,elem(ie),ttens(:,:,k,ie),var_coef=var_coef1)
115 :
116 3606800400 : tmp=elem(ie)%state%dp3d(:,:,k,nt)-elem(ie)%derived%dp_ref(:,:,k)
117 171752400 : call laplace_sphere_wk(tmp,deriv,elem(ie),dptens(:,:,k,ie),var_coef=var_coef1)
118 :
119 171752400 : call vlaplace_sphere_wk(elem(ie)%state%v(:,:,:,k,nt),deriv,elem(ie),.true.,vtens(:,:,:,k,ie), &
120 345351600 : var_coef=var_coef1,nu_ratio=nu_ratio1)
121 : enddo
122 :
123 1846800 : kptr = kbeg - 1
124 1846800 : call edgeVpack(edge3,ttens(:,:,kbeg:kend,ie),kblk,kptr,ie)
125 :
126 1846800 : kptr = kbeg - 1 + nlev
127 3608647200 : call edgeVpack(edge3,vtens(:,:,1,kbeg:kend,ie),kblk,kptr,ie)
128 :
129 1846800 : kptr = kbeg - 1 + 2*nlev
130 3608647200 : call edgeVpack(edge3,vtens(:,:,2,kbeg:kend,ie),kblk,kptr,ie)
131 :
132 1846800 : kptr = kbeg - 1 + 3*nlev
133 2109456 : call edgeVpack(edge3,dptens(:,:,kbeg:kend,ie),kblk,kptr,ie)
134 : enddo
135 :
136 262656 : call bndry_exchange(hybrid,edge3,location='biharmonic_wk_dp3d')
137 :
138 2109456 : do ie=nets,nete
139 : !CLEAN rspheremv => elem(ie)%rspheremp(:,:)
140 :
141 1846800 : kptr = kbeg - 1
142 1846800 : call edgeVunpack(edge3,ttens(:,:,kbeg:kend,ie),kblk,kptr,ie)
143 :
144 1846800 : kptr = kbeg - 1 + nlev
145 7215447600 : call edgeVunpack(edge3,vtens(:,:,1,kbeg:kend,ie),kblk,kptr,ie)
146 :
147 1846800 : kptr = kbeg - 1 + 2*nlev
148 7215447600 : call edgeVunpack(edge3,vtens(:,:,2,kbeg:kend,ie),kblk,kptr,ie)
149 :
150 1846800 : kptr = kbeg - 1 + 3*nlev
151 1846800 : call edgeVunpack(edge3,dptens(:,:,kbeg:kend,ie),kblk,kptr,ie)
152 :
153 1846800 : if (use_cslam) then
154 173599200 : do k=1,nlev
155 : !CLEAN tmp(:,:)= rspheremv(:,:)*dptens(:,:,k,ie)
156 3606800400 : tmp(:,:)= elem(ie)%rspheremp(:,:)*dptens(:,:,k,ie)
157 173599200 : call subcell_Laplace_fluxes(tmp, deriv, elem(ie), np, nc,dpflux(:,:,:,k,ie))
158 : enddo
159 : endif
160 :
161 : ! apply inverse mass matrix, then apply laplace again
162 : !$omp parallel do num_threads(vert_num_threads) private(k,v,tmp,tmp2)
163 173861856 : do k=kbeg,kend
164 : !CLEAN tmp(:,:)=rspheremv(:,:)*ttens(:,:,k,ie)
165 3606800400 : tmp(:,:)=elem(ie)%rspheremp(:,:)*ttens(:,:,k,ie)
166 171752400 : call laplace_sphere_wk(tmp,deriv,elem(ie),ttens(:,:,k,ie),var_coef=.true.)
167 : !CLEAN tmp2(:,:)=rspheremv(:,:)*dptens(:,:,k,ie)
168 3606800400 : tmp2(:,:)=elem(ie)%rspheremp(:,:)*dptens(:,:,k,ie)
169 171752400 : call laplace_sphere_wk(tmp2,deriv,elem(ie),dptens(:,:,k,ie),var_coef=.true.)
170 : !CLEAN v(:,:,1)=rspheremv(:,:)*vtens(:,:,1,k,ie)
171 : !CLEAN v(:,:,2)=rspheremv(:,:)*vtens(:,:,2,k,ie)
172 :
173 3606800400 : v(:,:,1)=elem(ie)%rspheremp(:,:)*vtens(:,:,1,k,ie)
174 3606800400 : v(:,:,2)=elem(ie)%rspheremp(:,:)*vtens(:,:,2,k,ie)
175 : call vlaplace_sphere_wk(v(:,:,:),deriv,elem(ie),.true.,vtens(:,:,:,k,ie), &
176 173599200 : var_coef=.true.,nu_ratio=nu_ratio2)
177 :
178 : enddo
179 : enddo
180 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
181 262656 : end subroutine biharmonic_wk_dp3d
182 :
183 :
184 48384 : subroutine biharmonic_wk_omega(elem,ptens,deriv,edge3,hybrid,nets,nete,kbeg,kend)
185 : type (hybrid_t) , intent(in) :: hybrid
186 : type (element_t) , intent(inout), target :: elem(:)
187 : integer , intent(in) :: nets,nete
188 : integer , intent(in) :: kbeg, kend
189 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: ptens
190 : type (EdgeBuffer_t) , intent(inout) :: edge3
191 : type (derivative_t) , intent(in) :: deriv
192 :
193 : ! local
194 : integer :: i,j,k,kptr,ie,kblk
195 48384 : real (kind=r8), dimension(:,:), pointer :: rspheremv
196 : real (kind=r8), dimension(np,np) :: tmp
197 : real (kind=r8), dimension(np,np) :: tmp2
198 : real (kind=r8), dimension(np,np,2) :: v
199 : real (kind=r8) :: nu_ratio1, nu_ratio2
200 : logical var_coef1
201 :
202 48384 : kblk = kend - kbeg + 1
203 :
204 : !if tensor hyperviscosity with tensor V is used, then biharmonic operator is (\grad\cdot V\grad) (\grad \cdot \grad)
205 : !so tensor is only used on second call to laplace_sphere_wk
206 48384 : var_coef1 = .true.
207 0 : if(hypervis_scaling > 0) var_coef1 = .false.
208 :
209 48384 : nu_ratio1=1
210 48384 : nu_ratio2=1
211 :
212 388584 : do ie=nets,nete
213 :
214 : !$omp parallel do num_threads(vert_num_threads) private(k,tmp)
215 31978800 : do k=kbeg,kend
216 664410600 : tmp=elem(ie)%derived%omega(:,:,k)
217 31978800 : call laplace_sphere_wk(tmp,deriv,elem(ie),ptens(:,:,k,ie),var_coef=var_coef1)
218 : enddo
219 :
220 340200 : kptr = kbeg - 1
221 388584 : call edgeVpack(edge3,ptens(:,:,kbeg:kend,ie),kblk,kptr,ie)
222 : enddo
223 :
224 48384 : call bndry_exchange(hybrid,edge3,location='biharmonic_wk_omega')
225 :
226 388584 : do ie=nets,nete
227 340200 : rspheremv => elem(ie)%rspheremp(:,:)
228 :
229 340200 : kptr = kbeg - 1
230 340200 : call edgeVunpack(edge3,ptens(:,:,kbeg:kend,ie),kblk,kptr,ie)
231 :
232 : ! apply inverse mass matrix, then apply laplace again
233 : !$omp parallel do num_threads(vert_num_threads) private(k,tmp)
234 32027184 : do k=kbeg,kend
235 664410600 : tmp(:,:)=rspheremv(:,:)*ptens(:,:,k,ie)
236 31978800 : call laplace_sphere_wk(tmp,deriv,elem(ie),ptens(:,:,k,ie),var_coef=.true.)
237 : enddo
238 : enddo
239 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
240 311040 : end subroutine biharmonic_wk_omega
241 :
242 :
243 0 : subroutine biharmonic_wk_scalar(elem,qtens,deriv,edgeq,hybrid,nets,nete)
244 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
245 : ! compute weak biharmonic operator
246 : ! input: qtens = Q
247 : ! output: qtens = weak biharmonic of Q
248 : !
249 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
250 : type (hybrid_t) , intent(in) :: hybrid
251 : type (element_t) , intent(inout), target :: elem(:)
252 : integer :: nets,nete
253 : real (kind=r8), dimension(np,np,nlev,qsize,nets:nete) :: qtens
254 : type (EdgeBuffer_t) , intent(inout) :: edgeq
255 : type (derivative_t) , intent(in) :: deriv
256 :
257 : ! local
258 : integer :: k,kptr,i,j,ie,ic,q
259 : integer :: kbeg,kend,qbeg,qend
260 : real (kind=r8), dimension(np,np) :: lap_p
261 : logical var_coef1
262 : integer :: kblk,qblk ! The per thead size of the vertical and tracers
263 :
264 0 : call get_loop_ranges(hybrid,kbeg=kbeg,kend=kend,qbeg=qbeg,qend=qend)
265 :
266 : !if tensor hyperviscosity with tensor V is used, then biharmonic operator is (\grad\cdot V\grad) (\grad \cdot \grad)
267 : !so tensor is only used on second call to laplace_sphere_wk
268 0 : var_coef1 = .true.
269 0 : if(hypervis_scaling > 0) var_coef1 = .false.
270 :
271 :
272 0 : kblk = kend - kbeg + 1 ! calculate size of the block of vertical levels
273 0 : qblk = qend - qbeg + 1 ! calculate size of the block of tracers
274 :
275 0 : do ie=nets,nete
276 0 : do q=qbeg,qend
277 0 : do k=kbeg,kend
278 0 : lap_p(:,:)=qtens(:,:,k,q,ie)
279 0 : call laplace_sphere_wk(lap_p,deriv,elem(ie),qtens(:,:,k,q,ie),var_coef=var_coef1)
280 : enddo
281 0 : kptr = nlev*(q-1) + kbeg - 1
282 0 : call edgeVpack(edgeq, qtens(:,:,kbeg:kend,q,ie),kblk,kptr,ie)
283 : enddo
284 : enddo
285 :
286 :
287 0 : call bndry_exchange(hybrid,edgeq,location='biharmonic_wk_scalar')
288 :
289 0 : do ie=nets,nete
290 :
291 : ! apply inverse mass matrix, then apply laplace again
292 0 : do q=qbeg,qend
293 0 : kptr = nlev*(q-1) + kbeg - 1
294 0 : call edgeVunpack(edgeq, qtens(:,:,kbeg:kend,q,ie),kblk,kptr,ie)
295 0 : do k=kbeg,kend
296 0 : lap_p(:,:)=elem(ie)%rspheremp(:,:)*qtens(:,:,k,q,ie)
297 0 : call laplace_sphere_wk(lap_p,deriv,elem(ie),qtens(:,:,k,q,ie),var_coef=.true.)
298 : enddo
299 : enddo
300 : enddo
301 :
302 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
303 0 : end subroutine biharmonic_wk_scalar
304 :
305 :
306 0 : subroutine make_C0(zeta,elem,hybrid,nets,nete)
307 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
308 : ! apply DSS (aka assembly procedure) to zeta.
309 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
310 :
311 : type (hybrid_t) , intent(in) :: hybrid
312 : type (element_t) , intent(in), target :: elem(:)
313 : integer :: nets,nete
314 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: zeta
315 :
316 : ! local
317 : integer :: k,i,j,ie,ic,kptr,nthread_save
318 :
319 :
320 0 : call initEdgeBuffer(hybrid%par,edge1,elem,nlev)
321 :
322 0 : do ie=nets,nete
323 : #if (defined COLUMN_OPENMP)
324 : !$omp parallel do num_threads(vert_num_threads) private(k)
325 : #endif
326 0 : do k=1,nlev
327 0 : zeta(:,:,k,ie)=zeta(:,:,k,ie)*elem(ie)%spheremp(:,:)
328 : enddo
329 0 : kptr=0
330 0 : call edgeVpack(edge1, zeta(1,1,1,ie),nlev,kptr,ie)
331 : enddo
332 0 : call bndry_exchange(hybrid,edge1,location='make_C0')
333 0 : do ie=nets,nete
334 0 : kptr=0
335 0 : call edgeVunpack(edge1, zeta(1,1,1,ie),nlev,kptr, ie)
336 : #if (defined COLUMN_OPENMP)
337 : !$omp parallel do num_threads(vert_num_threads) private(k)
338 : #endif
339 0 : do k=1,nlev
340 0 : zeta(:,:,k,ie)=zeta(:,:,k,ie)*elem(ie)%rspheremp(:,:)
341 : enddo
342 : enddo
343 :
344 0 : call FreeEdgeBuffer(edge1)
345 :
346 0 : end subroutine
347 :
348 :
349 0 : subroutine make_C0_vector(v,elem,hybrid,nets,nete)
350 : #if 1
351 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
352 : ! apply DSS to a velocity vector
353 : ! this is a low-performance routine used for I/O and analysis.
354 : ! no need to optimize
355 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
356 : type (hybrid_t) , intent(in) :: hybrid
357 : type (element_t) , intent(in), target :: elem(:)
358 : integer :: nets,nete
359 : real (kind=r8), dimension(np,np,2,nlev,nets:nete) :: v
360 :
361 : ! local
362 : integer :: k,i,j,ie,ic,kptr
363 : type (EdgeBuffer_t) :: edge2
364 0 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: v1
365 :
366 0 : v1(:,:,:,:) = v(:,:,1,:,:)
367 0 : call make_C0(v1,elem,hybrid,nets,nete)
368 0 : v(:,:,1,:,:) = v1(:,:,:,:)
369 :
370 0 : v1(:,:,:,:) = v(:,:,2,:,:)
371 0 : call make_C0(v1,elem,hybrid,nets,nete)
372 0 : v(:,:,2,:,:) = v1(:,:,:,:)
373 : #else
374 : type (hybrid_t) , intent(in) :: hybrid
375 : type (element_t) , intent(in), target :: elem(:)
376 : integer :: nets,nete
377 : real (kind=r8), dimension(np,np,2,nlev,nets:nete) :: v
378 :
379 : ! local
380 : integer :: k,i,j,ie,ic,kptr
381 : type (EdgeBuffer_t) :: edge2
382 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: v1
383 :
384 :
385 :
386 : call initEdgeBuffer(hybrid%par,edge2,elem,2*nlev)
387 :
388 : do ie=nets,nete
389 : #if (defined COLUMN_OPENMP)
390 : !$omp parallel do num_threads(vert_num_threads) private(k)
391 : #endif
392 : do k=1,nlev
393 : v(:,:,1,k,ie)=v(:,:,1,k,ie)*elem(ie)%spheremp(:,:)
394 : v(:,:,2,k,ie)=v(:,:,2,k,ie)*elem(ie)%spheremp(:,:)
395 : enddo
396 : kptr=0
397 : call edgeVpack(edge2, v(1,1,1,1,ie),2*nlev,kptr,ie)
398 : enddo
399 : call bndry_exchange(hybrid,edge2,location='make_C0_vector')
400 : do ie=nets,nete
401 : kptr=0
402 : call edgeVunpack(edge2, v(1,1,1,1,ie),2*nlev,kptr,ie)
403 : #if (defined COLUMN_OPENMP)
404 : !$omp parallel do num_threads(vert_num_threads) private(k)
405 : #endif
406 : do k=1,nlev
407 : v(:,:,1,k,ie)=v(:,:,1,k,ie)*elem(ie)%rspheremp(:,:)
408 : v(:,:,2,k,ie)=v(:,:,2,k,ie)*elem(ie)%rspheremp(:,:)
409 : enddo
410 : enddo
411 :
412 : call FreeEdgeBuffer(edge2)
413 : #endif
414 0 : end subroutine
415 :
416 :
417 :
418 :
419 :
420 :
421 0 : subroutine compute_zeta_C0_contra(zeta,elem,hybrid,nets,nete,nt)
422 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
423 : ! compute C0 vorticity. That is, solve:
424 : ! < PHI, zeta > = <PHI, curl(elem%state%v >
425 : !
426 : ! input: v (stored in elem()%, in contra-variant coordinates)
427 : ! output: zeta(:,:,:,:)
428 : !
429 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
430 :
431 : type (hybrid_t) , intent(in) :: hybrid
432 : type (element_t) , intent(in), target :: elem(:)
433 : integer :: nt,nets,nete
434 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: zeta
435 : real (kind=r8), dimension(np,np,2) :: ulatlon
436 : real (kind=r8), dimension(np,np) :: v1,v2
437 :
438 : ! local
439 : integer :: k,ie
440 : type (derivative_t) :: deriv
441 :
442 0 : call derivinit(deriv)
443 :
444 0 : do k=1,nlev
445 0 : do ie=nets,nete
446 0 : v1 = elem(ie)%state%v(:,:,1,k,nt)
447 0 : v2 = elem(ie)%state%v(:,:,2,k,nt)
448 0 : ulatlon(:,:,1) = elem(ie)%D(:,:,1,1)*v1 + elem(ie)%D(:,:,1,2)*v2
449 0 : ulatlon(:,:,2) = elem(ie)%D(:,:,2,1)*v1 + elem(ie)%D(:,:,2,2)*v2
450 0 : call vorticity_sphere(ulatlon,deriv,elem(ie),zeta(:,:,k,ie))
451 : enddo
452 : enddo
453 :
454 0 : call make_C0(zeta,elem,hybrid,nets,nete)
455 :
456 0 : end subroutine
457 :
458 :
459 :
460 0 : subroutine compute_div_C0_contra(zeta,elem,hybrid,nets,nete,nt)
461 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
462 : ! compute C0 divergence. That is, solve:
463 : ! < PHI, zeta > = <PHI, div(elem%state%v >
464 : !
465 : ! input: v (stored in elem()%, in contra-variant coordinates)
466 : ! output: zeta(:,:,:,:)
467 : !
468 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
469 :
470 : type (hybrid_t) , intent(in) :: hybrid
471 : type (element_t) , intent(in), target :: elem(:)
472 : integer :: nt,nets,nete
473 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: zeta
474 : real (kind=r8), dimension(np,np,2) :: ulatlon
475 : real (kind=r8), dimension(np,np) :: v1,v2
476 :
477 : ! local
478 : integer :: k,ie
479 : type (derivative_t) :: deriv
480 :
481 0 : call derivinit(deriv)
482 :
483 0 : do k=1,nlev
484 0 : do ie=nets,nete
485 0 : v1 = elem(ie)%state%v(:,:,1,k,nt)
486 0 : v2 = elem(ie)%state%v(:,:,2,k,nt)
487 0 : ulatlon(:,:,1) = elem(ie)%D(:,:,1,1)*v1 + elem(ie)%D(:,:,1,2)*v2
488 0 : ulatlon(:,:,2) = elem(ie)%D(:,:,2,1)*v1 + elem(ie)%D(:,:,2,2)*v2
489 0 : call divergence_sphere(ulatlon,deriv,elem(ie),zeta(:,:,k,ie))
490 : enddo
491 : enddo
492 :
493 0 : call make_C0(zeta,elem,hybrid,nets,nete)
494 :
495 0 : end subroutine
496 :
497 0 : subroutine compute_zeta_C0_par(zeta,elem,par,nt)
498 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
499 : ! compute C0 vorticity. That is, solve:
500 : ! < PHI, zeta > = <PHI, curl(elem%state%v >
501 : !
502 : ! input: v (stored in elem()%, in lat-lon coordinates)
503 : ! output: zeta(:,:,:,:)
504 : !
505 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
506 : type (parallel_t) :: par
507 : type (element_t) , intent(in), target :: elem(:)
508 : real (kind=r8), dimension(np,np,nlev,nelemd) :: zeta
509 : integer :: nt
510 :
511 : ! local
512 : type (hybrid_t) :: hybrid
513 : integer :: k,i,j,ie,ic
514 : type (derivative_t) :: deriv
515 :
516 : ! single thread
517 0 : hybrid = config_thread_region(par,'serial')
518 :
519 0 : call compute_zeta_C0_hybrid(zeta,elem,hybrid,1,nelemd,nt)
520 :
521 0 : end subroutine
522 :
523 :
524 0 : subroutine compute_div_C0_par(zeta,elem,par,nt)
525 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
526 : ! compute C0 divergence. That is, solve:
527 : ! < PHI, zeta > = <PHI, div(elem%state%v >
528 : !
529 : ! input: v (stored in elem()%, in lat-lon coordinates)
530 : ! output: zeta(:,:,:,:)
531 : !
532 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
533 :
534 : type (parallel_t) :: par
535 : type (element_t) , intent(in), target :: elem(:)
536 : real (kind=r8), dimension(np,np,nlev,nelemd) :: zeta
537 : integer :: nt
538 :
539 : ! local
540 : type (hybrid_t) :: hybrid
541 : integer :: k,i,j,ie,ic
542 : type (derivative_t) :: deriv
543 :
544 : ! single thread
545 0 : hybrid = config_thread_region(par,'serial')
546 :
547 0 : call compute_div_C0_hybrid(zeta,elem,hybrid,1,nelemd,nt)
548 :
549 0 : end subroutine
550 :
551 :
552 :
553 0 : subroutine compute_zeta_C0_hybrid(zeta,elem,hybrid,nets,nete,nt)
554 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
555 : ! compute C0 vorticity. That is, solve:
556 : ! < PHI, zeta > = <PHI, curl(elem%state%v >
557 : !
558 : ! input: v (stored in elem()%, in lat-lon coordinates)
559 : ! output: zeta(:,:,:,:)
560 : !
561 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
562 :
563 : type (hybrid_t) , intent(in) :: hybrid
564 : type (element_t) , intent(in), target :: elem(:)
565 : integer :: nt,nets,nete
566 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: zeta
567 :
568 : ! local
569 : integer :: k,i,j,ie,ic
570 : type (derivative_t) :: deriv
571 :
572 0 : call derivinit(deriv)
573 :
574 0 : do ie=nets,nete
575 : #if (defined COLUMN_OPENMP)
576 : !$omp parallel do num_threads(vert_num_threads) private(k)
577 : #endif
578 0 : do k=1,nlev
579 0 : call vorticity_sphere(elem(ie)%state%v(:,:,:,k,nt),deriv,elem(ie),zeta(:,:,k,ie))
580 : enddo
581 : enddo
582 :
583 0 : call make_C0(zeta,elem,hybrid,nets,nete)
584 :
585 0 : end subroutine
586 :
587 :
588 0 : subroutine compute_div_C0_hybrid(zeta,elem,hybrid,nets,nete,nt)
589 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
590 : ! compute C0 divergence. That is, solve:
591 : ! < PHI, zeta > = <PHI, div(elem%state%v >
592 : !
593 : ! input: v (stored in elem()%, in lat-lon coordinates)
594 : ! output: zeta(:,:,:,:)
595 : !
596 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
597 :
598 : type (hybrid_t) , intent(in) :: hybrid
599 : type (element_t) , intent(in), target :: elem(:)
600 : integer :: nt,nets,nete
601 : real (kind=r8), dimension(np,np,nlev,nets:nete) :: zeta
602 :
603 : ! local
604 : integer :: k,i,j,ie,ic
605 : type (derivative_t) :: deriv
606 :
607 0 : call derivinit(deriv)
608 :
609 0 : do ie=nets,nete
610 : #if (defined COLUMN_OPENMP)
611 : !$omp parallel do num_threads(vert_num_threads) private(k)
612 : #endif
613 0 : do k=1,nlev
614 0 : call divergence_sphere(elem(ie)%state%v(:,:,:,k,nt),deriv,elem(ie),zeta(:,:,k,ie))
615 : enddo
616 : enddo
617 :
618 0 : call make_C0(zeta,elem,hybrid,nets,nete)
619 :
620 0 : end subroutine
621 :
622 :
623 :
624 :
625 :
626 :
627 :
628 :
629 0 : subroutine neighbor_minmax(hybrid,edgeMinMax,nets,nete,min_neigh,max_neigh)
630 :
631 : type (hybrid_t) , intent(in) :: hybrid
632 : type (EdgeBuffer_t) , intent(inout) :: edgeMinMax
633 : integer :: nets,nete
634 : real (kind=r8) :: min_neigh(nlev,qsize,nets:nete)
635 : real (kind=r8) :: max_neigh(nlev,qsize,nets:nete)
636 : integer :: kblk, qblk
637 : ! local
638 : integer:: ie, q, k, kptr
639 : integer:: kbeg, kend, qbeg, qend
640 :
641 0 : call get_loop_ranges(hybrid,kbeg=kbeg,kend=kend,qbeg=qbeg,qend=qend)
642 :
643 0 : kblk = kend - kbeg + 1 ! calculate size of the block of vertical levels
644 0 : qblk = qend - qbeg + 1 ! calculate size of the block of tracers
645 :
646 0 : do ie=nets,nete
647 0 : do q = qbeg, qend
648 0 : kptr = nlev*(q - 1) + kbeg - 1
649 0 : call edgeSpack(edgeMinMax,min_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
650 0 : kptr = qsize*nlev + nlev*(q - 1) + kbeg - 1
651 0 : call edgeSpack(edgeMinMax,max_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
652 : enddo
653 : enddo
654 :
655 0 : call bndry_exchange(hybrid,edgeMinMax,location='neighbor_minmax')
656 :
657 0 : do ie=nets,nete
658 0 : do q=qbeg,qend
659 0 : kptr = nlev*(q - 1) + kbeg - 1
660 0 : call edgeSunpackMIN(edgeMinMax,min_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
661 0 : kptr = qsize*nlev + nlev*(q - 1) + kbeg - 1
662 0 : call edgeSunpackMAX(edgeMinMax,max_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
663 0 : do k=kbeg,kend
664 0 : min_neigh(k,q,ie) = max(min_neigh(k,q,ie),0.0_r8)
665 : enddo
666 : enddo
667 : enddo
668 :
669 0 : end subroutine neighbor_minmax
670 :
671 :
672 0 : subroutine neighbor_minmax_start(hybrid,edgeMinMax,nets,nete,min_neigh,max_neigh)
673 :
674 : type (hybrid_t) , intent(in) :: hybrid
675 : type (EdgeBuffer_t) , intent(inout) :: edgeMinMax
676 : integer :: nets,nete
677 : real (kind=r8) :: min_neigh(nlev,qsize,nets:nete)
678 : real (kind=r8) :: max_neigh(nlev,qsize,nets:nete)
679 : integer :: kblk, qblk
680 : integer :: kbeg, kend, qbeg, qend
681 :
682 : ! local
683 : integer :: ie,q, k,kptr
684 :
685 0 : call get_loop_ranges(hybrid,kbeg=kbeg,kend=kend,qbeg=qbeg,qend=qend)
686 :
687 0 : kblk = kend - kbeg + 1 ! calculate size of the block of vertical levels
688 0 : qblk = qend - qbeg + 1 ! calculate size of the block of tracers
689 :
690 0 : do ie=nets,nete
691 0 : do q=qbeg, qend
692 0 : kptr = nlev*(q - 1) + kbeg - 1
693 0 : call edgeSpack(edgeMinMax,min_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
694 0 : kptr = qsize*nlev + nlev*(q - 1) + kbeg - 1
695 0 : call edgeSpack(edgeMinMax,max_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
696 : enddo
697 : enddo
698 :
699 0 : call bndry_exchange_start(hybrid,edgeMinMax,location='neighbor_minmax_start')
700 :
701 0 : end subroutine neighbor_minmax_start
702 :
703 0 : subroutine neighbor_minmax_finish(hybrid,edgeMinMax,nets,nete,min_neigh,max_neigh)
704 :
705 : type (hybrid_t) , intent(in) :: hybrid
706 : type (EdgeBuffer_t) , intent(inout) :: edgeMinMax
707 : integer :: nets,nete
708 : real (kind=r8) :: min_neigh(nlev,qsize,nets:nete)
709 : real (kind=r8) :: max_neigh(nlev,qsize,nets:nete)
710 : integer :: kblk, qblk
711 : integer :: ie,q, k,kptr
712 : integer :: kbeg, kend, qbeg, qend
713 :
714 0 : call get_loop_ranges(hybrid,kbeg=kbeg,kend=kend,qbeg=qbeg,qend=qend)
715 :
716 0 : kblk = kend - kbeg + 1 ! calculate size of the block of vertical levels
717 0 : qblk = qend - qbeg + 1 ! calculate size of the block of tracers
718 :
719 0 : call bndry_exchange_finish(hybrid,edgeMinMax,location='neighbor_minmax_finish')
720 :
721 0 : do ie=nets,nete
722 0 : do q=qbeg, qend
723 0 : kptr = nlev*(q - 1) + kbeg - 1
724 0 : call edgeSunpackMIN(edgeMinMax,min_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
725 0 : kptr = qsize*nlev + nlev*(q - 1) + kbeg - 1
726 0 : call edgeSunpackMAX(edgeMinMax,max_neigh(kbeg:kend,q,ie),kblk,kptr,ie)
727 0 : do k=kbeg,kend
728 0 : min_neigh(k,q,ie) = max(min_neigh(k,q,ie),0.0_r8)
729 : enddo
730 : enddo
731 : enddo
732 :
733 0 : end subroutine neighbor_minmax_finish
734 :
735 : end module viscosity_mod
|