62 real,
parameter::
r3 = 1./3.
69 real,
parameter::
b1 = 0.0375
70 real,
parameter::
b2 = -7./30.
71 real,
parameter::
b3 = -23./120.
72 real,
parameter::
b4 = 13./30.
73 real,
parameter::
b5 = -11./240.
76 real,
parameter::
b1 = 1./30.
77 real,
parameter::
b2 = -13./60.
78 real,
parameter::
b3 = -13./60.
79 real,
parameter::
b4 = 0.45
80 real,
parameter::
b5 = -0.05
82 real,
parameter::
t11 = 27./28.,
t12 = -13./28.,
t13=3./7.
83 real,
parameter::
s11 = 11./14.,
s14 = 4./7.,
s15=3./14.
88 real,
parameter::
c1 = -2./14.
89 real,
parameter::
c2 = 11./14.
90 real,
parameter::
c3 = 5./14.
94 real,
parameter::
p1 = 7./12.
95 real,
parameter::
p2 = -1./12.
108 subroutine fv_tp_2d(q, crx, cry, npx, npy, hord, fx, fy, xfx, yfx, &
109 gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx, mfy, mass, nord, damp_c)
111 integer,
intent(in):: npx, npy
112 integer,
intent(in)::hord
114 real,
intent(in):: crx(bd%is:bd%ie+1,bd%jsd:bd%jed)
115 real,
intent(in):: xfx(bd%is:bd%ie+1,bd%jsd:bd%jed)
116 real,
intent(in):: cry(bd%isd:bd%ied,bd%js:bd%je+1 )
117 real,
intent(in):: yfx(bd%isd:bd%ied,bd%js:bd%je+1 )
118 real,
intent(in):: ra_x(bd%is:bd%ie,bd%jsd:bd%jed)
119 real,
intent(in):: ra_y(bd%isd:bd%ied,bd%js:bd%je)
120 real,
intent(inout):: q(bd%isd:bd%ied,bd%jsd:bd%jed)
121 real,
intent(out)::fx(bd%is:bd%ie+1 ,bd%js:bd%je)
122 real,
intent(out)::fy(bd%is:bd%ie, bd%js:bd%je+1 )
126 real,
intent(in):: lim_fac
127 logical,
intent(in):: regional
129 real,
OPTIONAL,
intent(in):: mfx(bd%is:bd%ie+1,bd%js:bd%je )
130 real,
OPTIONAL,
intent(in):: mfy(bd%is:bd%ie ,bd%js:bd%je+1)
131 real,
OPTIONAL,
intent(in):: mass(bd%isd:bd%ied,bd%jsd:bd%jed)
132 real,
OPTIONAL,
intent(in):: damp_c
133 integer,
OPTIONAL,
intent(in):: nord
135 integer ord_ou, ord_in
136 real q_i(bd%isd:bd%ied,bd%js:bd%je)
137 real q_j(bd%is:bd%ie,bd%jsd:bd%jed)
138 real fx2(bd%is:bd%ie+1,bd%jsd:bd%jed)
139 real fy2(bd%isd:bd%ied,bd%js:bd%je+1)
140 real fyy(bd%isd:bd%ied,bd%js:bd%je+1)
141 real fx1(bd%is:bd%ie+1)
145 integer:: is, ie, js, je, isd, ied, jsd, jed
156 if ( hord == 10 )
then 163 if (.not. (regional))
call copy_corners(q, npx, npy, 2, gridstruct%nested, bd, &
164 gridstruct%sw_corner, gridstruct%se_corner, gridstruct%nw_corner, gridstruct%ne_corner)
166 call yppm(fy2, q, cry, ord_in, isd,ied,isd,ied, js,je,jsd,jed, npx,npy, gridstruct%dya, gridstruct%nested, gridstruct%grid_type, lim_fac,regional)
170 fyy(i,j) = yfx(i,j) * fy2(i,j)
175 q_i(i,j) = (q(i,j)*gridstruct%area(i,j) + fyy(i,j)-fyy(i,j+1))/ra_y(i,j)
179 call xppm(fx, q_i, crx(is,js), ord_ou, is,ie,isd,ied, js,je,jsd,jed, npx,npy, gridstruct%dxa, gridstruct%nested, gridstruct%grid_type, lim_fac,regional)
181 if (.not. (regional))
call copy_corners(q, npx, npy, 1, gridstruct%nested, bd, &
182 gridstruct%sw_corner, gridstruct%se_corner, gridstruct%nw_corner, gridstruct%ne_corner)
184 call xppm(fx2, q, crx, ord_in, is,ie,isd,ied, jsd,jed,jsd,jed, npx,npy, gridstruct%dxa, gridstruct%nested, gridstruct%grid_type, lim_fac,regional)
188 fx1(i) = xfx(i,j) * fx2(i,j)
191 q_j(i,j) = (q(i,j)*gridstruct%area(i,j) + fx1(i)-fx1(i+1))/ra_x(i,j)
195 call yppm(fy, q_j, cry, ord_ou, is,ie,isd,ied, js,je,jsd,jed, npx, npy, gridstruct%dya, gridstruct%nested, gridstruct%grid_type, lim_fac,regional)
201 if (
present(mfx) .and.
present(mfy) )
then 207 fx(i,j) = 0.5*(fx(i,j) + fx2(i,j)) * mfx(i,j)
212 fy(i,j) = 0.5*(fy(i,j) + fy2(i,j)) * mfy(i,j)
215 if (
present(nord) .and.
present(damp_c) .and.
present(mass) )
then 216 if ( damp_c > 1.e-4 )
then 217 damp = (damp_c * gridstruct%da_min)**(nord+1)
218 call deln_flux(nord, is,ie,js,je, npx, npy, damp, q, fx, fy, gridstruct,regional, bd, mass)
227 fx(i,j) = 0.5*(fx(i,j) + fx2(i,j)) * xfx(i,j)
232 fy(i,j) = 0.5*(fy(i,j) + fy2(i,j)) * yfx(i,j)
235 if (
present(nord) .and.
present(damp_c) )
then 236 if ( damp_c > 1.e-4 )
then 237 damp = (damp_c * gridstruct%da_min)**(nord+1)
238 call deln_flux(nord, is,ie,js,je, npx, npy, damp, q, fx, fy, gridstruct, regional, bd)
247 sw_corner, se_corner, nw_corner, ne_corner)
249 integer,
intent(in):: npx, npy, dir
250 real,
intent(inout):: q(bd%isd:bd%ied,bd%jsd:bd%jed)
251 logical,
intent(IN) :: nested, sw_corner, se_corner, nw_corner, ne_corner
258 if ( sw_corner )
then 265 if ( se_corner )
then 268 q(i,j) = q(npy-j,i-npx+1)
272 if ( ne_corner )
then 275 q(i,j) = q(j,2*npx-1-i)
279 if ( nw_corner )
then 282 q(i,j) = q(npy-j,i-1+npx)
287 elseif ( dir == 2 )
then 290 if ( sw_corner )
then 297 if ( se_corner )
then 300 q(i,j) = q(npy+j-1,npx-i)
304 if ( ne_corner )
then 307 q(i,j) = q(2*npy-1-j,i)
311 if ( nw_corner )
then 314 q(i,j) = q(j+1-npx,npy-i)
323 subroutine xppm(flux, q, c, iord, is,ie,isd,ied, jfirst,jlast,jsd,jed, npx, npy, dxa, nested, grid_type, lim_fac,regional)
324 integer,
INTENT(IN) :: is, ie, isd, ied, jsd, jed
325 integer,
INTENT(IN) :: jfirst, jlast
326 integer,
INTENT(IN) :: iord
327 integer,
INTENT(IN) :: npx, npy
328 real ,
INTENT(IN) :: q(isd:ied,jfirst:jlast)
329 real ,
INTENT(IN) :: c(is:ie+1,jfirst:jlast)
330 real ,
intent(IN) :: dxa(isd:ied,jsd:jed)
331 logical,
intent(IN) :: nested,regional
332 integer,
intent(IN) :: grid_type
333 real ,
intent(IN) :: lim_fac
335 real ,
INTENT(OUT) :: flux(is:ie+1,jfirst:jlast)
337 real,
dimension(is-1:ie+1):: bl, br, b0
339 real,
dimension(is:ie+1):: fx0, fx1, xt1
340 logical,
dimension(is-1:ie+1):: smt5, smt6
341 logical,
dimension(is:ie+1):: hi5, hi6
345 integer:: i, j, ie3, is1, ie1, mord
346 real:: x0, x1, xt, qtmp, pmp_1, lac_1, pmp_2, lac_2
348 if ( .not. (nested .or. regional) .and. grid_type<3 )
then 349 is1 = max(3,is-1); ie3 = min(npx-2,ie+2)
350 ie1 = min(npx-3,ie+1)
352 is1 = is-1; ie3 = ie+2
357 do 666 j=jfirst,jlast
368 al(i) =
p1*(q1(i-1)+q1(i)) +
p2*(q1(i-2)+q1(i+1))
371 if ( .not. (nested .or. regional) .and. grid_type<3 )
then 373 al(0) =
c1*q1(-2) +
c2*q1(-1) +
c3*q1(0)
374 al(1) = 0.5*(((2.*dxa(0,j)+dxa(-1,j))*q1(0)-dxa(0,j)*q1(-1))/(dxa(-1,j)+dxa(0,j)) &
375 + ((2.*dxa(1,j)+dxa( 2,j))*q1(1)-dxa(1,j)*q1( 2))/(dxa(1, j)+dxa(2,j)))
376 al(2) =
c3*q1(1) +
c2*q1(2) +
c1*q1(3)
378 if ( (ie+1)==npx )
then 379 al(npx-1) =
c1*q1(npx-3) +
c2*q1(npx-2) +
c3*q1(npx-1)
380 al(npx) = 0.5*(((2.*dxa(npx-1,j)+dxa(npx-2,j))*q1(npx-1)-dxa(npx-1,j)*q1(npx-2))/(dxa(npx-2,j)+dxa(npx-1,j)) &
381 + ((2.*dxa(npx, j)+dxa(npx+1,j))*q1(npx )-dxa(npx, j)*q1(npx+1))/(dxa(npx, j)+dxa(npx+1,j)))
382 al(npx+1) =
c3*q1(npx) +
c2*q1(npx+1) +
c1*q1(npx+2)
388 al(i) = max(0., al(i))
394 bl(i) = al(i) - q1(i)
395 br(i) = al(i+1) - q1(i)
396 b0(i) = bl(i) + br(i)
397 smt5(i) = abs(lim_fac*b0(i)) < abs(bl(i)-br(i))
401 if ( c(i,j) > 0. )
then 402 fx1(i) = (1.-c(i,j))*(br(i-1) - c(i,j)*b0(i-1))
405 fx1(i) = (1.+c(i,j))*(bl(i) + c(i,j)*b0(i))
408 if (smt5(i-1).or.smt5(i)) flux(i,j) = flux(i,j) + fx1(i)
411 elseif ( mord==2 )
then 418 flux(i,j) = qtmp + (1.-xt)*(al(i)-qtmp-xt*(al(i-1)+al(i)-(qtmp+qtmp)))
421 flux(i,j) = qtmp + (1.+xt)*(al(i)-qtmp+xt*(al(i)+al(i+1)-(qtmp+qtmp)))
429 elseif ( mord==3 )
then 432 bl(i) = al(i) - q1(i)
433 br(i) = al(i+1) - q1(i)
434 b0(i) = bl(i) + br(i)
436 xt = abs(bl(i)-br(i))
443 hi5(i) = smt5(i-1) .and. smt5(i)
444 hi6(i) = smt6(i-1) .or. smt6(i)
447 if ( xt1(i) > 0. )
then 449 fx1(i) = br(i-1) - xt1(i)*b0(i-1)
450 elseif ( hi5(i) )
then 451 fx1(i) = sign(min(abs(bl(i-1)),abs(br(i-1))), br(i-1))
453 flux(i,j) = q1(i-1) + (1.-xt1(i))*fx1(i)
456 fx1(i) = bl(i) + xt1(i)*b0(i)
457 elseif ( hi5(i) )
then 458 fx1(i) = sign(min(abs(bl(i)), abs(br(i))), bl(i))
460 flux(i,j) = q1(i) + (1.+xt1(i))*fx1(i)
464 elseif ( mord==4 )
then 467 bl(i) = al(i) - q1(i)
468 br(i) = al(i+1) - q1(i)
469 b0(i) = bl(i) + br(i)
471 xt = abs(bl(i)-br(i))
477 hi5(i) = smt5(i-1) .and. smt5(i)
478 hi6(i) = smt6(i-1) .or. smt6(i)
479 hi5(i) = hi5(i) .or. hi6(i)
483 if ( xt1(i) > 0. )
then 484 fx1(i) = (1.-xt1(i))*(br(i-1) - xt1(i)*b0(i-1))
487 fx1(i) = (1.+xt1(i))*(bl(i) + xt1(i)*b0(i))
490 if ( hi5(i) ) flux(i,j) = flux(i,j) + fx1(i)
497 bl(i) = al(i) - q1(i)
498 br(i) = al(i+1) - q1(i)
499 b0(i) = bl(i) + br(i)
500 smt5(i) = bl(i)*br(i) < 0.
504 bl(i) = al(i) - q1(i)
505 br(i) = al(i+1) - q1(i)
506 b0(i) = bl(i) + br(i)
507 smt5(i) = 3.*abs(b0(i)) < abs(bl(i)-br(i))
513 if ( c(i,j) > 0. )
then 514 fx1(i) = (1.-c(i,j))*(br(i-1) - c(i,j)*b0(i-1))
517 fx1(i) = (1.+c(i,j))*(bl(i) + c(i,j)*b0(i))
520 if (smt5(i-1).or.smt5(i)) flux(i,j) = flux(i,j) + fx1(i)
534 xt = 0.25*(q1(i+1) - q1(i-1))
535 dm(i) = sign(min(abs(xt), max(q1(i-1), q1(i), q1(i+1)) - q1(i), &
536 q1(i) - min(q1(i-1), q1(i), q1(i+1))), xt)
539 al(i) = 0.5*(q1(i-1)+q1(i)) +
r3*(dm(i-1)-dm(i))
545 bl(i) = -sign(min(abs(xt), abs(al(i )-q1(i))), xt)
546 br(i) = sign(min(abs(xt), abs(al(i+1)-q1(i))), xt)
548 elseif ( iord==11 )
then 552 bl(i) = -sign(min(abs(xt), abs(al(i )-q1(i))), xt)
553 br(i) = sign(min(abs(xt), abs(al(i+1)-q1(i))), xt)
557 dq(i) = 2.*(q1(i+1) - q1(i))
560 bl(i) = al(i ) - q1(i)
561 br(i) = al(i+1) - q1(i)
562 if ( abs(dm(i-1))+abs(dm(i))+abs(dm(i+1)) <
near_zero )
then 565 elseif( abs(3.*(bl(i)+br(i))) > abs(bl(i)-br(i)) )
then 567 lac_2 = pmp_2 - 0.75*dq(i-2)
568 br(i) = min( max(0., pmp_2, lac_2), max(br(i), min(0., pmp_2, lac_2)) )
570 lac_1 = pmp_1 + 0.75*dq(i+1)
571 bl(i) = min( max(0., pmp_1, lac_1), max(bl(i), min(0., pmp_1, lac_1)) )
576 if(iord==9 .or. iord==13)
call pert_ppm(ie1-is1+1, q1(is1), bl(is1), br(is1), 0)
578 if (.not. (nested .or. regional) .and. grid_type<3)
then 580 bl(0) =
s14*dm(-1) +
s11*(q1(-1)-q1(0))
582 xt = 0.5*(((2.*dxa(0,j)+dxa(-1,j))*q1(0)-dxa(0,j)*q1(-1))/(dxa(-1,j)+dxa(0,j)) &
583 + ((2.*dxa(1,j)+dxa( 2,j))*q1(1)-dxa(1,j)*q1( 2))/(dxa(1, j)+dxa(2,j)))
585 xt = max(xt, min(q1(-1),q1(0),q1(1),q1(2)))
586 xt = min(xt, max(q1(-1),q1(0),q1(1),q1(2)))
594 br(2) = al(3) - q1(2)
595 call pert_ppm(3, q1(0), bl(0), br(0), 1)
597 if ( (ie+1)==npx )
then 598 bl(npx-2) = al(npx-2) - q1(npx-2)
600 xt =
s15*q1(npx-1) +
s11*q1(npx-2) +
s14*dm(npx-2)
601 br(npx-2) = xt - q1(npx-2)
602 bl(npx-1) = xt - q1(npx-1)
604 xt = 0.5*(((2.*dxa(npx-1,j)+dxa(npx-2,j))*q1(npx-1)-dxa(npx-1,j)*q1(npx-2))/(dxa(npx-2,j)+dxa(npx-1,j)) &
605 + ((2.*dxa(npx, j)+dxa(npx+1,j))*q1(npx )-dxa(npx, j)*q1(npx+1))/(dxa(npx, j)+dxa(npx+1,j)))
607 xt = max(xt, min(q1(npx-2),q1(npx-1),q1(npx),q1(npx+1)))
608 xt = min(xt, max(q1(npx-2),q1(npx-1),q1(npx),q1(npx+1)))
610 br(npx-1) = xt - q1(npx-1)
611 bl(npx ) = xt - q1(npx )
613 br(npx) =
s11*(q1(npx+1)-q1(npx)) -
s14*dm(npx+1)
614 call pert_ppm(3, q1(npx-2), bl(npx-2), br(npx-2), 1)
622 flux(i,j) = q1(i-1) + (1.-c(i,j))*(br(i-1)-c(i,j)*(bl(i-1)+br(i-1)))
624 flux(i,j) = q1(i ) + (1.+c(i,j))*(bl(i )+c(i,j)*(bl(i)+br(i)))
632 subroutine yppm(flux, q, c, jord, ifirst,ilast, isd,ied, js,je,jsd,jed, npx, npy, dya, nested, grid_type, lim_fac,regional)
633 integer,
INTENT(IN) :: ifirst,ilast
634 integer,
INTENT(IN) :: isd,ied, js,je,jsd,jed
635 integer,
INTENT(IN) :: jord
636 integer,
INTENT(IN) :: npx, npy
637 real ,
INTENT(IN) :: q(ifirst:ilast,jsd:jed)
638 real ,
intent(in) :: c(isd:ied,js:je+1 )
639 real ,
INTENT(OUT):: flux(ifirst:ilast,js:je+1)
640 real ,
intent(IN) :: dya(isd:ied,jsd:jed)
641 logical,
intent(IN) :: nested,regional
642 integer,
intent(IN) :: grid_type
643 real ,
intent(IN) :: lim_fac
645 real:: dm(ifirst:ilast,js-2:je+2)
646 real:: al(ifirst:ilast,js-1:je+2)
647 real,
dimension(ifirst:ilast,js-1:je+1):: bl, br, b0
648 real:: dq(ifirst:ilast,js-3:je+2)
649 real,
dimension(ifirst:ilast):: fx0, fx1, xt1
650 logical,
dimension(ifirst:ilast,js-1:je+1):: smt5, smt6
651 logical,
dimension(ifirst:ilast):: hi5, hi6
652 real:: x0, xt, qtmp, pmp_1, lac_1, pmp_2, lac_2, r1
653 integer:: i, j, js1, je3, je1, mord
655 if ( .not. (nested .or. regional) .and. grid_type < 3 )
then 657 js1 = max(3,js-1); je3 = min(npy-2,je+2)
658 je1 = min(npy-3,je+1)
661 js1 = js-1; je3 = je+2
671 al(i,j) =
p1*(q(i,j-1)+q(i,j)) +
p2*(q(i,j-2)+q(i,j+1))
675 if ( .not. (nested .or. regional) .and. grid_type<3 )
then 678 al(i,0) =
c1*q(i,-2) +
c2*q(i,-1) +
c3*q(i,0)
679 al(i,1) = 0.5*(((2.*dya(i,0)+dya(i,-1))*q(i,0)-dya(i,0)*q(i,-1))/(dya(i,-1)+dya(i,0)) &
680 + ((2.*dya(i,1)+dya(i,2))*q(i,1)-dya(i,1)*q(i,2))/(dya(i,1)+dya(i,2)))
681 al(i,2) =
c3*q(i,1) +
c2*q(i,2) +
c1*q(i,3)
684 if( (je+1)==npy )
then 686 al(i,npy-1) =
c1*q(i,npy-3) +
c2*q(i,npy-2) +
c3*q(i,npy-1)
687 al(i,npy) = 0.5*(((2.*dya(i,npy-1)+dya(i,npy-2))*q(i,npy-1)-dya(i,npy-1)*q(i,npy-2))/(dya(i,npy-2)+dya(i,npy-1)) &
688 + ((2.*dya(i,npy)+dya(i,npy+1))*q(i,npy)-dya(i,npy)*q(i,npy+1))/(dya(i,npy)+dya(i,npy+1)))
689 al(i,npy+1) =
c3*q(i,npy) +
c2*q(i,npy+1) +
c1*q(i,npy+2)
697 al(i,j) = max(0., al(i,j))
705 bl(i,j) = al(i,j ) - q(i,j)
706 br(i,j) = al(i,j+1) - q(i,j)
707 b0(i,j) = bl(i,j) + br(i,j)
708 smt5(i,j) = abs(lim_fac*b0(i,j)) < abs(bl(i,j)-br(i,j))
714 if ( c(i,j) > 0. )
then 715 fx1(i) = (1.-c(i,j))*(br(i,j-1) - c(i,j)*b0(i,j-1))
718 fx1(i) = (1.+c(i,j))*(bl(i,j) + c(i,j)*b0(i,j))
721 if (smt5(i,j-1).or.smt5(i,j)) flux(i,j) = flux(i,j) + fx1(i)
725 elseif ( mord==2 )
then 734 flux(i,j) = qtmp + (1.-xt)*(al(i,j)-qtmp-xt*(al(i,j-1)+al(i,j)-(qtmp+qtmp)))
737 flux(i,j) = qtmp + (1.+xt)*(al(i,j)-qtmp+xt*(al(i,j)+al(i,j+1)-(qtmp+qtmp)))
742 elseif ( mord==3 )
then 746 bl(i,j) = al(i,j ) - q(i,j)
747 br(i,j) = al(i,j+1) - q(i,j)
748 b0(i,j) = bl(i,j) + br(i,j)
750 xt = abs(bl(i,j)-br(i,j))
752 smt6(i,j) = 3.*x0 < xt
759 hi5(i) = smt5(i,j-1) .and. smt5(i,j)
760 hi6(i) = smt6(i,j-1) .or. smt6(i,j)
763 if ( xt1(i) > 0. )
then 765 fx1(i) = br(i,j-1) - xt1(i)*b0(i,j-1)
766 elseif ( hi5(i) )
then 767 fx1(i) = sign(min(abs(bl(i,j-1)),abs(br(i,j-1))),br(i,j-1))
769 flux(i,j) = q(i,j-1) + (1.-xt1(i))*fx1(i)
772 fx1(i) = bl(i,j) + xt1(i)*b0(i,j)
773 elseif ( hi5(i) )
then 774 fx1(i) = sign(min(abs(bl(i,j)),abs(br(i,j))), bl(i,j))
776 flux(i,j) = q(i,j) + (1.+xt1(i))*fx1(i)
781 elseif ( mord==4 )
then 785 bl(i,j) = al(i,j ) - q(i,j)
786 br(i,j) = al(i,j+1) - q(i,j)
787 b0(i,j) = bl(i,j) + br(i,j)
789 xt = abs(bl(i,j)-br(i,j))
791 smt6(i,j) = 3.*x0 < xt
797 hi5(i) = smt5(i,j-1) .and. smt5(i,j)
798 hi6(i) = smt6(i,j-1) .or. smt6(i,j)
799 hi5(i) = hi5(i) .or. hi6(i)
803 if ( xt1(i) > 0. )
then 804 fx1(i) = (1.-xt1(i))*(br(i,j-1) - xt1(i)*b0(i,j-1))
807 fx1(i) = (1.+xt1(i))*(bl(i,j) + xt1(i)*b0(i,j))
810 if ( hi5(i) ) flux(i,j) = flux(i,j) + fx1(i)
818 bl(i,j) = al(i,j ) - q(i,j)
819 br(i,j) = al(i,j+1) - q(i,j)
820 b0(i,j) = bl(i,j) + br(i,j)
821 smt5(i,j) = bl(i,j)*br(i,j) < 0.
827 bl(i,j) = al(i,j ) - q(i,j)
828 br(i,j) = al(i,j+1) - q(i,j)
829 b0(i,j) = bl(i,j) + br(i,j)
830 smt5(i,j) = 3.*abs(b0(i,j)) < abs(bl(i,j)-br(i,j))
838 if ( c(i,j) > 0. )
then 839 fx1(i) = (1.-c(i,j))*(br(i,j-1) - c(i,j)*b0(i,j-1))
842 fx1(i) = (1.+c(i,j))*(bl(i,j) + c(i,j)*b0(i,j))
845 if (smt5(i,j-1).or.smt5(i,j)) flux(i,j) = flux(i,j) + fx1(i)
859 xt = 0.25*(q(i,j+1) - q(i,j-1))
860 dm(i,j) = sign(min(abs(xt), max(q(i,j-1), q(i,j), q(i,j+1)) - q(i,j), &
861 q(i,j) - min(q(i,j-1), q(i,j), q(i,j+1))), xt)
866 al(i,j) = 0.5*(q(i,j-1)+q(i,j)) +
r3*(dm(i,j-1) - dm(i,j))
874 bl(i,j) = -sign(min(abs(xt), abs(al(i,j)-q(i,j))), xt)
875 br(i,j) = sign(min(abs(xt), abs(al(i,j+1)-q(i,j))), xt)
878 elseif ( jord==11 )
then 882 bl(i,j) = -sign(min(abs(xt), abs(al(i,j)-q(i,j))), xt)
883 br(i,j) = sign(min(abs(xt), abs(al(i,j+1)-q(i,j))), xt)
889 dq(i,j) = 2.*(q(i,j+1) - q(i,j))
894 bl(i,j) = al(i,j ) - q(i,j)
895 br(i,j) = al(i,j+1) - q(i,j)
896 if ( abs(dm(i,j-1))+abs(dm(i,j))+abs(dm(i,j+1)) <
near_zero )
then 899 elseif( abs(3.*(bl(i,j)+br(i,j))) > abs(bl(i,j)-br(i,j)) )
then 901 lac_2 = pmp_2 - 0.75*dq(i,j-2)
902 br(i,j) = min(max(0.,pmp_2,lac_2), max(br(i,j), min(0.,pmp_2,lac_2)))
904 lac_1 = pmp_1 + 0.75*dq(i,j+1)
905 bl(i,j) = min(max(0.,pmp_1,lac_1), max(bl(i,j), min(0.,pmp_1,lac_1)))
910 if ( jord==9 .or. jord==13 )
then 913 call pert_ppm(ilast-ifirst+1, q(ifirst,j), bl(ifirst,j), br(ifirst,j), 0)
917 if (.not. (nested .or. regional) .and. grid_type<3)
then 920 bl(i,0) =
s14*dm(i,-1) +
s11*(q(i,-1)-q(i,0))
922 xt = 0.5*(((2.*dya(i,0)+dya(i,-1))*q(i,0)-dya(i,0)*q(i,-1))/(dya(i,-1)+dya(i,0)) &
923 + ((2.*dya(i,1)+dya(i,2))*q(i,1)-dya(i,1)*q(i,2))/(dya(i,1)+dya(i,2)))
925 xt = max(xt, min(q(i,-1),q(i,0),q(i,1),q(i,2)))
926 xt = min(xt, max(q(i,-1),q(i,0),q(i,1),q(i,2)))
928 br(i,0) = xt - q(i,0)
929 bl(i,1) = xt - q(i,1)
931 xt =
s15*q(i,1) +
s11*q(i,2) -
s14*dm(i,2)
932 br(i,1) = xt - q(i,1)
933 bl(i,2) = xt - q(i,2)
935 br(i,2) = al(i,3) - q(i,2)
937 call pert_ppm(3*(ilast-ifirst+1), q(ifirst,0), bl(ifirst,0), br(ifirst,0), 1)
939 if( (je+1)==npy )
then 941 bl(i,npy-2) = al(i,npy-2) - q(i,npy-2)
943 xt =
s15*q(i,npy-1) +
s11*q(i,npy-2) +
s14*dm(i,npy-2)
944 br(i,npy-2) = xt - q(i,npy-2)
945 bl(i,npy-1) = xt - q(i,npy-1)
947 xt = 0.5*(((2.*dya(i,npy-1)+dya(i,npy-2))*q(i,npy-1)-dya(i,npy-1)*q(i,npy-2))/(dya(i,npy-2)+dya(i,npy-1)) &
948 + ((2.*dya(i,npy)+dya(i,npy+1))*q(i,npy)-dya(i,npy)*q(i,npy+1))/(dya(i,npy)+dya(i,npy+1)))
950 xt = max(xt, min(q(i,npy-2),q(i,npy-1),q(i,npy),q(i,npy+1)))
951 xt = min(xt, max(q(i,npy-2),q(i,npy-1),q(i,npy),q(i,npy+1)))
953 br(i,npy-1) = xt - q(i,npy-1)
954 bl(i,npy ) = xt - q(i,npy)
956 br(i,npy) =
s11*(q(i,npy+1)-q(i,npy)) -
s14*dm(i,npy+1)
958 call pert_ppm(3*(ilast-ifirst+1), q(ifirst,npy-2), bl(ifirst,npy-2), br(ifirst,npy-2), 1)
967 flux(i,j) = q(i,j-1) + (1.-c(i,j))*(br(i,j-1)-c(i,j)*(bl(i,j-1)+br(i,j-1)))
969 flux(i,j) = q(i,j ) + (1.+c(i,j))*(bl(i,j )+c(i,j)*(bl(i,j)+br(i,j)))
977 subroutine mp_ghost_ew(im, jm, km, nq, ifirst, ilast, jfirst, jlast, &
978 kfirst, klast, ng_w, ng_e, ng_s, ng_n, q_ghst, q)
981 integer,
intent(in):: im, jm, km, nq
982 integer,
intent(in):: ifirst, ilast
983 integer,
intent(in):: jfirst, jlast
984 integer,
intent(in):: kfirst, klast
985 integer,
intent(in):: ng_e
986 integer,
intent(in):: ng_w
987 integer,
intent(in):: ng_s
988 integer,
intent(in):: ng_n
989 real,
intent(inout):: q_ghst(ifirst-ng_w:ilast+ng_e,jfirst-ng_s:jlast+ng_n,kfirst:klast,nq)
990 real,
optional,
intent(in):: q(ifirst:ilast,jfirst:jlast,kfirst:klast,nq)
1004 if (
present(q))
then 1005 q_ghst(ifirst:ilast,jfirst:jlast,kfirst:klast,1:nq) = &
1006 q(ifirst:ilast,jfirst:jlast,kfirst:klast,1:nq)
1012 do j=jfirst-ng_s,jlast+ng_n
1014 q_ghst(ifirst-i,j,k,n) = q_ghst(ilast-i+1,j,k,n)
1017 q_ghst(ilast+i,j,k,n) = q_ghst(ifirst+i-1,j,k,n)
1027 subroutine pert_ppm(im, a0, al, ar, iv)
1028 integer,
intent(in):: im
1029 integer,
intent(in):: iv
1030 real,
intent(in) :: a0(im)
1031 real,
intent(inout):: al(im), ar(im)
1033 real a4, da1, da2, a6da, fmin
1035 real,
parameter:: r12 = 1./12.
1044 if ( a0(i) <= 0. )
then 1048 a4 = -3.*(ar(i) + al(i))
1050 if( abs(da1) < -a4 )
then 1051 fmin = a0(i) + 0.25/a4*da1**2 + a4*r12
1052 if( fmin < 0. )
then 1053 if( ar(i)>0. .and. al(i)>0. )
then 1056 elseif( da1 > 0. )
then 1068 if ( al(i)*ar(i) < 0. )
then 1071 a6da = 3.*(al(i)+ar(i))*da1
1073 if( a6da < -da2 )
then 1075 elseif( a6da > da2 )
then 1089 subroutine deln_flux(nord,is,ie,js,je, npx, npy, damp, q, fx, fy, gridstruct,regional, bd, mass )
1097 type(fv_grid_bounds_type),
intent(IN) :: bd
1098 integer,
intent(in):: nord
1099 integer,
intent(in):: is,ie,js,je, npx, npy
1100 real,
intent(in):: damp
1101 real,
intent(in):: q(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng)
1102 type(fv_grid_type),
intent(IN),
target :: gridstruct
1103 logical,
intent(in):: regional
1104 real,
optional,
intent(in):: mass(bd%isd:bd%ied, bd%jsd:bd%jed)
1106 real,
intent(inout):: fx(bd%is:bd%ie+1,bd%js:bd%je), fy(bd%is:bd%ie,bd%js:bd%je+1)
1108 real fx2(bd%isd:bd%ied+1,bd%jsd:bd%jed), fy2(bd%isd:bd%ied,bd%jsd:bd%jed+1)
1109 real d2(bd%isd:bd%ied,bd%jsd:bd%jed)
1111 integer i,j, n, nt, i1, i2, j1, j2
1114 real,
pointer,
dimension(:,:) :: dx, dy, rdxc, rdyc
1115 real,
pointer,
dimension(:,:,:) :: sin_sg
1118 rdxc => gridstruct%rdxc
1119 rdyc => gridstruct%rdyc
1120 sin_sg => gridstruct%sin_sg
1123 i1 = is-1-nord; i2 = ie+1+nord
1124 j1 = js-1-nord; j2 = je+1+nord
1126 if ( .not.
present(mass) )
then 1129 d2(i,j) = damp*q(i,j)
1140 if( nord>0 .and. (.not. (regional)))
call copy_corners(d2, npx, npy, 1, gridstruct%nested, bd, &
1141 gridstruct%sw_corner, gridstruct%se_corner, gridstruct%nw_corner, gridstruct%ne_corner)
1143 do j=js-nord,je+nord
1144 do i=is-nord,ie+nord+1
1146 fx2(i,j) = 0.5*(sin_sg(i-1,j,3)+sin_sg(i,j,1))*dy(i,j)*(d2(i-1,j)-d2(i,j))*rdxc(i,j)
1148 fx2(i,j) = gridstruct%del6_v(i,j)*(d2(i-1,j)-d2(i,j))
1153 if( nord>0 .and. (.not. (regional)))
call copy_corners(d2, npx, npy, 2, gridstruct%nested, bd, &
1154 gridstruct%sw_corner, gridstruct%se_corner, gridstruct%nw_corner, gridstruct%ne_corner)
1155 do j=js-nord,je+nord+1
1156 do i=is-nord,ie+nord
1158 fy2(i,j) = 0.5*(sin_sg(i,j-1,4)+sin_sg(i,j,2))*dx(i,j)*(d2(i,j-1)-d2(i,j))*rdyc(i,j)
1160 fy2(i,j) = gridstruct%del6_u(i,j)*(d2(i,j-1)-d2(i,j))
1175 do j=js-nt-1,je+nt+1
1176 do i=is-nt-1,ie+nt+1
1177 d2(i,j) = (fx2(i,j)-fx2(i+1,j)+fy2(i,j)-fy2(i,j+1))*gridstruct%rarea(i,j)
1181 if (.not.(regional))
call copy_corners(d2, npx, npy, 1, gridstruct%nested, bd, &
1182 gridstruct%sw_corner, gridstruct%se_corner, gridstruct%nw_corner, gridstruct%ne_corner)
1186 fx2(i,j) = 0.5*(sin_sg(i-1,j,3)+sin_sg(i,j,1))*dy(i,j)*(d2(i,j)-d2(i-1,j))*rdxc(i,j)
1188 fx2(i,j) = gridstruct%del6_v(i,j)*(d2(i,j)-d2(i-1,j))
1193 if (.not.(regional))
call copy_corners(d2, npx, npy, 2, gridstruct%nested, bd, &
1194 gridstruct%sw_corner, gridstruct%se_corner, gridstruct%nw_corner, gridstruct%ne_corner)
1198 fy2(i,j) = 0.5*(sin_sg(i,j-1,4)+sin_sg(i,j,2))*dx(i,j)*(d2(i,j)-d2(i,j-1))*rdyc(i,j)
1200 fy2(i,j) = gridstruct%del6_u(i,j)*(d2(i,j)-d2(i,j-1))
1212 if (
present(mass) )
then 1217 fx(i,j) = fx(i,j) + damp2*(mass(i-1,j)+mass(i,j))*fx2(i,j)
1222 fy(i,j) = fy(i,j) + damp2*(mass(i,j-1)+mass(i,j))*fy2(i,j)
1228 fx(i,j) = fx(i,j) + fx2(i,j)
1233 fy(i,j) = fy(i,j) + fy2(i,j)
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
The type 'fv_grid_type' is made up of grid-dependent information from fv_grid_tools and fv_grid_utils...
real, parameter ppm_limiter
The module 'tp_core' is a collection of routines to support FV transport.
subroutine, public copy_corners(q, npx, npy, dir, nested, bd, sw_corner, se_corner, nw_corner, ne_corner)
real, parameter near_zero
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
real, parameter, public big_number
subroutine, public pert_ppm(im, a0, al, ar, iv)
subroutine yppm(flux, q, c, jord, ifirst, ilast, isd, ied, js, je, jsd, jed, npx, npy, dya, nested, grid_type, lim_fac, regional)
real, parameter ppm_fac
nonlinear scheme limiter: between 1 and 2
subroutine xppm(flux, q, c, iord, is, ie, isd, ied, jfirst, jlast, jsd, jed, npx, npy, dxa, nested, grid_type, lim_fac, regional)
integer, parameter, public ng
The module 'fv_grid_utils' contains routines for setting up and computing grid-related quantities...
subroutine, public fv_tp_2d(q, crx, cry, npx, npy, hord, fx, fy, xfx, yfx, gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx, mfy, mass, nord, damp_c)
The subroutine 'fv_tp_2d' contains the FV advection scheme .
subroutine mp_ghost_ew(im, jm, km, nq, ifirst, ilast, jfirst, jlast, kfirst, klast, ng_w, ng_e, ng_s, ng_n, q_ghst, q)
subroutine deln_flux(nord, is, ie, js, je, npx, npy, damp, q, fx, fy, gridstruct, regional, bd, mass)