120 qmax, rx, ry, rz, xc, yc, zc, &
121 x, y, z, lcmesh3D, elem, &
122 IntrpPolyOrder_h, IntrpPolyOrder_v, &
123 z_func_type, z_func_params, cosbell_exponent )
129 real(rp),
intent(out) :: q(elem%np,lcmesh3d%nea)
130 real(rp),
intent(in) :: qmax
131 real(rp),
intent(in) :: rx, ry, rz
132 real(rp),
intent(in) :: xc, yc, zc
133 real(rp),
intent(in) :: x(elem%np,lcmesh3d%ne)
134 real(rp),
intent(in) :: y(elem%np,lcmesh3d%ne)
135 real(rp),
intent(in) :: z(elem%np,lcmesh3d%ne)
136 integer,
intent(in) :: intrppolyorder_h
137 integer,
intent(in) :: intrppolyorder_v
138 character(len=*),
optional,
intent(in) :: z_func_type
139 real(rp),
optional,
intent(in) :: z_func_params(:)
140 integer,
intent(in),
optional :: cosbell_exponent
145 real(rp),
allocatable :: x_intrp(:,:), y_intrp(:,:), z_intrp(:,:)
146 real(rp),
allocatable :: r_intrp(:)
147 real(rp),
allocatable :: z_func(:,:)
148 real(rp) :: vx(elem%nv), vy(elem%nv), vz(elem%nv)
150 real(rp),
allocatable :: intrpmat(:,:)
151 real(rp),
allocatable :: q_intrp(:)
156 if (
present(cosbell_exponent) )
then
157 exponent = cosbell_exponent
162 call elem_intrp%Init( intrppolyorder_h, intrppolyorder_v, .false. )
164 allocate( intrpmat(elem%Np,elem_intrp%Np) )
167 allocate( x_intrp(elem_intrp%Np,lcmesh3d%Ne), y_intrp(elem_intrp%Np,lcmesh3d%Ne), z_intrp(elem_intrp%Np,lcmesh3d%Ne) )
168 allocate( z_func(elem_intrp%Np,lcmesh3d%Ne))
169 allocate( r_intrp(elem_intrp%Np) )
170 allocate( q_intrp(elem_intrp%Np) )
177 do ke=lcmesh3d%NeS, lcmesh3d%NeE
179 vx(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),1)
180 vy(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),2)
181 vz(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),3)
182 x_intrp(:,ke) = vx(1) + 0.5_rp * ( elem_intrp%x1(:) + 1.0_rp ) * ( vx(2) - vx(1) )
183 y_intrp(:,ke) = vy(1) + 0.5_rp * ( elem_intrp%x2(:) + 1.0_rp ) * ( vy(4) - vy(1) )
184 z_intrp(:,ke) = vz(1) + 0.5_rp * ( elem_intrp%x3(:) + 1.0_rp ) * ( vz(5) - vz(1) )
186 z_func(:,ke) = 1.0_rp
191 if (
present(z_func_type) )
then
192 select case(z_func_type)
195 do ke=lcmesh3d%NeS, lcmesh3d%NeE
196 z_func(:,ke) = sin( z_func_params(1) * pi * z_intrp(:,ke) / z_func_params(2) )
203 do ke=lcmesh3d%NeS, lcmesh3d%NeE
205 ( (x_intrp(:,ke) - xc) / rx )**2 &
206 + ( (y_intrp(:,ke) - yc) / ry )**2 &
207 + ( (z_intrp(:,ke) - zc) / rz )**2 )
209 where( r_intrp(:) <= 1.0_rp )
210 q_intrp(:) = qmax * ( 0.5_rp * (1.0_rp + cos( pi * r_intrp(:) ) ) )**exponent
214 q(:,ke) = matmul(intrpmat, q_intrp(:) * z_func(:,ke))
219 call elem_intrp%Final()
234 qmax, rh, lonc, latc, rplanet, &
235 x, y, z, lcmesh3D, elem, &
236 IntrpPolyOrder_h, IntrpPolyOrder_v, &
237 z_func_type, z_func_params, cosbell_exponent )
245 real(rp),
intent(out) :: q(elem%np,lcmesh3d%nea)
246 real(rp),
intent(in) :: qmax
247 real(rp),
intent(in) :: rh
248 real(rp),
intent(in) :: lonc, latc
249 real(rp),
intent(in) :: rplanet
250 real(rp),
intent(in) :: x(elem%np,lcmesh3d%ne)
251 real(rp),
intent(in) :: y(elem%np,lcmesh3d%ne)
252 real(rp),
intent(in) :: z(elem%np,lcmesh3d%ne)
253 integer,
intent(in) :: intrppolyorder_h
254 integer,
intent(in) :: intrppolyorder_v
255 character(len=*),
optional,
intent(in) :: z_func_type
256 real(rp),
optional,
intent(in) :: z_func_params(:)
257 integer,
intent(in),
optional :: cosbell_exponent
262 real(rp),
allocatable :: x_intrp(:,:), y_intrp(:,:), z_intrp(:,:)
263 real(rp),
allocatable :: gam_intrp(:,:)
264 real(rp),
allocatable :: lon_intrp(:,:), lat_intrp(:,:)
265 real(rp),
allocatable :: z_func(:,:)
266 real(rp),
allocatable :: r_intrp(:)
267 real(rp) :: vx(elem%nv), vy(elem%nv), vz(elem%nv)
269 real(rp),
allocatable :: intrpmat(:,:)
270 real(rp),
allocatable :: q_intrp(:)
276 if (
present(cosbell_exponent) )
then
277 exponent = cosbell_exponent
282 call elem_intrp%Init( intrppolyorder_h, intrppolyorder_v, .false. )
284 allocate( intrpmat(elem%Np,elem_intrp%Np) )
287 allocate( x_intrp(elem_intrp%Np,lcmesh3d%Ne), y_intrp(elem_intrp%Np,lcmesh3d%Ne), z_intrp(elem_intrp%Np,lcmesh3d%Ne) )
288 allocate( gam_intrp(elem_intrp%Np,lcmesh3d%Ne) )
289 allocate( lon_intrp(elem_intrp%Np,lcmesh3d%Ne), lat_intrp(elem_intrp%Np,lcmesh3d%Ne) )
290 allocate( z_func(elem_intrp%Np,lcmesh3d%Ne) )
291 allocate( r_intrp(elem_intrp%Np) )
292 allocate( q_intrp(elem_intrp%Np) )
296 do ke=lcmesh3d%NeS, lcmesh3d%NeE
297 vx(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),1)
298 vy(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),2)
299 vz(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),3)
300 x_intrp(:,ke) = vx(1) + 0.5_rp * ( elem_intrp%x1(:) + 1.0_rp ) * ( vx(2) - vx(1) )
301 y_intrp(:,ke) = vy(1) + 0.5_rp * ( elem_intrp%x2(:) + 1.0_rp ) * ( vy(4) - vy(1) )
302 z_intrp(:,ke) = vz(1) + 0.5_rp * ( elem_intrp%x3(:) + 1.0_rp ) * ( vz(5) - vz(1) )
303 gam_intrp(:,ke) = 1.0_rp
305 z_func(:,ke) = 1.0_rp
309 elem_intrp%Np * lcmesh3d%Ne, &
310 lon_intrp(:,:), lat_intrp(:,:) )
313 if (
present(z_func_type) )
then
314 select case(z_func_type)
317 do ke=lcmesh3d%NeS, lcmesh3d%NeE
318 z_func(:,ke) = sin( z_func_params(1) * pi * z_intrp(:,ke) / z_func_params(2) )
324 do ke=lcmesh3d%NeS, lcmesh3d%NeE
327 r_intrp(:) = rplanet / rh * acos( sin(latc) * sin(lat_intrp(:,ke)) + cos(latc) * cos(lat_intrp(:,ke)) * cos(lon_intrp(:,ke) - lonc) )
328 where( r_intrp(:) <= 1.0_rp )
329 q_intrp(:) = qmax * ( 0.5_rp * (1.0_rp + cos( pi * r_intrp(:) ) ) )**exponent
335 q(:,ke) = matmul(intrpmat, q_intrp(:) * z_func(:,ke))
338 call elem_intrp%Final()
349 func, IntrpPolyOrder_h, IntrpPolyOrder_v, &
355 real(rp),
intent(out) :: q(elem%np,lcmesh3d%nea)
356 integer,
intent(in) :: intrppolyorder_h
357 integer,
intent(in) :: intrppolyorder_v
360 subroutine func( q_intrp, &
361 x, y, z, elem_intrp )
365 real(rp),
intent(out) :: q_intrp(elem_intrp%np)
366 real(rp),
intent(in) :: x(elem_intrp%np)
367 real(rp),
intent(in) :: y(elem_intrp%np)
368 real(rp),
intent(in) :: z(elem_intrp%np)
373 real(rp),
allocatable :: x_intrp(:,:), y_intrp(:,:), z_intrp(:,:)
374 real(rp) :: vx(elem%nv), vy(elem%nv), vz(elem%nv)
376 real(rp),
allocatable :: intrpmat(:,:)
377 real(rp),
allocatable :: q_intrp(:)
382 call elem_intrp%Init( intrppolyorder_h, intrppolyorder_v, .false. )
384 allocate( intrpmat(elem%Np,elem_intrp%Np) )
387 allocate( x_intrp(elem_intrp%Np,lcmesh3d%Ne), y_intrp(elem_intrp%Np,lcmesh3d%Ne), z_intrp(elem_intrp%Np,lcmesh3d%Ne) )
388 allocate( q_intrp(elem_intrp%Np) )
391 do ke=lcmesh3d%NeS, lcmesh3d%NeE
392 vx(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),1)
393 vy(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),2)
394 vz(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),3)
395 x_intrp(:,ke) = vx(1) + 0.5_rp * ( elem_intrp%x1(:) + 1.0_rp ) * ( vx(2) - vx(1) )
396 y_intrp(:,ke) = vy(1) + 0.5_rp * ( elem_intrp%x2(:) + 1.0_rp ) * ( vy(4) - vy(1) )
397 z_intrp(:,ke) = vz(1) + 0.5_rp * ( elem_intrp%x3(:) + 1.0_rp ) * ( vz(5) - vz(1) )
401 do ke=lcmesh3d%NeS, lcmesh3d%NeE
403 call func( q_intrp, &
404 x_intrp(:,ke), y_intrp(:,ke), z_intrp(:,ke), &
408 q(:,ke) = matmul( intrpmat, q_intrp )
411 call elem_intrp%Final()
421 func, IntrpPolyOrder_h, IntrpPolyOrder_v, &
422 lcmesh3D, elem, rplanet )
430 real(rp),
intent(out) :: q(elem%np,lcmesh3d%nea)
431 integer,
intent(in) :: intrppolyorder_h
432 integer,
intent(in) :: intrppolyorder_v
433 real(rp),
intent(in) :: rplanet
436 subroutine func( q_intrp, &
437 lon, lat, z, elem_intrp, rplanet_ )
441 real(rp),
intent(out) :: q_intrp(elem_intrp%np)
442 real(rp),
intent(in) :: lon(elem_intrp%np)
443 real(rp),
intent(in) :: lat(elem_intrp%np)
444 real(rp),
intent(in) :: z(elem_intrp%np)
445 real(rp),
intent(in) :: rplanet_
450 real(rp),
allocatable :: x_intrp(:,:), y_intrp(:,:), z_intrp(:,:)
451 real(rp),
allocatable :: gam_intrp(:,:)
452 real(rp),
allocatable :: lon_intrp(:,:), lat_intrp(:,:)
453 real(rp) :: vx(elem%nv), vy(elem%nv), vz(elem%nv)
455 real(rp),
allocatable :: intrpmat(:,:)
456 real(rp),
allocatable :: q_intrp(:)
461 call elem_intrp%Init( intrppolyorder_h, intrppolyorder_v, .false. )
463 allocate( intrpmat(elem%Np,elem_intrp%Np) )
466 allocate( x_intrp(elem_intrp%Np,lcmesh3d%Ne), y_intrp(elem_intrp%Np,lcmesh3d%Ne), z_intrp(elem_intrp%Np,lcmesh3d%Ne) )
467 allocate( gam_intrp(elem_intrp%Np,lcmesh3d%Ne) )
468 allocate( lon_intrp(elem_intrp%Np,lcmesh3d%Ne), lat_intrp(elem_intrp%Np,lcmesh3d%Ne) )
469 allocate( q_intrp(elem_intrp%Np) )
472 do ke=lcmesh3d%NeS, lcmesh3d%NeE
473 vx(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),1)
474 vy(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),2)
475 vz(:) = lcmesh3d%pos_ev(lcmesh3d%EToV(ke,:),3)
476 x_intrp(:,ke) = vx(1) + 0.5_rp * ( elem_intrp%x1(:) + 1.0_rp ) * ( vx(2) - vx(1) )
477 y_intrp(:,ke) = vy(1) + 0.5_rp * ( elem_intrp%x2(:) + 1.0_rp ) * ( vy(4) - vy(1) )
478 z_intrp(:,ke) = vz(1) + 0.5_rp * ( elem_intrp%x3(:) + 1.0_rp ) * ( vz(5) - vz(1) )
480 gam_intrp(:,ke) = 1.0_rp
484 elem_intrp%Np * lcmesh3d%Ne, &
485 lon_intrp(:,:), lat_intrp(:,:) )
488 do ke=lcmesh3d%NeS, lcmesh3d%NeE
490 call func( q_intrp, &
491 lon_intrp(:,ke), lat_intrp(:,ke), z_intrp(:,ke), &
492 elem_intrp, rplanet )
495 q(:,ke) = matmul( intrpmat, q_intrp )
498 call elem_intrp%Final()
subroutine, public mkinitutil_calc_cosinebell(q, qmax, rx, ry, rz, xc, yc, zc, x, y, z, lcmesh3d, elem, intrppolyorder_h, intrppolyorder_v, z_func_type, z_func_params, cosbell_exponent)
Calculate the distribution function of a cosine bell in regional domain.
subroutine, public mkinitutil_calc_cosinebell_global(q, qmax, rh, lonc, latc, rplanet, x, y, z, lcmesh3d, elem, intrppolyorder_h, intrppolyorder_v, z_func_type, z_func_params, cosbell_exponent)
Calculate the distribution function of a cosine bell in global domain.