Dataset S1. Second-Order Moments (SOM) advection FORTRAN code C---(QVECT.f)----SOM (Prather, 1986) advection core subroutine c--- new INTRINSIC/VECTOR formulation of CTM advection core no conditionals c--- >>>NB<<< this abstract code is correct and easier to follow, BUT c--- it runs much more slowly than the same code with do-loops on typical c--- computers (e.g., IBM p-series 8- and 32-proc boards) c C---from UCIrvine (p-code 5.x, prather 8/2003, minor rev 4/2007) c--- Note that all advectionis piped flow and must be cyclical in new code c--- This is easily achieved for E-W, for N-S with over-the-pole flow, c--- and for vertical by concatenating many pipes with QU=0 between. c---subroutines: QVECT c----------------------------------------------------------------------- subroutine QVECT (QM,QU,NQ,LIMTR, & QTT,QXT,QXX,QXY,QXZ,QYT,QYY,QZT,QZZ,QYZ) c----------------------------------------------------------------------- C Computes SOM advection for a single 'piped flow' of tracer in the X-direction C-----that is cyclic = connect QM(NQ) to QM(1), assume QU cyclic: QU(NQ+1)=QU(1) C----- =QU(1)=> BOX[1] =QU(2)=> ... =QU(NQ)=> BOX[NQ] =QU(1)=> BOX[1] ... c----------------------------------------------------------------------- c Different flux limiters are allowed: C----LIMTR=0 no limits, neg. tracer allowed - OK C----LIMTR=1 pos. definite, simple SOM fix (per Prather 1986):adjust QXT & QXX C----LIMTR=2 pos. definite & monotonic (ca. 1994-95): adjust QXT,QXX,QXY,QXZ C----LIMTR=3 min/max limiters (revised 2003) c----------------------------------------------------------------------- c Flux diagnostic, Q0F, tracer flux from [I-1] to [I] is no longer computed c----------------------------------------------------------------------- implicit none integer, intent(in) :: NQ ! no. of boxes integer, intent(in) :: LIMTR ! flux limiter parameter real*8, intent(in) :: QU(NQ+1) ! air mass to be moved [I-1]->[I] in adv step real*8, intent(inout):: QM(NQ+1) ! air mass in box at start/finish real*8, intent(inout):: QTT(NQ+1) ! tracer mass in box [I] at start/finish real*8, intent(inout):: QXT(NQ+1),QXX(NQ+1) ! 1st & 2nd moments of tracer in X real*8, intent(inout):: QXY(NQ+1),QXZ(NQ+1) ! cross-moments with X component real*8, intent(inout):: QYT(NQ+1),QYY(NQ+1) ! other SOMs without X real*8, intent(inout):: QZT(NQ+1),QZZ(NQ+1) real*8, intent(inout):: QYZ(NQ+1) C local parameters integer J,NJ real*8, allocatable, dimension(:):: FTT,FXX,FYY,FZZ,FXT,FYT,FZT real*8, allocatable, dimension(:):: FXY,FXZ,FYZ, EPS real*8, allocatable, dimension(:):: A0M, A0TT,A0XT,A0XX,A0YT,A0YY real*8, allocatable, dimension(:):: A0ZT,A0ZZ,A0XY,A0XZ,A0YZ real*8, allocatable, dimension(:):: A1M, A1TT,A1XT,A1XX,A1YT,A1YY real*8, allocatable, dimension(:):: A1ZT,A1ZZ,A1XY,A1XZ,A1YZ real*8, allocatable, dimension(:):: B0M, B0TT,B0XT,B0XX,B0YT,B0YY real*8, allocatable, dimension(:):: B0ZT,B0ZZ,B0XY,B0XZ,B0YZ real*8, allocatable, dimension(:):: B1M, B1TT,B1XT,B1XX,B1YT,B1YY real*8, allocatable, dimension(:):: B1ZT,B1ZZ,B1XY,B1XZ,B1YZ real*8, allocatable, dimension(:):: UAB logical, allocatable, dimension(:):: LGC_AB logical, allocatable, dimension(:):: LGC_L2, LGC_LX real*8, allocatable, dimension(:):: A3XT,A3XX,A4XT,A4XX real*8, allocatable, dimension(:):: ABSCAL real*8, allocatable, dimension(:):: ABXMX,ABXMN,ABYMX real*8, allocatable, dimension(:):: ABYMN,ABZMX,ABZMN real*8, allocatable, dimension(:):: XXAMX,XYAMX,XZAMX real*8, allocatable, dimension(:):: XXAMN,XYAMN,XZAMN real*8, allocatable, dimension(:):: XXBMX,XYBMX,XZBMX real*8, allocatable, dimension(:):: XXBMN,XYBMN,XZBMN real*8 C033,C034,C149,C000, CXYZ c---C000 = effective zero mass, needed for min to avoid zero divide C000 = 1.d-30 c---CXYZ used to limit cross moments, could be = 1.0 or 1.5 CXYZ = 1.0d0 C033 = 1.0d0/3.0d0 C149 = 1.49999d0 C034 = 0.33334d0 c----------------------------------------------------------------------- C---REQUIRE that NQ be even (easy if adding non-cyclic columns, but not if true C--- odd-numbered, cyclic solution) - cannot do separartion otherwise. if (mod(NQ,2) .ne. 0) STOP C--Assume that always cyclical and add box #1 as #NQ+1. C-- this is easy even if LCYCLE = .F., just zero connecting velocity. C-- if we are running a window-on-the-world, C--- then create pseudo boxes 1 and NQ+1 ahead of time C--Think of the sequence of boxes (begin with I=1) as: A B A B A B A B A B A B (A) C--Now first prepare a column of 'advects' A => B C-- if U>0: A1 & B1 have all moments stored C-- if U<0: switch A & B, reverse sign of X moments c-- this must be done in the calling program: QU(NQ+1) = QU(1) NJ = NQ/2 allocate (A0M(NJ),A0TT(NJ),A0XT(NJ),A0XX(NJ),A0YT(NJ),A0YY(NJ)) allocate (A0ZT(NJ),A0ZZ(NJ),A0XY(NJ),A0YZ(NJ),A0XZ(NJ)) allocate (A1M(NJ),A1TT(NJ),A1XT(NJ),A1XX(NJ),A1YT(NJ),A1YY(NJ)) allocate (A1ZT(NJ),A1ZZ(NJ),A1XY(NJ),A1YZ(NJ),A1XZ(NJ)) allocate (B0M(NJ),B0TT(NJ),B0XT(NJ),B0XX(NJ),B0YT(NJ),B0YY(NJ)) allocate (B0ZT(NJ),B0ZZ(NJ),B0XY(NJ),B0YZ(NJ),B0XZ(NJ)) allocate (B1M(NJ),B1TT(NJ),B1XT(NJ),B1XX(NJ),B1YT(NJ),B1YY(NJ)) allocate (B1ZT(NJ),B1ZZ(NJ),B1XY(NJ),B1YZ(NJ),B1XZ(NJ)) allocate (FTT(NJ),FXX(NJ),FYY(NJ),FZZ(NJ),FXT(NJ),FYT(NJ),FZT(NJ)) allocate (FXY(NJ),FXZ(NJ),FYZ(NJ), EPS(NJ)) allocate (UAB(NJ), LGC_AB(NJ)) allocate (LGC_L2(NJ),LGC_LX(NJ)) allocate (A3XT(NJ),A3XX(NJ),A4XT(NJ),A4XX(NJ)) allocate (ABSCAL(NJ)) allocate (ABXMX(NJ),ABXMN(NJ),ABYMX(NJ)) allocate (ABYMN(NJ),ABZMX(NJ),ABZMN(NJ)) allocate (XXAMX(NJ),XYAMX(NJ),XZAMX(NJ)) allocate (XXAMN(NJ),XYAMN(NJ),XZAMN(NJ)) allocate (XXBMX(NJ),XYBMX(NJ),XZBMX(NJ)) allocate (XXBMN(NJ),XYBMN(NJ),XZBMN(NJ)) C---#1 NOW load the EVEN pairs: do J=1,NJ LGC_AB(J) = QU(J+J).ge.0.d0 UAB(J) = abs(QU(J+J)) enddo do J=1,NJ C--- load the EVEN pair of Q--'s into A--'s and B--'s A1M(J) = QM(J+J-1) A1TT(J) = QTT(J+J-1) A1XX(J) = QXX(J+J-1) A1YT(J) = QYT(J+J-1) A1YY(J) = QYY(J+J-1) A1ZT(J) = QZT(J+J-1) A1ZZ(J) = QZZ(J+J-1) A1YZ(J) = QYZ(J+J-1) A1XT(J) = QXT(J+J-1) A1XY(J) = QXY(J+J-1) A1XZ(J) = QXZ(J+J-1) C--- same for the B's B1M(J) = QM(J+J) B1TT(J) = QTT(J+J) B1XX(J) = QXX(J+J) B1YT(J) = QYT(J+J) B1YY(J) = QYY(J+J) B1ZT(J) = QZT(J+J) B1ZZ(J) = QZZ(J+J) B1YZ(J) = QYZ(J+J) B1XT(J) = QXT(J+J) B1XY(J) = QXY(J+J) B1XZ(J) = QXZ(J+J) enddo C---Now create vector set of A0 & B0 quantities reflecting sign of QU C--- NB first-order X-moments switch sign for reverse advection A0M = merge (A1M, B1M, LGC_AB) A0TT = merge (A1TT, B1TT, LGC_AB) A0XX = merge (A1XX, B1XX, LGC_AB) A0YT = merge (A1YT, B1YT, LGC_AB) A0YY = merge (A1YY, B1YY, LGC_AB) A0ZT = merge (A1ZT, B1ZT, LGC_AB) A0ZZ = merge (A1ZZ, B1ZZ, LGC_AB) A0YZ = merge (A1YZ, B1YZ, LGC_AB) A0XT = merge (A1XT,-B1XT, LGC_AB) A0XY = merge (A1XY,-B1XY, LGC_AB) A0XZ = merge (A1XZ,-B1XZ, LGC_AB) C--- B0M = merge (B1M, A1M, LGC_AB) B0TT = merge (B1TT, A1TT, LGC_AB) B0XX = merge (B1XX, A1XX, LGC_AB) B0YT = merge (B1YT, A1YT, LGC_AB) B0YY = merge (B1YY, A1YY, LGC_AB) B0ZT = merge (B1ZT, A1ZT, LGC_AB) B0ZZ = merge (B1ZZ, A1ZZ, LGC_AB) B0YZ = merge (B1YZ, A1YZ, LGC_AB) B0XT = merge (B1XT,-A1XT, LGC_AB) B0XY = merge (B1XY,-A1XY, LGC_AB) B0XZ = merge (B1XZ,-A1XZ, LGC_AB) C---LIMTR = 1 C---positive definite, reset Sxx to allow maximum Sx (1.5*So) if (LIMTR.eq.1) then A0XT = min(+C149*A0TT, max(-C149*A0TT, A0XT)) A0XX = min((A0TT+A0TT-C034*abs(A0XT)), max(abs(A0XT)-A0TT, A0XX)) A0XY = min(+CXYZ*A0TT, max(-CXYZ*A0TT,A0XY)) A0XZ = min(+CXYZ*A0TT, max(-CXYZ*A0TT,A0XZ)) C--------DO Not impose limits on downwind B since it is not transported C B0XT = min(+C149*B0TT, max(-C149*B0TT,B0XT)) C B0XX = min((B0TT+B0TT-C034*abs(B0XT)), max(abs(B0XT)-B0TT, B0XX)) C B0XY = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XY)) C B0XZ = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XZ)) endif C---LIMTR = 2 C-------- positive and monotonic in box (-dim only) if (LIMTR.ge.2) then LGC_L2 = A0XX .gt. 0.d0 LGC_LX = A0XT .gt. 0.D0 C3---do not allow reversal of curvature: A0XX(J) > 0 A3XX = min(A0XX, abs(A0XT)*C033, 0.5d0*A0TT) A3XT = min(A0TT+A3XX, abs(A0XT)) A3XT = merge(A3XT, -A3XT, LGC_LX) C4---A0XX(J) < 0 = curved down at ends, allow if A0XT(J) < A0TT(J) A4XT = min(+A0TT, max(-A0TT, A0XT)) A4XX = max(A0XX, -A0TT+abs(A4XT), -abs(A4XT)*C033) A0XT = merge (A3XT, A4XT, LGC_L2) A0XX = merge (A3XX, A4XX, LGC_L2) A0XY = min(+CXYZ*A0TT,max(-CXYZ*A0TT,A0XY)) A0XZ = min(+CXYZ*A0TT,max(-CXYZ*A0TT,A0XZ)) C-----IF LIMTR=3, then force monotonic on downwind B also if (LIMTR.eq.3) then LGC_L2 = B0XX .gt. 0.d0 LGC_LX = B0XT .gt. 0.D0 A3XX = min(B0XX, abs(B0XT)*C033, 0.5d0*B0TT) A3XT = min(B0TT+A3XX, abs(B0XT)) A3XT = merge(A3XT, -A3XT, LGC_LX) A4XT = min(+B0TT, max(-B0TT, B0XT)) A4XX = max(B0XX, -B0TT+abs(A4XT), -abs(A4XT)*C033) B0XT = merge (A3XT, A4XT, LGC_L2) B0XX = merge (A3XX, A4XX, LGC_L2) B0XY = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XY)) B0XZ = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XZ)) endif endif C---LIMTR = 3 if (LIMTR.eq.3) then c---min/max in A-- & B-- boxes occurs at end points if monotonic XXAMX = (A0TT + abs(A0XT) + A0XX)/A0M XYAMX = (A0TT + abs(A0XT) + A0XX + abs(A0XY))/A0M XZAMX = (A0TT + abs(A0XT) + A0XX + abs(A0XZ))/A0M XXAMN = (A0TT - abs(A0XT) + A0XX)/A0M XYAMN = (A0TT - abs(A0XT) + A0XX - abs(A0XY))/A0M XZAMN = (A0TT - abs(A0XT) + A0XX - abs(A0XZ))/A0M XXBMX = (B0TT + abs(B0XT) + B0XX)/B0M XYBMX = (B0TT + abs(B0XT) + B0XX + abs(B0XY))/B0M XZBMX = (B0TT + abs(B0XT) + B0XX + abs(B0XZ))/B0M XXBMN = (B0TT - abs(B0XT) + B0XX)/B0M XYBMN = (B0TT - abs(B0XT) + B0XX - abs(B0XY))/B0M XZBMN = (B0TT - abs(B0XT) + B0XX - abs(B0XZ))/B0M ABXMX = max(XXAMX,XXBMX) ABXMN = min(XXAMN,XXBMN) ABYMX = max(XYAMX,XYBMX) ABYMN = min(XYAMN,XYBMN) ABZMX = max(XZAMX,XZBMX) ABZMN = min(XZAMN,XZBMN) endif C---Advect from A0 ==[UAB]==> to B0 for even pairs [1:2] [3:4] [5:6] EPS = UAB/A0M c------Seperate left A0 box to 2 parts: F's = section tranferred to B0 box FYY = EPS*A0YY FZZ = EPS*A0ZZ FYZ = EPS*A0YZ FYT = EPS*(A0YT + (1.d0-EPS)*A0XY) FZT = EPS*(A0ZT + (1.d0-EPS)*A0XZ) FXY = A0XY*EPS**2 FXZ = A0XZ*EPS**2 FXX = A0XX*EPS**3 FXT = (A0XT + 3.d0*(1.d0-EPS)*A0XX)*EPS**2 FTT = EPS*(A0TT + (1.d0-EPS)*A0XT & + (1.d0-EPS)*(1.d0-EPS-EPS)*A0XX) C---ORDER is important for A0XT! A0YY = A0YY - FYY A0ZZ = A0ZZ - FZZ A0YZ = A0YZ - FYZ A0YT = A0YT - FYT A0ZT = A0ZT - FZT A0XY = A0XY*(1.d0-EPS)**2 A0XZ = A0XZ*(1.d0-EPS)**2 A0XT = (A0XT - 3.d0*EPS*A0XX)*(1.d0-EPS)**2 A0XX = A0XX*(1.d0-EPS)**3 A0TT = A0TT - FTT A0M = A0M - UAB B0M = B0M + UAB EPS = UAB/B0M C combine transferred moments (F's) into box B0 - ORDER important ! B0ZZ = B0ZZ + FZZ B0YY = B0YY + FYY B0YZ = B0YZ + FYZ B0XY = EPS*FXY + (1.d0-EPS)*B0XY & +3.d0*(EPS*B0YT - (1.d0-EPS)*FYT) B0YT = B0YT + FYT B0XZ = EPS*FXZ + (1.d0-EPS)*B0XZ & +3.d0*(EPS*B0ZT - (1.d0-EPS)*FZT) B0ZT = B0ZT + FZT B0XX = FXX*EPS**2 + B0XX*(1.d0-EPS)**2 & +5.d0*( EPS*(1.d0-EPS)*(B0XT - FXT) & +(1.d0-EPS-EPS)*((1.d0-EPS)*FTT-EPS*B0TT)) B0XT = EPS*FXT + (1.d0-EPS)*B0XT & +3.d0*(EPS*B0TT- (1.d0-EPS)*FTT) B0TT = B0TT + FTT C---LIMTR = 3 if (LIMTR.eq.3) then c------------force daughter box B to be monotonic(L=2) first LGC_L2 = B0XX .gt. 0.d0 LGC_LX = B0XT .gt. 0.D0 A3XX = min(B0XX, abs(B0XT)*C033, 0.5d0*B0TT) A3XT = min(B0TT+A3XX, abs(B0XT)) A3XT = merge(A3XT, -A3XT, LGC_LX) A4XT = min(+B0TT, max(-B0TT, B0XT)) A4XX = max(B0XX, -B0TT+abs(A4XT), -abs(A4XT)*C033) B0XT = merge (A3XT, A4XT, LGC_L2) B0XX = merge (A3XX, A4XX, LGC_L2) B0XY = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XY)) B0XZ = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XZ)) c------------limit daughter box B by min/max of parent A & B boxes C***NB have had to assume a near-zero level to avoid zero divide c---rescale for X-max ABSCAL = max(0.d0, min(1.d0, & (ABXMX*B0M-B0TT)/max(abs(B0XT)+B0XX, C000))) B0XT = ABSCAL*B0XT B0XX = ABSCAL*B0XX c---rescale for X-min ABSCAL = max(0.d0, min(1.d0, & (B0TT-ABXMN*B0M)/max(abs(B0XT)-B0XX, C000))) B0XT = ABSCAL*B0XT B0XX = ABSCAL*B0XX c---rescale for XY moment ABSCAL = max(0.d0, min(1.d0, & (ABYMX*B0M-B0TT)/max(abs(B0XT)+B0XX+abs(B0XY), C000))) B0XY = ABSCAL*B0XY ABSCAL = max(0.d0, min(1.d0, & (B0TT-ABYMN*B0M)/max(abs(B0XT)-B0XX+abs(B0XY), C000))) B0XY = ABSCAL*B0XY c---rescale for XZ moment ABSCAL = max(0.d0, min(1.d0, & (ABZMX*B0M-B0TT)/max(abs(B0XT)+B0XX+abs(B0XZ), C000))) B0XZ = ABSCAL*B0XZ ABSCAL = max(0.d0, min(1.d0, & (B0TT-ABZMN*B0M)/max(abs(B0XT)-B0XX+abs(B0XZ), C000))) B0XZ = ABSCAL*B0XZ endif C---Now re-mask back onto A1-- AND B1-- arrays A1M = merge (A0M, B0M, LGC_AB) A1TT = merge (A0TT, B0TT, LGC_AB) A1XX = merge (A0XX, B0XX, LGC_AB) A1YT = merge (A0YT, B0YT, LGC_AB) A1YY = merge (A0YY, B0YY, LGC_AB) A1ZT = merge (A0ZT, B0ZT, LGC_AB) A1ZZ = merge (A0ZZ, B0ZZ, LGC_AB) A1YZ = merge (A0YZ, B0YZ, LGC_AB) A1XT = merge (A0XT,-B0XT, LGC_AB) A1XY = merge (A0XY,-B0XY, LGC_AB) A1XZ = merge (A0XZ,-B0XZ, LGC_AB) C--- B1M = merge (B0M, A0M, LGC_AB) B1TT = merge (B0TT, A0TT, LGC_AB) B1XX = merge (B0XX, A0XX, LGC_AB) B1YT = merge (B0YT, A0YT, LGC_AB) B1YY = merge (B0YY, A0YY, LGC_AB) B1ZT = merge (B0ZT, A0ZT, LGC_AB) B1ZZ = merge (B0ZZ, A0ZZ, LGC_AB) B1YZ = merge (B0YZ, A0YZ, LGC_AB) B1XT = merge (B0XT,-A0XT, LGC_AB) B1XY = merge (B0XY,-A0XY, LGC_AB) B1XZ = merge (B0XZ,-A0XZ, LGC_AB) C---#2 NOW we need to repeat the whole process with ODD pairs: C--- this wraps box(NQ) to box(1), need to have QU(NQ+1) defined do J=1,NJ LGC_AB(J) = QU(J+J+1).ge.0.d0 UAB(J) = abs(QU(J+J+1)) enddo C---Now shift A,B down by one, B(1)=>A(1), A(2)=>B(1),...B(N)=>A(N), A(1)=>B(N) A0M = B1M ! temporary store A0TT = B1TT A0XT = B1XT A0XX = B1XX A0YT = B1YT A0YY = B1YY A0ZT = B1ZT A0ZZ = B1ZZ A0XY = B1XY A0XZ = B1XZ A0YZ = B1YZ B1M = cshift(A1M, 1) B1TT = cshift(A1TT, 1) B1XT = cshift(A1XT, 1) B1XX = cshift(A1XX, 1) B1YT = cshift(A1YT, 1) B1YY = cshift(A1YY, 1) B1ZT = cshift(A1ZT, 1) B1ZZ = cshift(A1ZZ, 1) B1XY = cshift(A1XY, 1) B1XZ = cshift(A1XZ, 1) B1YZ = cshift(A1YZ, 1) A1M = A0M A1TT = A0TT A1XT = A0XT A1XX = A0XX A1YT = A0YT A1YY = A0YY A1ZT = A0ZT A1ZZ = A0ZZ A1XY = A0XY A1XZ = A0XZ A1YZ = A0YZ C---Now create vector set of A0 & B0 quantities reflecting sign of QU C--- NB first-order X-moments switch sign for reverse advection A0M = merge (A1M, B1M, LGC_AB) A0TT = merge (A1TT, B1TT, LGC_AB) A0XX = merge (A1XX, B1XX, LGC_AB) A0YT = merge (A1YT, B1YT, LGC_AB) A0YY = merge (A1YY, B1YY, LGC_AB) A0ZT = merge (A1ZT, B1ZT, LGC_AB) A0ZZ = merge (A1ZZ, B1ZZ, LGC_AB) A0YZ = merge (A1YZ, B1YZ, LGC_AB) A0XT = merge (A1XT,-B1XT, LGC_AB) A0XY = merge (A1XY,-B1XY, LGC_AB) A0XZ = merge (A1XZ,-B1XZ, LGC_AB) C--- B0M = merge (B1M, A1M, LGC_AB) B0TT = merge (B1TT, A1TT, LGC_AB) B0XX = merge (B1XX, A1XX, LGC_AB) B0YT = merge (B1YT, A1YT, LGC_AB) B0YY = merge (B1YY, A1YY, LGC_AB) B0ZT = merge (B1ZT, A1ZT, LGC_AB) B0ZZ = merge (B1ZZ, A1ZZ, LGC_AB) B0YZ = merge (B1YZ, A1YZ, LGC_AB) B0XT = merge (B1XT,-A1XT, LGC_AB) B0XY = merge (B1XY,-A1XY, LGC_AB) B0XZ = merge (B1XZ,-A1XZ, LGC_AB) C---LIMTR = 1 C---positive definite, reset Sxx to allow maximum Sx (1.5*So) if (LIMTR.eq.1) then A0XT = min(+C149*A0TT, max(-C149*A0TT, A0XT)) A0XX = min((A0TT+A0TT-C034*abs(A0XT)), max(abs(A0XT)-A0TT, A0XX)) A0XY = min(+CXYZ*A0TT, max(-CXYZ*A0TT,A0XY)) A0XZ = min(+CXYZ*A0TT, max(-CXYZ*A0TT,A0XZ)) C--------DO Not impose limits on downwind B since it is not transported C B0XT = min(+C149*B0TT, max(-C149*B0TT,B0XT)) C B0XX = min((B0TT+B0TT-C034*abs(B0XT)), max(abs(B0XT)-B0TT, B0XX)) C B0XY = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XY)) C B0XZ = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XZ)) endif C---LIMTR = 2 C-------- positive and monotonic in box (-dim only) if (LIMTR.ge.2) then LGC_L2 = A0XX .gt. 0.d0 LGC_LX = A0XT .gt. 0.D0 C3---do not allow reversal of curvature: A0XX(J) > 0 A3XX = min(A0XX, abs(A0XT)*C033, 0.5d0*A0TT) A3XT = min(A0TT+A3XX, abs(A0XT)) A3XT = merge(A3XT, -A3XT, LGC_LX) C4---A0XX(J) < 0 = curved down at ends, allow if A0XT(J) < A0TT(J) A4XT = min(+A0TT, max(-A0TT, A0XT)) A4XX = max(A0XX, -A0TT+abs(A4XT), -abs(A4XT)*C033) A0XT = merge (A3XT, A4XT, LGC_L2) A0XX = merge (A3XX, A4XX, LGC_L2) A0XY = min(+CXYZ*A0TT,max(-CXYZ*A0TT,A0XY)) A0XZ = min(+CXYZ*A0TT,max(-CXYZ*A0TT,A0XZ)) C-----IF LIMTR=3, then force monotonic on downwind B also if (LIMTR.eq.3) then LGC_L2 = B0XX .gt. 0.d0 LGC_LX = B0XT .gt. 0.D0 A3XX = min(B0XX, abs(B0XT)*C033, 0.5d0*B0TT) A3XT = min(B0TT+A3XX, abs(B0XT)) A3XT = merge(A3XT, -A3XT, LGC_LX) A4XT = min(+B0TT, max(-B0TT, B0XT)) A4XX = max(B0XX, -B0TT+abs(A4XT), -abs(A4XT)*C033) B0XT = merge (A3XT, A4XT, LGC_L2) B0XX = merge (A3XX, A4XX, LGC_L2) B0XY = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XY)) B0XZ = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XZ)) endif endif C---LIMTR = 3 if (LIMTR.eq.3) then c---min/max in A-- & B-- boxes occurs at end points if monotonic XXAMX = (A0TT + abs(A0XT) + A0XX)/A0M XYAMX = (A0TT + abs(A0XT) + A0XX + abs(A0XY))/A0M XZAMX = (A0TT + abs(A0XT) + A0XX + abs(A0XZ))/A0M XXAMN = (A0TT - abs(A0XT) + A0XX)/A0M XYAMN = (A0TT - abs(A0XT) + A0XX - abs(A0XY))/A0M XZAMN = (A0TT - abs(A0XT) + A0XX - abs(A0XZ))/A0M XXBMX = (B0TT + abs(B0XT) + B0XX)/B0M XYBMX = (B0TT + abs(B0XT) + B0XX + abs(B0XY))/B0M XZBMX = (B0TT + abs(B0XT) + B0XX + abs(B0XZ))/B0M XXBMN = (B0TT - abs(B0XT) + B0XX)/B0M XYBMN = (B0TT - abs(B0XT) + B0XX - abs(B0XY))/B0M XZBMN = (B0TT - abs(B0XT) + B0XX - abs(B0XZ))/B0M ABXMX = max(XXAMX,XXBMX) ABXMN = min(XXAMN,XXBMN) ABYMX = max(XYAMX,XYBMX) ABYMN = min(XYAMN,XYBMN) ABZMX = max(XZAMX,XZBMX) ABZMN = min(XZAMN,XZBMN) endif C---Advect from A0 ==[UAB]==> to B0 for odd pairs [2:3] [4:5] [6:7] [N:1] EPS = UAB/A0M c------Seperate left A0 box to 2 parts: F's = section tranferred to B0 box FYY = EPS*A0YY FZZ = EPS*A0ZZ FYZ = EPS*A0YZ FYT = EPS*(A0YT + (1.d0-EPS)*A0XY) FZT = EPS*(A0ZT + (1.d0-EPS)*A0XZ) FXY = A0XY*EPS**2 FXZ = A0XZ*EPS**2 FXX = A0XX*EPS**3 FXT = (A0XT + 3.d0*(1.d0-EPS)*A0XX)*EPS**2 FTT = EPS*(A0TT + (1.d0-EPS)*A0XT & + (1.d0-EPS)*(1.d0-EPS-EPS)*A0XX) C---ORDER is important for A0XT! A0YY = A0YY - FYY A0ZZ = A0ZZ - FZZ A0YZ = A0YZ - FYZ A0YT = A0YT - FYT A0ZT = A0ZT - FZT A0XY = A0XY*(1.d0-EPS)**2 A0XZ = A0XZ*(1.d0-EPS)**2 A0XT = (A0XT - 3.d0*EPS*A0XX)*(1.d0-EPS)**2 A0XX = A0XX*(1.d0-EPS)**3 A0TT = A0TT - FTT A0M = A0M - UAB B0M = B0M + UAB EPS = UAB/B0M C combine transferred moments (F's) into box B0 - ORDER important ! B0ZZ = B0ZZ + FZZ B0YY = B0YY + FYY B0YZ = B0YZ + FYZ B0XY = EPS*FXY + (1.d0-EPS)*B0XY & +3.d0*(EPS*B0YT - (1.d0-EPS)*FYT) B0YT = B0YT + FYT B0XZ = EPS*FXZ + (1.d0-EPS)*B0XZ & +3.d0*(EPS*B0ZT - (1.d0-EPS)*FZT) B0ZT = B0ZT + FZT B0XX = FXX*EPS**2 + B0XX*(1.d0-EPS)**2 & +5.d0*( EPS*(1.d0-EPS)*(B0XT - FXT) & +(1.d0-EPS-EPS)*((1.d0-EPS)*FTT-EPS*B0TT)) B0XT = EPS*FXT + (1.d0-EPS)*B0XT & +3.d0*(EPS*B0TT- (1.d0-EPS)*FTT) B0TT = B0TT + FTT C---LIMTR = 3 if (LIMTR.eq.3) then c------------force daughter box B to be monotonic(L=2) first LGC_L2 = B0XX .gt. 0.d0 LGC_LX = B0XT .gt. 0.D0 A3XX = min(B0XX, abs(B0XT)*C033, 0.5d0*B0TT) A3XT = min(B0TT+A3XX, abs(B0XT)) A3XT = merge(A3XT, -A3XT, LGC_LX) A4XT = min(+B0TT, max(-B0TT, B0XT)) A4XX = max(B0XX, -B0TT+abs(A4XT), -abs(A4XT)*C033) B0XT = merge (A3XT, A4XT, LGC_L2) B0XX = merge (A3XX, A4XX, LGC_L2) B0XY = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XY)) B0XZ = min(+CXYZ*B0TT,max(-CXYZ*B0TT,B0XZ)) c------------limit daughter box B by min/max of parent A & B boxes C***NB have had to assume a near-zero level to avoid zero divide c---rescale for X-max ABSCAL = max(0.d0, min(1.d0, & (ABXMX*B0M-B0TT)/max(abs(B0XT)+B0XX, C000))) B0XT = ABSCAL*B0XT B0XX = ABSCAL*B0XX c---rescale for X-min ABSCAL = max(0.d0, min(1.d0, & (B0TT-ABXMN*B0M)/max(abs(B0XT)-B0XX, C000))) B0XT = ABSCAL*B0XT B0XX = ABSCAL*B0XX c---rescale for XY moment ABSCAL = max(0.d0, min(1.d0, & (ABYMX*B0M-B0TT)/max(abs(B0XT)+B0XX+abs(B0XY), C000))) B0XY = ABSCAL*B0XY ABSCAL = max(0.d0, min(1.d0, & (B0TT-ABYMN*B0M)/max(abs(B0XT)-B0XX+abs(B0XY), C000))) B0XY = ABSCAL*B0XY c---rescale for XZ moment ABSCAL = max(0.d0, min(1.d0, & (ABZMX*B0M-B0TT)/max(abs(B0XT)+B0XX+abs(B0XZ), C000))) B0XZ = ABSCAL*B0XZ ABSCAL = max(0.d0, min(1.d0, & (B0TT-ABZMN*B0M)/max(abs(B0XT)-B0XX+abs(B0XZ), C000))) B0XZ = ABSCAL*B0XZ endif C---Now re-mask back onto A1-- AND B1-- arrays A1M = merge (A0M, B0M, LGC_AB) A1TT = merge (A0TT, B0TT, LGC_AB) A1XX = merge (A0XX, B0XX, LGC_AB) A1YT = merge (A0YT, B0YT, LGC_AB) A1YY = merge (A0YY, B0YY, LGC_AB) A1ZT = merge (A0ZT, B0ZT, LGC_AB) A1ZZ = merge (A0ZZ, B0ZZ, LGC_AB) A1YZ = merge (A0YZ, B0YZ, LGC_AB) A1XT = merge (A0XT,-B0XT, LGC_AB) A1XY = merge (A0XY,-B0XY, LGC_AB) A1XZ = merge (A0XZ,-B0XZ, LGC_AB) C--- B1M = merge (B0M, A0M, LGC_AB) B1TT = merge (B0TT, A0TT, LGC_AB) B1XX = merge (B0XX, A0XX, LGC_AB) B1YT = merge (B0YT, A0YT, LGC_AB) B1YY = merge (B0YY, A0YY, LGC_AB) B1ZT = merge (B0ZT, A0ZT, LGC_AB) B1ZZ = merge (B0ZZ, A0ZZ, LGC_AB) B1YZ = merge (B0YZ, A0YZ, LGC_AB) B1XT = merge (B0XT,-A0XT, LGC_AB) B1XY = merge (B0XY,-A0XY, LGC_AB) B1XZ = merge (B0XZ,-A0XZ, LGC_AB) C---re-load to Q--'s: B1M = cshift(B1M, -1) B1TT = cshift(B1TT, -1) B1XT = cshift(B1XT, -1) B1XX = cshift(B1XX, -1) B1YT = cshift(B1YT, -1) B1YY = cshift(B1YY, -1) B1ZT = cshift(B1ZT, -1) B1ZZ = cshift(B1ZZ, -1) B1XY = cshift(B1XY, -1) B1XZ = cshift(B1XZ, -1) B1YZ = cshift(B1YZ, -1) do J=1,NJ C--- B's are 1,3,5,7 now QM(J+J-1) = B1M(J) QTT(J+J-1) = B1TT(J) QXX(J+J-1) = B1XX(J) QYT(J+J-1) = B1YT(J) QYY(J+J-1) = B1YY(J) QZT(J+J-1) = B1ZT(J) QZZ(J+J-1) = B1ZZ(J) QYZ(J+J-1) = B1YZ(J) QXT(J+J-1) = B1XT(J) QXY(J+J-1) = B1XY(J) QXZ(J+J-1) = B1XZ(J) C--- A's are 2,4,6,8,.. QM(J+J) = A1M(J) QTT(J+J) = A1TT(J) QXX(J+J) = A1XX(J) QYT(J+J) = A1YT(J) QYY(J+J) = A1YY(J) QZT(J+J) = A1ZT(J) QZZ(J+J) = A1ZZ(J) QYZ(J+J) = A1YZ(J) QXT(J+J) = A1XT(J) QXY(J+J) = A1XY(J) QXZ(J+J) = A1XZ(J) enddo deallocate (A0M,A0TT,A0XT,A0XX,A0YT,A0YY,A0ZT,A0ZZ,A0XY,A0YZ,A0XZ) deallocate (A1M,A1TT,A1XT,A1XX,A1YT,A1YY,A1ZT,A1ZZ,A1XY,A1YZ,A1XZ) deallocate (B0M,B0TT,B0XT,B0XX,B0YT,B0YY,B0ZT,B0ZZ,B0XY,B0YZ,B0XZ) deallocate (B1M,B1TT,B1XT,B1XX,B1YT,B1YY,B1ZT,B1ZZ,B1XY,B1YZ,B1XZ) deallocate (FTT,FXX,FYY,FZZ,FXT,FYT,FZT,FXY,FXZ,FYZ, EPS) deallocate (UAB, LGC_AB) deallocate (A3XT,A3XX,A4XT,A4XX,LGC_L2,LGC_LX) deallocate (ABSCAL, ABXMX,ABXMN,ABYMX,ABYMN,ABZMX,ABZMN) deallocate (XXAMX,XYAMX,XZAMX,XXAMN,XYAMN,XZAMN) deallocate (XXBMX,XYBMX,XZBMX,XXBMN,XYBMN,XZBMN) return end