source: LMDZ6/trunk/libf/dyn3d_common/ppm3d.f90 @ 5456

Last change on this file since 5456 was 5246, checked in by abarral, 2 months ago

Convert fixed-form to free-form sources .F -> .{f,F}90
(WIP: some .F remain, will be handled in subsequent commits)

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 50.9 KB
Line 
1!
2! $Id: ppm3d.f90 5246 2024-10-21 12:58:45Z aborella $
3!
4
5!From lin@explorer.gsfc.nasa.gov Wed Apr 15 17:44:44 1998
6!Date: Wed, 15 Apr 1998 11:37:03 -0400
7!From: lin@explorer.gsfc.nasa.gov
8!To: Frederic.Hourdin@lmd.jussieu.fr
9!Subject: 3D transport module of the GSFC CTM and GEOS GCM
10
11
12!This code is sent to you by S-J Lin, DAO, NASA-GSFC
13
14!Note: this version is intended for machines like CRAY
15!-90. No multitasking directives implemented.
16
17
18! ********************************************************************
19!
20! TransPort Core for Goddard Chemistry Transport Model (G-CTM), Goddard
21! Earth Observing System General Circulation Model (GEOS-GCM), and Data
22! Assimilation System (GEOS-DAS).
23!
24! ********************************************************************
25!
26! Purpose: given horizontal winds on  a hybrid sigma-p surfaces,
27!      one call to tpcore updates the 3-D mixing ratio
28!      fields one time step (NDT). [vertical mass flux is computed
29!      internally consistent with the discretized hydrostatic mass
30!      continuity equation of the C-Grid GEOS-GCM (for IGD=1)].
31!
32! Schemes: Multi-dimensional Flux Form Semi-Lagrangian (FFSL) scheme based
33!      on the van Leer or PPM.
34!      (see Lin and Rood 1996).
35! Version 4.5
36! Last modified: Dec. 5, 1996
37! Major changes from version 4.0: a more general vertical hybrid sigma-
38! pressure coordinate.
39! Subroutines modified: xtp, ytp, fzppm, qckxyz
40! Subroutines deleted: vanz
41!
42! Author: Shian-Jiann Lin
43! mail address:
44!             Shian-Jiann Lin*
45!             Code 910.3, NASA/GSFC, Greenbelt, MD 20771
46!             Phone: 301-286-9540
47!             E-mail: lin@dao.gsfc.nasa.gov
48!
49! *affiliation:
50!             Joint Center for Earth Systems Technology
51!             The University of Maryland Baltimore County
52!             NASA - Goddard Space Flight Center
53! References:
54!
55! 1. Lin, S.-J., and R. B. Rood, 1996: Multidimensional flux form semi-
56!    Lagrangian transport schemes. Mon. Wea. Rev., 124, 2046-2070.
57!
58! 2. Lin, S.-J., W. C. Chao, Y. C. Sud, and G. K. Walker, 1994: A class of
59!    the van Leer-type transport schemes and its applications to the moist-
60!    ure transport in a General Circulation Model. Mon. Wea. Rev., 122,
61!    1575-1593.
62!
63! ****6***0*********0*********0*********0*********0*********0**********72
64!
65subroutine ppm3d(IGD,Q,PS1,PS2,U,V,W,NDT,IORD,JORD,KORD,NC,IMR, &
66        JNP,j1,NLAY,AP,BP,PT,AE,fill,dum,Umax)
67
68  implicit none
69
70  ! rajout de d�clarations
71  !  integer Jmax,kmax,ndt0,nstep,k,j,i,ic,l,js,jn,imh,iad,jad,krd
72  !  integer iu,iiu,j2,jmr,js0,jt
73  !  real dtdy,dtdy5,rcap,iml,jn0,imjm,pi,dl,dp
74  !  real dt,cr1,maxdt,ztc,d5,sum1,sum2,ru
75  !
76  ! ********************************************************************
77  !
78  ! =============
79  ! INPUT:
80  ! =============
81  !
82  ! Q(IMR,JNP,NLAY,NC): mixing ratios at current time (t)
83  ! NC: total # of constituents
84  ! IMR: first dimension (E-W); # of Grid intervals in E-W is IMR
85  ! JNP: 2nd dimension (N-S); # of Grid intervals in N-S is JNP-1
86  ! NLAY: 3rd dimension (# of layers); vertical index increases from 1 at
87  !   the model top to NLAY near the surface (see fig. below).
88  !   It is assumed that 6 <= NLAY <= JNP (for dynamic memory allocation)
89  !
90  ! PS1(IMR,JNP): surface pressure at current time (t)
91  ! PS2(IMR,JNP): surface pressure at mid-time-level (t+NDT/2)
92  ! PS2 is replaced by the predicted PS (at t+NDT) on output.
93  ! Note: surface pressure can have any unit or can be multiplied by any
94  !   const.
95  !
96  ! The pressure at layer edges are defined as follows:
97  !
98  !    p(i,j,k) = AP(k)*PT  +  BP(k)*PS(i,j)          (1)
99  !
100  ! Where PT is a constant having the same unit as PS.
101  ! AP and BP are unitless constants given at layer edges
102  ! defining the vertical coordinate.
103  ! BP(1) = 0., BP(NLAY+1) = 1.
104  ! The pressure at the model top is PTOP = AP(1)*PT
105  !
106  ! For pure sigma system set AP(k) = 1 for all k, PT = PTOP,
107  ! BP(k) = sige(k) (sigma at edges), PS = Psfc - PTOP.
108  !
109  ! Note: the sigma-P coordinate is a subset of Eq. 1, which in turn
110  ! is a subset of the following even more general sigma-P-thelta coord.
111  ! currently under development.
112  !  p(i,j,k) = (AP(k)*PT + BP(k)*PS(i,j))/(D(k)-C(k)*TE**(-1/kapa))
113  !
114  !              /////////////////////////////////
115  !          / \ ------------- PTOP --------------  AP(1), BP(1)
116  !           |
117  !    delp(1)    |  ........... Q(i,j,1) ............
118  !           |
119  !  W(1)    \ / ---------------------------------  AP(2), BP(2)
120  !
121  !
122  !
123  ! W(k-1)   / \ ---------------------------------  AP(k), BP(k)
124  !           |
125  !    delp(K)    |  ........... Q(i,j,k) ............
126  !           |
127  !  W(k)    \ / ---------------------------------  AP(k+1), BP(k+1)
128  !
129  !
130  !
131  !          / \ ---------------------------------  AP(NLAY), BP(NLAY)
132  !           |
133  !  delp(NLAY)   |  ........... Q(i,j,NLAY) .........
134  !           |
135  !   W(NLAY)=0  \ / ------------- surface ----------- AP(NLAY+1), BP(NLAY+1)
136  !             //////////////////////////////////
137  !
138  ! U(IMR,JNP,NLAY) & V(IMR,JNP,NLAY):winds (m/s) at mid-time-level (t+NDT/2)
139  ! U and V may need to be polar filtered in advance in some cases.
140  !
141  ! IGD:      grid type on which winds are defined.
142  ! IGD = 0:  A-Grid  [all variables defined at the same point from south
143  !               pole (j=1) to north pole (j=JNP) ]
144  !
145  ! IGD = 1  GEOS-GCM C-Grid
146  !                                  [North]
147  !
148  !                                   V(i,j)
149  !                                      |
150  !                                      |
151  !                                      |
152  !                         U(i-1,j)---Q(i,j)---U(i,j) [EAST]
153  !                                      |
154  !                                      |
155  !                                      |
156  !                                   V(i,j-1)
157  !
158  !     U(i,  1) is defined at South Pole.
159  !     V(i,  1) is half grid north of the South Pole.
160  !     V(i,JMR) is half grid south of the North Pole.
161  !
162  !     V must be defined at j=1 and j=JMR if IGD=1
163  !     V at JNP need not be given.
164  !
165  ! NDT: time step in seconds (need not be constant during the course of
166  !  the integration). Suggested value: 30 min. for 4x5, 15 min. for 2x2.5
167  !  (Lat-Lon) resolution. Smaller values are recommanded if the model
168  !  has a well-resolved stratosphere.
169  !
170  ! J1 defines the size of the polar cap:
171  ! South polar cap edge is located at -90 + (j1-1.5)*180/(JNP-1) deg.
172  ! North polar cap edge is located at  90 - (j1-1.5)*180/(JNP-1) deg.
173  ! There are currently only two choices (j1=2 or 3).
174  ! IMR must be an even integer if j1 = 2. Recommended value: J1=3.
175  !
176  ! IORD, JORD, and KORD are integers controlling various options in E-W, N-S,
177  ! and vertical transport, respectively. Recommended values for positive
178  ! definite scalars: IORD=JORD=3, KORD=5. Use KORD=3 for non-
179  ! positive definite scalars or when linear correlation between constituents
180  ! is to be maintained.
181  !
182  !  _ORD=
183  !    1: 1st order upstream scheme (too diffusive, not a useful option; it
184  !       can be used for debugging purposes; this is THE only known "linear"
185  !       monotonic advection scheme.).
186  !    2: 2nd order van Leer (full monotonicity constraint;
187  !       see Lin et al 1994, MWR)
188  !    3: monotonic PPM* (slightly improved PPM of Collela & Woodward 1984)
189  !    4: semi-monotonic PPM (same as 3, but overshoots are allowed)
190  !    5: positive-definite PPM (constraint on the subgrid distribution is
191  !       only strong enough to prevent generation of negative values;
192  !       both overshoots & undershoots are possible).
193  !    6: un-constrained PPM (nearly diffusion free; slightly faster but
194  !       positivity not quaranteed. Use this option only when the fields
195  !       and winds are very smooth).
196  !
197  ! *PPM: Piece-wise Parabolic Method
198  !
199  ! Note that KORD <=2 options are no longer supported. DO not use option 4 or 5.
200  ! for non-positive definite scalars (such as Ertel Potential Vorticity).
201  !
202  ! The implicit numerical diffusion decreases as _ORD increases.
203  ! The last two options (ORDER=5, 6) should only be used when there is
204  ! significant explicit diffusion (such as a turbulence parameterization). You
205  ! might get dispersive results otherwise.
206  ! No filter of any kind is applied to the constituent fields here.
207  !
208  ! AE: Radius of the sphere (meters).
209  ! Recommended value for the planet earth: 6.371E6
210  !
211  ! fill(logical):   flag to do filling for negatives (see note below).
212  !
213  ! Umax: Estimate (upper limit) of the maximum U-wind speed (m/s).
214  ! (220 m/s is a good value for troposphere model; 280 m/s otherwise)
215  !
216  ! =============
217  ! Output
218  ! =============
219  !
220  ! Q: mixing ratios at future time (t+NDT) (original values are over-written)
221  ! W(NLAY): large-scale vertical mass flux as diagnosed from the hydrostatic
222  !      relationship. W will have the same unit as PS1 and PS2 (eg, mb).
223  !      W must be divided by NDT to get the correct mass-flux unit.
224  !      The vertical Courant number C = W/delp_UPWIND, where delp_UPWIND
225  !      is the pressure thickness in the "upwind" direction. For example,
226  !      C(k) = W(k)/delp(k)   if W(k) > 0;
227  !      C(k) = W(k)/delp(k+1) if W(k) < 0.
228  !          ( W > 0 is downward, ie, toward surface)
229  ! PS2: predicted PS at t+NDT (original values are over-written)
230  !
231  ! ********************************************************************
232  ! NOTES:
233  ! This forward-in-time upstream-biased transport scheme reduces to
234  ! the 2nd order center-in-time center-in-space mass continuity eqn.
235  ! if Q = 1 (constant fields will remain constant). This also ensures
236  ! that the computed vertical velocity to be identical to GEOS-1 GCM
237  ! for on-line transport.
238  !
239  ! A larger polar cap is used if j1=3 (recommended for C-Grid winds or when
240  ! winds are noisy near poles).
241  !
242  ! Flux-Form Semi-Lagrangian transport in the East-West direction is used
243  ! when and where Courant # is greater than one.
244  !
245  ! The user needs to change the parameter Jmax or Kmax if the resolution
246  ! is greater than 0.5 deg in N-S or 150 layers in the vertical direction.
247  ! (this TransPort Core is otherwise resolution independent and can be used
248  ! as a library routine).
249  !
250  ! PPM is 4th order accurate when grid spacing is uniform (x & y); 3rd
251  ! order accurate for non-uniform grid (vertical sigma coord.).
252  !
253  ! Time step is limitted only by transport in the meridional direction.
254  ! (the FFSL scheme is not implemented in the meridional direction).
255  !
256  ! Since only 1-D limiters are applied, negative values could
257  ! potentially be generated when large time step is used and when the
258  ! initial fields contain discontinuities.
259  ! This does not necessarily imply the integration is unstable.
260  ! These negatives are typically very small. A filling algorithm is
261  ! activated if the user set "fill" to be true.
262  !
263  ! The van Leer scheme used here is nearly as accurate as the original PPM
264  ! due to the use of a 4th order accurate reference slope. The PPM imple-
265  ! mented here is an improvement over the original and is also based on
266  ! the 4th order reference slope.
267  !
268  ! ****6***0*********0*********0*********0*********0*********0**********72
269  !
270  ! User modifiable parameters
271  !
272  integer,parameter :: Jmax = 361, kmax = 150
273  !
274  ! ****6***0*********0*********0*********0*********0*********0**********72
275  !
276  ! Input-Output arrays
277  !
278
279  real :: Q(IMR,JNP,NLAY,NC),PS1(IMR,JNP),PS2(IMR,JNP), &
280        U(IMR,JNP,NLAY),V(IMR,JNP,NLAY),AP(NLAY+1), &
281        BP(NLAY+1),W(IMR,JNP,NLAY),NDT,val(NLAY),Umax
282  integer :: IGD,IORD,JORD,KORD,NC,IMR,JNP,j1,NLAY,AE
283  integer :: IMRD2
284  real :: PT
285  logical :: cross, fill, dum
286  !
287  ! Local dynamic arrays
288  !
289  real :: CRX(IMR,JNP),CRY(IMR,JNP),xmass(IMR,JNP),ymass(IMR,JNP), &
290        fx1(IMR+1),DPI(IMR,JNP,NLAY),delp1(IMR,JNP,NLAY), &
291        WK1(IMR,JNP,NLAY),PU(IMR,JNP),PV(IMR,JNP),DC2(IMR,JNP), &
292        delp2(IMR,JNP,NLAY),DQ(IMR,JNP,NLAY,NC),VA(IMR,JNP), &
293        UA(IMR,JNP),qtmp(-IMR:2*IMR)
294  !
295  ! Local static  arrays
296  !
297  real :: DTDX(Jmax), DTDX5(Jmax), acosp(Jmax), &
298        cosp(Jmax), cose(Jmax), DAP(kmax),DBK(Kmax)
299  data NDT0, NSTEP /0, 0/
300  data cross /.true./
301  REAL :: DTDY, DTDY5, RCAP
302  INTEGER :: JS0, JN0, IML, JMR, IMJM
303  SAVE DTDY, DTDY5, RCAP, JS0, JN0, IML, &
304        DTDX, DTDX5, ACOSP, COSP, COSE, DAP,DBK
305  !
306  INTEGER :: NDT0, NSTEP, j2, k,j,i,ic,l,JS,JN,IMH
307  INTEGER :: IU,IIU,JT,iad,jad,krd
308  REAL :: r23,r3,PI,DL,DP,DT,CR1,MAXDT,ZTC,D5
309  REAL :: sum1,sum2,ru
310
311  JMR = JNP -1
312  IMJM  = IMR*JNP
313  j2 = JNP - j1 + 1
314  NSTEP = NSTEP + 1
315  !
316  ! *********** Initialization **********************
317  if(NSTEP.eq.1) then
318  !
319  write(6,*) '------------------------------------ '
320  write(6,*) 'NASA/GSFC Transport Core Version 4.5'
321  write(6,*) '------------------------------------ '
322  !
323  WRITE(6,*) 'IMR=',IMR,' JNP=',JNP,' NLAY=',NLAY,' j1=',j1
324  WRITE(6,*) 'NC=',NC,IORD,JORD,KORD,NDT
325  !
326  ! controles sur les parametres
327  if(NLAY.LT.6) then
328    write(6,*) 'NLAY must be >= 6'
329    stop
330  endif
331  if (JNP.LT.NLAY) then
332     write(6,*) 'JNP must be >= NLAY'
333    stop
334  endif
335  IMRD2=mod(IMR,2)
336  if (j1.eq.2.and.IMRD2.NE.0) then
337     write(6,*) 'if j1=2 IMR must be an even integer'
338    stop
339  endif
340
341  !
342  if(Jmax.lt.JNP .or. Kmax.lt.NLAY) then
343    write(6,*) 'Jmax or Kmax is too small'
344    stop
345  endif
346  !
347  DO k=1,NLAY
348  DAP(k) = (AP(k+1) - AP(k))*PT
349  DBK(k) =  BP(k+1) - BP(k)
350  ENDDO
351  !
352  PI = 4. * ATAN(1.)
353  DL = 2.*PI / REAL(IMR)
354  DP =    PI / REAL(JMR)
355  !
356  if(IGD.eq.0) then
357  ! Compute analytic cosine at cell edges
358        call cosa(cosp,cose,JNP,PI,DP)
359  else
360  ! Define cosine consistent with GEOS-GCM (using dycore2.0 or later)
361        call cosc(cosp,cose,JNP,PI,DP)
362  endif
363  !
364  do J=2,JMR
365    acosp(j) = 1. / cosp(j)
366  end do
367  !
368  ! Inverse of the Scaled polar cap area.
369  !
370  RCAP  = DP / (IMR*(1.-COS((j1-1.5)*DP)))
371  acosp(1)   = RCAP
372  acosp(JNP) = RCAP
373  endif
374  !
375  if(NDT0 .ne. NDT) then
376  DT   = NDT
377  NDT0 = NDT
378
379    if(Umax .lt. 180.) then
380     write(6,*) 'Umax may be too small!'
381    endif
382  CR1  = abs(Umax*DT)/(DL*AE)
383  MaxDT = DP*AE / abs(Umax) + 0.5
384  write(6,*)'Largest time step for max(V)=',Umax,' is ',MaxDT
385  if(MaxDT .lt. abs(NDT)) then
386        write(6,*) 'Warning!!! NDT maybe too large!'
387  endif
388  !
389  if(CR1.ge.0.95) then
390  JS0 = 0
391  JN0 = 0
392  IML = IMR-2
393  ZTC = 0.
394  else
395  ZTC  = acos(CR1) * (180./PI)
396  !
397  JS0 = REAL(JMR)*(90.-ZTC)/180. + 2
398  JS0 = max(JS0, J1+1)
399  IML = min(6*JS0/(J1-1)+2, 4*IMR/5)
400  JN0 = JNP-JS0+1
401  endif
402  !
403  !
404  do J=2,JMR
405  DTDX(j)  = DT / ( DL*AE*COSP(J) )
406
407  ! print*,'dtdx=',dtdx(j)
408  DTDX5(j) = 0.5*DTDX(j)
409  enddo
410  !
411
412  DTDY  = DT /(AE*DP)
413   ! print*,'dtdy=',dtdy
414  DTDY5 = 0.5*DTDY
415  !
416  !  write(6,*) 'J1=',J1,' J2=', J2
417  endif
418  !
419  ! *********** End Initialization **********************
420  !
421  ! delp = pressure thickness: the psudo-density in a hydrostatic system.
422  do  k=1,NLAY
423     do  j=1,JNP
424        do  i=1,IMR
425           delp1(i,j,k)=DAP(k)+DBK(k)*PS1(i,j)
426           delp2(i,j,k)=DAP(k)+DBK(k)*PS2(i,j)
427        enddo
428     enddo
429  enddo
430
431  !
432  if(j1.ne.2) then
433  DO IC=1,NC
434  DO L=1,NLAY
435  DO I=1,IMR
436  Q(I,  2,L,IC) = Q(I,  1,L,IC)
437    Q(I,JMR,L,IC) = Q(I,JNP,L,IC)
438  END DO
439  END DO
440  END DO
441  endif
442  !
443  ! Compute "tracer density"
444  DO IC=1,NC
445  DO k=1,NLAY
446  DO j=1,JNP
447  DO i=1,IMR
448    DQ(i,j,k,IC) = Q(i,j,k,IC)*delp1(i,j,k)
449  END DO
450  END DO
451  END DO
452  END DO
453  !
454  do k=1,NLAY
455  !
456  if(IGD.eq.0) then
457  ! Convert winds on A-Grid to Courant # on C-Grid.
458  call A2C(U(1,1,k),V(1,1,k),IMR,JMR,j1,j2,CRX,CRY,dtdx5,DTDY5)
459  else
460  ! Convert winds on C-grid to Courant #
461  do j=j1,j2
462  do i=2,IMR
463    CRX(i,J) = dtdx(j)*U(i-1,j,k)
464  end do
465  end do
466
467  !
468  do j=j1,j2
469    CRX(1,J) = dtdx(j)*U(IMR,j,k)
47050  continue
471  end do
472  !
473  do i=1,IMR*JMR
474    CRY(i,2) = DTDY*V(i,1,k)
475  end do
476  endif
477  !
478  ! Determine JS and JN
479  JS = j1
480  JN = j2
481  !
482  do j=JS0,j1+1,-1
483  do i=1,IMR
484  if(abs(CRX(i,j)).GT.1.) then
485        JS = j
486        go to 2222
487  endif
488  enddo
489  enddo
490  !
4912222   continue
492  do j=JN0,j2-1
493  do i=1,IMR
494  if(abs(CRX(i,j)).GT.1.) then
495        JN = j
496        go to 2233
497  endif
498  enddo
499  enddo
5002233   continue
501  !
502  if(j1.ne.2) then           ! Enlarged polar cap.
503  do i=1,IMR
504  DPI(i,  2,k) = 0.
505  DPI(i,JMR,k) = 0.
506  enddo
507  endif
508  !
509  ! ******* Compute horizontal mass fluxes ************
510  !
511  ! N-S component
512  do j=j1,j2+1
513  D5 = 0.5 * COSE(j)
514  do i=1,IMR
515  ymass(i,j) = CRY(i,j)*D5*(delp2(i,j,k) + delp2(i,j-1,k))
516  enddo
517  enddo
518  !
519  do j=j1,j2
520  DO i=1,IMR
521    DPI(i,j,k) = (ymass(i,j) - ymass(i,j+1)) * acosp(j)
522  END DO
523  end do
524  !
525  ! Poles
526  sum1 = ymass(IMR,j1  )
527  sum2 = ymass(IMR,J2+1)
528  do i=1,IMR-1
529  sum1 = sum1 + ymass(i,j1  )
530  sum2 = sum2 + ymass(i,J2+1)
531  enddo
532  !
533  sum1 = - sum1 * RCAP
534  sum2 =   sum2 * RCAP
535  do i=1,IMR
536  DPI(i,  1,k) = sum1
537  DPI(i,JNP,k) = sum2
538  enddo
539  !
540  ! E-W component
541  !
542  do j=j1,j2
543  do i=2,IMR
544  PU(i,j) = 0.5 * (delp2(i,j,k) + delp2(i-1,j,k))
545  enddo
546  enddo
547  !
548  do j=j1,j2
549  PU(1,j) = 0.5 * (delp2(1,j,k) + delp2(IMR,j,k))
550  enddo
551  !
552  do j=j1,j2
553  DO i=1,IMR
554    xmass(i,j) = PU(i,j)*CRX(i,j)
555  END DO
556  end do
557  !
558  DO j=j1,j2
559  DO i=1,IMR-1
560    DPI(i,j,k) = DPI(i,j,k) + xmass(i,j) - xmass(i+1,j)
561  END DO
562  END DO
563  !
564  DO j=j1,j2
565    DPI(IMR,j,k) = DPI(IMR,j,k) + xmass(IMR,j) - xmass(1,j)
566  END DO
567  !
568  DO j=j1,j2
569  do i=1,IMR-1
570  UA(i,j) = 0.5 * (CRX(i,j)+CRX(i+1,j))
571  enddo
572  enddo
573  !
574  DO j=j1,j2
575  UA(imr,j) = 0.5 * (CRX(imr,j)+CRX(1,j))
576  enddo
577  !cccccccccccccccccccccccccccccccccccccccccccccccccccccc
578  ! Rajouts pour LMDZ.3.3
579  !cccccccccccccccccccccccccccccccccccccccccccccccccccccc
580  do i=1,IMR
581     do j=1,JNP
582         VA(i,j)=0.
583     enddo
584  enddo
585
586  do i=1,imr*(JMR-1)
587  VA(i,2) = 0.5*(CRY(i,2)+CRY(i,3))
588  enddo
589  !
590  if(j1.eq.2) then
591    IMH = IMR/2
592  do i=1,IMH
593  VA(i,      1) = 0.5*(CRY(i,2)-CRY(i+IMH,2))
594  VA(i+IMH,  1) = -VA(i,1)
595  VA(i,    JNP) = 0.5*(CRY(i,JNP)-CRY(i+IMH,JMR))
596  VA(i+IMH,JNP) = -VA(i,JNP)
597  enddo
598  VA(IMR,1)=VA(1,1)
599  VA(IMR,JNP)=VA(1,JNP)
600  endif
601  !
602  ! ****6***0*********0*********0*********0*********0*********0**********72
603  do IC=1,NC
604  !
605  do i=1,IMJM
606  wk1(i,1,1) = 0.
607  wk1(i,1,2) = 0.
608  enddo
609  !
610  ! E-W advective cross term
611  do j=J1,J2
612  if(J.GT.JS  .and. J.LT.JN) GO TO 250
613  !
614  do i=1,IMR
615  qtmp(i) = q(i,j,k,IC)
616  enddo
617  !
618  do i=-IML,0
619  qtmp(i)       = q(IMR+i,j,k,IC)
620  qtmp(IMR+1-i) = q(1-i,j,k,IC)
621  enddo
622  !
623  DO i=1,IMR
624  iu = UA(i,j)
625  ru = UA(i,j) - iu
626  iiu = i-iu
627  if(UA(i,j).GE.0.) then
628  wk1(i,j,1) = qtmp(iiu)+ru*(qtmp(iiu-1)-qtmp(iiu))
629  else
630  wk1(i,j,1) = qtmp(iiu)+ru*(qtmp(iiu)-qtmp(iiu+1))
631  endif
632  wk1(i,j,1) = wk1(i,j,1) - qtmp(i)
633  END DO
634250 continue
635  end do
636  !
637  if(JN.ne.0) then
638  do j=JS+1,JN-1
639  !
640  do i=1,IMR
641  qtmp(i) = q(i,j,k,IC)
642  enddo
643  !
644  qtmp(0)     = q(IMR,J,k,IC)
645  qtmp(IMR+1) = q(  1,J,k,IC)
646  !
647  do i=1,imr
648  iu = i - UA(i,j)
649  wk1(i,j,1) = UA(i,j)*(qtmp(iu) - qtmp(iu+1))
650  enddo
651  enddo
652  endif
653  ! ****6***0*********0*********0*********0*********0*********0**********72
654  ! Contribution from the N-S advection
655  do i=1,imr*(j2-j1+1)
656  JT = REAL(J1) - VA(i,j1)
657  wk1(i,j1,2) = VA(i,j1) * (q(i,jt,k,IC) - q(i,jt+1,k,IC))
658  enddo
659  !
660  do i=1,IMJM
661  wk1(i,1,1) = q(i,1,k,IC) + 0.5*wk1(i,1,1)
662  wk1(i,1,2) = q(i,1,k,IC) + 0.5*wk1(i,1,2)
663  enddo
664  !
665    if(cross) then
666  ! Add cross terms in the vertical direction.
667    if(IORD .GE. 2) then
668            iad = 2
669    else
670            iad = 1
671    endif
672  !
673    if(JORD .GE. 2) then
674            jad = 2
675    else
676            jad = 1
677    endif
678  call xadv(IMR,JNP,j1,j2,wk1(1,1,2),UA,JS,JN,IML,DC2,iad)
679  call yadv(IMR,JNP,j1,j2,wk1(1,1,1),VA,PV,W,jad)
680  do j=1,JNP
681  do i=1,IMR
682  q(i,j,k,IC) = q(i,j,k,IC) + DC2(i,j) + PV(i,j)
683  enddo
684  enddo
685  endif
686  !
687  call xtp(IMR,JNP,IML,j1,j2,JN,JS,PU,DQ(1,1,k,IC),wk1(1,1,2) &
688        ,CRX,fx1,xmass,IORD)
689
690  call ytp(IMR,JNP,j1,j2,acosp,RCAP,DQ(1,1,k,IC),wk1(1,1,1),CRY, &
691        DC2,ymass,WK1(1,1,3),wk1(1,1,4),WK1(1,1,5),WK1(1,1,6),JORD)
692  !
693  end do
694  end do
695  !
696  ! ******* Compute vertical mass flux (same unit as PS) ***********
697  !
698  ! 1st step: compute total column mass CONVERGENCE.
699  !
700  do j=1,JNP
701  do i=1,IMR
702    CRY(i,j) = DPI(i,j,1)
703  end do
704  end do
705  !
706  do k=2,NLAY
707  do j=1,JNP
708  do i=1,IMR
709  CRY(i,j)  = CRY(i,j) + DPI(i,j,k)
710  end do
711  end do
712  end do
713  !
714  do j=1,JNP
715  do i=1,IMR
716  !
717  ! 2nd step: compute PS2 (PS at n+1) using the hydrostatic assumption.
718  ! Changes (increases) to surface pressure = total column mass convergence
719  !
720  PS2(i,j)  = PS1(i,j) + CRY(i,j)
721  !
722  ! 3rd step: compute vertical mass flux from mass conservation principle.
723  !
724  W(i,j,1) = DPI(i,j,1) - DBK(1)*CRY(i,j)
725  W(i,j,NLAY) = 0.
726  end do
727  end do
728  !
729  do k=2,NLAY-1
730  do j=1,JNP
731  do i=1,IMR
732  W(i,j,k) = W(i,j,k-1) + DPI(i,j,k) - DBK(k)*CRY(i,j)
733  end do
734  end do
735  end do
736  !
737  DO k=1,NLAY
738  DO j=1,JNP
739  DO i=1,IMR
740  delp2(i,j,k) = DAP(k) + DBK(k)*PS2(i,j)
741  END DO
742  END DO
743  END DO
744  !
745    KRD = max(3, KORD)
746  do IC=1,NC
747  !
748  !****6***0*********0*********0*********0*********0*********0**********72
749
750  call FZPPM(IMR,JNP,NLAY,j1,DQ(1,1,1,IC),W,Q(1,1,1,IC),WK1,DPI, &
751        DC2,CRX,CRY,PU,PV,xmass,ymass,delp1,KRD)
752  !
753
754  if(fill) call qckxyz(DQ(1,1,1,IC),DC2,IMR,JNP,NLAY,j1,j2, &
755        cosp,acosp,.false.,IC,NSTEP)
756  !
757  ! Recover tracer mixing ratio from "density" using predicted
758  ! "air density" (pressure thickness) at time-level n+1
759  !
760  DO k=1,NLAY
761  DO j=1,JNP
762  DO i=1,IMR
763        Q(i,j,k,IC) = DQ(i,j,k,IC) / delp2(i,j,k)
764         ! print*,'i=',i,'j=',j,'k=',k,'Q(i,j,k,IC)=',Q(i,j,k,IC)
765  enddo
766  enddo
767  enddo
768  !
769  if(j1.ne.2) then
770  DO k=1,NLAY
771  DO I=1,IMR
772  ! j=1 c'est le p�le Sud, j=JNP c'est le p�le Nord
773  Q(I,  2,k,IC) = Q(I,  1,k,IC)
774  Q(I,JMR,k,IC) = Q(I,JNP,k,IC)
775  END DO
776  END DO
777  endif
778  end do
779  !
780  if(j1.ne.2) then
781  DO k=1,NLAY
782  DO i=1,IMR
783  W(i,  2,k) = W(i,  1,k)
784  W(i,JMR,k) = W(i,JNP,k)
785  END DO
786  END DO
787  endif
788  !
789  RETURN
790END SUBROUTINE ppm3d
791!
792!****6***0*********0*********0*********0*********0*********0**********72
793subroutine FZPPM(IMR,JNP,NLAY,j1,DQ,WZ,P,DC,DQDT,AR,AL,A6, &
794        flux,wk1,wk2,wz2,delp,KORD)
795  implicit none
796  integer,parameter :: kmax = 150
797  real,parameter :: R23 = 2./3., R3 = 1./3.
798  integer :: IMR,JNP,NLAY,J1,KORD
799  real :: WZ(IMR,JNP,NLAY),P(IMR,JNP,NLAY),DC(IMR,JNP,NLAY), &
800        wk1(IMR,*),delp(IMR,JNP,NLAY),DQ(IMR,JNP,NLAY), &
801        DQDT(IMR,JNP,NLAY)
802  ! Assuming JNP >= NLAY
803  real :: AR(IMR,*),AL(IMR,*),A6(IMR,*),flux(IMR,*),wk2(IMR,*), &
804        wz2(IMR,*)
805  integer :: JMR,IMJM,NLAYM1,LMT,K,I,J
806  real :: c0,c1,c2,tmp,qmax,qmin,a,b,fct,a1,a2,cm,cp
807  !
808  JMR = JNP - 1
809  IMJM = IMR*JNP
810  NLAYM1 = NLAY - 1
811  !
812  LMT = KORD - 3
813  !
814  ! ****6***0*********0*********0*********0*********0*********0**********72
815  ! Compute DC for PPM
816  ! ****6***0*********0*********0*********0*********0*********0**********72
817  !
818  do k=1,NLAYM1
819  do i=1,IMJM
820  DQDT(i,1,k) = P(i,1,k+1) - P(i,1,k)
821  end do
822  end do
823  !
824  DO k=2,NLAYM1
825  DO I=1,IMJM
826   c0 =  delp(i,1,k) / (delp(i,1,k-1)+delp(i,1,k)+delp(i,1,k+1))
827   c1 = (delp(i,1,k-1)+0.5*delp(i,1,k))/(delp(i,1,k+1)+delp(i,1,k))
828   c2 = (delp(i,1,k+1)+0.5*delp(i,1,k))/(delp(i,1,k-1)+delp(i,1,k))
829  tmp = c0*(c1*DQDT(i,1,k) + c2*DQDT(i,1,k-1))
830  Qmax = max(P(i,1,k-1),P(i,1,k),P(i,1,k+1)) - P(i,1,k)
831  Qmin = P(i,1,k) - min(P(i,1,k-1),P(i,1,k),P(i,1,k+1))
832  DC(i,1,k) = sign(min(abs(tmp),Qmax,Qmin), tmp)
833  END DO
834  END DO
835
836  !
837  ! ****6***0*********0*********0*********0*********0*********0**********72
838  ! Loop over latitudes  (to save memory)
839  ! ****6***0*********0*********0*********0*********0*********0**********72
840  !
841  DO j=1,JNP
842  if((j.eq.2 .or. j.eq.JMR) .and. j1.ne.2) goto 2000
843  !
844  DO k=1,NLAY
845  DO i=1,IMR
846  wz2(i,k) =   WZ(i,j,k)
847  wk1(i,k) =    P(i,j,k)
848  wk2(i,k) = delp(i,j,k)
849  flux(i,k) = DC(i,j,k)  !this flux is actually the monotone slope
850  enddo
851  enddo
852  !
853  !****6***0*********0*********0*********0*********0*********0**********72
854  ! Compute first guesses at cell interfaces
855  ! First guesses are required to be continuous.
856  ! ****6***0*********0*********0*********0*********0*********0**********72
857  !
858  ! three-cell parabolic subgrid distribution at model top
859  ! two-cell parabolic with zero gradient subgrid distribution
860  ! at the surface.
861  !
862  ! First guess top edge value
863  DO i=1,IMR
864  ! three-cell PPM
865  ! Compute a,b, and c of q = aP**2 + bP + c using cell averages and delp
866  a = 3.*( DQDT(i,j,2) - DQDT(i,j,1)*(wk2(i,2)+wk2(i,3))/ &
867        (wk2(i,1)+wk2(i,2)) ) / &
868        ( (wk2(i,2)+wk2(i,3))*(wk2(i,1)+wk2(i,2)+wk2(i,3)) )
869  b = 2.*DQDT(i,j,1)/(wk2(i,1)+wk2(i,2)) - &
870        R23*a*(2.*wk2(i,1)+wk2(i,2))
871  AL(i,1) =  wk1(i,1) - wk2(i,1)*(R3*a*wk2(i,1) + 0.5*b)
872  AL(i,2) =  wk2(i,1)*(a*wk2(i,1) + b) + AL(i,1)
873  !
874  ! Check if change sign
875  if(wk1(i,1)*AL(i,1).le.0.) then
876             AL(i,1) = 0.
877         flux(i,1) = 0.
878    else
879         flux(i,1) =  wk1(i,1) - AL(i,1)
880    endif
881  END DO
882  !
883  ! Bottom
884  DO i=1,IMR
885  ! 2-cell PPM with zero gradient right at the surface
886  !
887  fct = DQDT(i,j,NLAYM1)*wk2(i,NLAY)**2 / &
888        ( (wk2(i,NLAY)+wk2(i,NLAYM1))*(2.*wk2(i,NLAY)+wk2(i,NLAYM1)))
889  AR(i,NLAY) = wk1(i,NLAY) + fct
890  AL(i,NLAY) = wk1(i,NLAY) - (fct+fct)
891  if(wk1(i,NLAY)*AR(i,NLAY).le.0.) AR(i,NLAY) = 0.
892  flux(i,NLAY) = AR(i,NLAY) -  wk1(i,NLAY)
893  END DO
894
895  !
896  !****6***0*********0*********0*********0*********0*********0**********72
897  ! 4th order interpolation in the interior.
898  !****6***0*********0*********0*********0*********0*********0**********72
899  !
900  DO k=3,NLAYM1
901  DO i=1,IMR
902  c1 =  DQDT(i,j,k-1)*wk2(i,k-1) / (wk2(i,k-1)+wk2(i,k))
903  c2 =  2. / (wk2(i,k-2)+wk2(i,k-1)+wk2(i,k)+wk2(i,k+1))
904  A1   =  (wk2(i,k-2)+wk2(i,k-1)) / (2.*wk2(i,k-1)+wk2(i,k))
905  A2   =  (wk2(i,k  )+wk2(i,k+1)) / (2.*wk2(i,k)+wk2(i,k-1))
906  AL(i,k) = wk1(i,k-1) + c1 + c2 * &
907        ( wk2(i,k  )*(c1*(A1 - A2)+A2*flux(i,k-1)) - &
908        wk2(i,k-1)*A1*flux(i,k)  )
909   ! print *,'AL1',i,k, AL(i,k)
910  END DO
911  END DO
912  !
913  do i=1,IMR*NLAYM1
914  AR(i,1) = AL(i,2)
915   ! print *,'AR1',i,AR(i,1)
916  end do
917  !
918  do i=1,IMR*NLAY
919  A6(i,1) = 3.*(wk1(i,1)+wk1(i,1) - (AL(i,1)+AR(i,1)))
920   ! print *,'A61',i,A6(i,1)
921  end do
922  !
923  !****6***0*********0*********0*********0*********0*********0**********72
924  ! Top & Bot always monotonic
925  call lmtppm(flux(1,1),A6(1,1),AR(1,1),AL(1,1),wk1(1,1),IMR,0)
926  call lmtppm(flux(1,NLAY),A6(1,NLAY),AR(1,NLAY),AL(1,NLAY), &
927        wk1(1,NLAY),IMR,0)
928  !
929  ! Interior depending on KORD
930  if(LMT.LE.2) &
931        call lmtppm(flux(1,2),A6(1,2),AR(1,2),AL(1,2),wk1(1,2), &
932        IMR*(NLAY-2),LMT)
933  !
934  !****6***0*********0*********0*********0*********0*********0**********72
935  !
936  DO i=1,IMR*NLAYM1
937  IF(wz2(i,1).GT.0.) then
938    CM = wz2(i,1) / wk2(i,1)
939    flux(i,2) = AR(i,1)+0.5*CM*(AL(i,1)-AR(i,1)+A6(i,1)*(1.-R23*CM))
940  else
941     ! print *,'test2-0',i,j,wz2(i,1),wk2(i,2)
942    CP= wz2(i,1) / wk2(i,2)
943     ! print *,'testCP',CP
944    flux(i,2) = AL(i,2)+0.5*CP*(AL(i,2)-AR(i,2)-A6(i,2)*(1.+R23*CP))
945     ! print *,'test2',i, AL(i,2),AR(i,2),A6(i,2),R23
946  endif
947  END DO
948  !
949  DO i=1,IMR*NLAYM1
950  flux(i,2) = wz2(i,1) * flux(i,2)
951250 CONTINUE
952  END DO
953  !
954  do i=1,IMR
955  DQ(i,j,   1) = DQ(i,j,   1) - flux(i,   2)
956  DQ(i,j,NLAY) = DQ(i,j,NLAY) + flux(i,NLAY)
957  end do
958  !
959  do k=2,NLAYM1
960  do i=1,IMR
961    DQ(i,j,k) = DQ(i,j,k) + flux(i,k) - flux(i,k+1)
962  end do
963  end do
9642000 CONTINUE
965  END DO
966  return
967end subroutine fzppm
968!
969subroutine xtp(IMR,JNP,IML,j1,j2,JN,JS,PU,DQ,Q,UC, &
970        fx1,xmass,IORD)
971  implicit none
972  integer :: IMR,JNP,IML,j1,j2,JN,JS,IORD
973  real :: PU,DQ,Q,UC,fx1,xmass
974  real :: dc,qtmp
975  integer :: ISAVE(IMR)
976  dimension UC(IMR,*),DC(-IML:IMR+IML+1),xmass(IMR,JNP) &
977        ,fx1(IMR+1),DQ(IMR,JNP),qtmp(-IML:IMR+1+IML)
978  dimension PU(IMR,JNP),Q(IMR,JNP)
979  integer :: jvan,j1vl,j2vl,j,i,iu,itmp,ist,imp
980  real :: rut
981  !
982  IMP = IMR + 1
983  !
984  ! van Leer at high latitudes
985  jvan = max(1,JNP/18)
986  j1vl = j1+jvan
987  j2vl = j2-jvan
988  !
989  do j=j1,j2
990  !
991  do i=1,IMR
992  qtmp(i) = q(i,j)
993  enddo
994  !
995  if(j.ge.JN .or. j.le.JS) goto 2222
996  ! ************* Eulerian **********
997  !
998  qtmp(0)     = q(IMR,J)
999  qtmp(-1)    = q(IMR-1,J)
1000  qtmp(IMP)   = q(1,J)
1001  qtmp(IMP+1) = q(2,J)
1002  !
1003  IF(IORD.eq.1 .or. j.eq.j1.OR. j.eq.j2) THEN
1004  DO i=1,IMR
1005  iu = REAL(i) - uc(i,j)
1006    fx1(i) = qtmp(iu)
1007  END DO
1008  ELSE
1009  call xmist(IMR,IML,Qtmp,DC)
1010  DC(0) = DC(IMR)
1011  !
1012  if(IORD.eq.2 .or. j.le.j1vl .or. j.ge.j2vl) then
1013  DO i=1,IMR
1014  iu = REAL(i) - uc(i,j)
1015    fx1(i) = qtmp(iu) + DC(iu)*(sign(1.,uc(i,j))-uc(i,j))
1016  END DO
1017  else
1018  call fxppm(IMR,IML,UC(1,j),Qtmp,DC,fx1,IORD)
1019  endif
1020  !
1021  ENDIF
1022  !
1023  DO i=1,IMR
1024    fx1(i) = fx1(i)*xmass(i,j)
1025  END DO
1026  !
1027  goto 1309
1028  !
1029  ! ***** Conservative (flux-form) Semi-Lagrangian transport *****
1030  !
10312222   continue
1032  !
1033  do i=-IML,0
1034  qtmp(i)     = q(IMR+i,j)
1035  qtmp(IMP-i) = q(1-i,j)
1036  enddo
1037  !
1038  IF(IORD.eq.1 .or. j.eq.j1.OR. j.eq.j2) THEN
1039  DO i=1,IMR
1040  itmp = INT(uc(i,j))
1041  ISAVE(i) = i - itmp
1042  iu = i - uc(i,j)
1043    fx1(i) = (uc(i,j) - itmp)*qtmp(iu)
1044  END DO
1045  ELSE
1046  call xmist(IMR,IML,Qtmp,DC)
1047  !
1048  do i=-IML,0
1049  DC(i)     = DC(IMR+i)
1050  DC(IMP-i) = DC(1-i)
1051  enddo
1052  !
1053  DO i=1,IMR
1054  itmp = INT(uc(i,j))
1055  rut  = uc(i,j) - itmp
1056  ISAVE(i) = i - itmp
1057  iu = i - uc(i,j)
1058    fx1(i) = rut*(qtmp(iu) + DC(iu)*(sign(1.,rut) - rut))
1059  END DO
1060  ENDIF
1061  !
1062  do i=1,IMR
1063  IF(uc(i,j).GT.1.) then
1064  !DIR$ NOVECTOR
1065    do ist = ISAVE(i),i-1
1066    fx1(i) = fx1(i) + qtmp(ist)
1067    enddo
1068  elseIF(uc(i,j).LT.-1.) then
1069    do ist = i,ISAVE(i)-1
1070    fx1(i) = fx1(i) - qtmp(ist)
1071    enddo
1072  !DIR$ VECTOR
1073  endif
1074  end do
1075  do i=1,IMR
1076  fx1(i) = PU(i,j)*fx1(i)
1077  enddo
1078  !
1079  ! ***************************************
1080  !
10811309   fx1(IMP) = fx1(1)
1082  DO i=1,IMR
1083    DQ(i,j) =  DQ(i,j) + fx1(i)-fx1(i+1)
1084  END DO
1085  !
1086  ! ***************************************
1087  !
1088  end do
1089  return
1090end subroutine xtp
1091!
1092subroutine fxppm(IMR,IML,UT,P,DC,flux,IORD)
1093  implicit none
1094  integer :: IMR,IML,IORD
1095  real :: UT,P,DC,flux
1096  real,parameter ::  R3 = 1./3., R23 = 2./3.
1097  DIMENSION UT(*),flux(*),P(-IML:IMR+IML+1),DC(-IML:IMR+IML+1)
1098  REAL :: AR(0:IMR),AL(0:IMR),A6(0:IMR)
1099  integer :: LMT,IMP,JLVL,i
1100   ! logical first
1101   ! data first /.true./
1102   ! SAVE LMT
1103   ! if(first) then
1104  !
1105  ! correction calcul de LMT a chaque passage pour pouvoir choisir
1106  ! plusieurs schemas PPM pour differents traceurs
1107  !  IF (IORD.LE.0) then
1108  !        if(IMR.GE.144) then
1109  !              LMT = 0
1110  !        elseif(IMR.GE.72) then
1111  !              LMT = 1
1112  !        else
1113  !              LMT = 2
1114  !        endif
1115  !  else
1116  !        LMT = IORD - 3
1117  !  endif
1118  !
1119  LMT = IORD - 3
1120   ! write(6,*) 'PPM option in E-W direction = ', LMT
1121   ! first = .false.
1122   ! endif
1123  !
1124  DO i=1,IMR
1125    AL(i) = 0.5*(p(i-1)+p(i)) + (DC(i-1) - DC(i))*R3
1126  END DO
1127  !
1128  do i=1,IMR-1
1129    AR(i) = AL(i+1)
1130  end do
1131  AR(IMR) = AL(1)
1132  !
1133  do i=1,IMR
1134    A6(i) = 3.*(p(i)+p(i)  - (AL(i)+AR(i)))
1135  end do
1136  !
1137  if(LMT.LE.2) call lmtppm(DC(1),A6(1),AR(1),AL(1),P(1),IMR,LMT)
1138  !
1139  AL(0) = AL(IMR)
1140  AR(0) = AR(IMR)
1141  A6(0) = A6(IMR)
1142  !
1143  DO i=1,IMR
1144  IF(UT(i).GT.0.) then
1145  flux(i) = AR(i-1) + 0.5*UT(i)*(AL(i-1) - AR(i-1) + &
1146        A6(i-1)*(1.-R23*UT(i)) )
1147  else
1148  flux(i) = AL(i) - 0.5*UT(i)*(AR(i) - AL(i) + &
1149        A6(i)*(1.+R23*UT(i)))
1150  endif
1151  enddo
1152  return
1153end subroutine fxppm
1154!
1155subroutine xmist(IMR,IML,P,DC)
1156  implicit none
1157  integer :: IMR,IML
1158  real,parameter :: R24 = 1./24.
1159  real :: P(-IML:IMR+1+IML),DC(-IML:IMR+1+IML)
1160  integer :: i
1161  real :: tmp,pmax,pmin
1162  !
1163  do  i=1,IMR
1164  tmp = R24*(8.*(p(i+1) - p(i-1)) + p(i-2) - p(i+2))
1165  Pmax = max(P(i-1), p(i), p(i+1)) - p(i)
1166  Pmin = p(i) - min(P(i-1), p(i), p(i+1))
1167    DC(i) = sign(min(abs(tmp),Pmax,Pmin), tmp)
1168  end do
1169  return
1170end subroutine xmist
1171!
1172subroutine ytp(IMR,JNP,j1,j2,acosp,RCAP,DQ,P,VC,DC2 &
1173        ,ymass,fx,A6,AR,AL,JORD)
1174  implicit none
1175  integer :: IMR,JNP,j1,j2,JORD
1176  real :: acosp,RCAP,DQ,P,VC,DC2,ymass,fx,A6,AR,AL
1177  dimension P(IMR,JNP),VC(IMR,JNP),ymass(IMR,JNP) &
1178        ,DC2(IMR,JNP),DQ(IMR,JNP),acosp(JNP)
1179  ! Work array
1180  DIMENSION fx(IMR,JNP),AR(IMR,JNP),AL(IMR,JNP),A6(IMR,JNP)
1181  integer :: JMR,len,i,jt,j
1182  real :: sum1,sum2
1183  !
1184  JMR = JNP - 1
1185  len = IMR*(J2-J1+2)
1186  !
1187  if(JORD.eq.1) then
1188  DO i=1,len
1189  JT = REAL(J1) - VC(i,J1)
1190    fx(i,j1) = p(i,JT)
1191  END DO
1192  else
1193
1194  call ymist(IMR,JNP,j1,P,DC2,4)
1195  !
1196  if(JORD.LE.0 .or. JORD.GE.3) then
1197
1198  call fyppm(VC,P,DC2,fx,IMR,JNP,j1,j2,A6,AR,AL,JORD)
1199
1200  else
1201  DO i=1,len
1202  JT = REAL(J1) - VC(i,J1)
1203    fx(i,j1) = p(i,JT) + (sign(1.,VC(i,j1))-VC(i,j1))*DC2(i,JT)
1204  END DO
1205  endif
1206  endif
1207  !
1208  DO i=1,len
1209    fx(i,j1) = fx(i,j1)*ymass(i,j1)
1210  END DO
1211  !
1212  DO j=j1,j2
1213  DO i=1,IMR
1214    DQ(i,j) = DQ(i,j) + (fx(i,j) - fx(i,j+1)) * acosp(j)
1215  END DO
1216  END DO
1217  !
1218  ! Poles
1219  sum1 = fx(IMR,j1  )
1220  sum2 = fx(IMR,J2+1)
1221  do i=1,IMR-1
1222  sum1 = sum1 + fx(i,j1  )
1223  sum2 = sum2 + fx(i,J2+1)
1224  enddo
1225  !
1226  sum1 = DQ(1,  1) - sum1 * RCAP
1227  sum2 = DQ(1,JNP) + sum2 * RCAP
1228  do i=1,IMR
1229  DQ(i,  1) = sum1
1230  DQ(i,JNP) = sum2
1231  enddo
1232  !
1233  if(j1.ne.2) then
1234  do i=1,IMR
1235  DQ(i,  2) = sum1
1236  DQ(i,JMR) = sum2
1237  enddo
1238  endif
1239  !
1240  return
1241end subroutine ytp
1242!
1243subroutine  ymist(IMR,JNP,j1,P,DC,ID)
1244  implicit none
1245  integer :: IMR,JNP,j1,ID
1246  real,parameter :: R24 = 1./24.
1247  real :: P(IMR,JNP),DC(IMR,JNP)
1248  integer :: iimh,jmr,ijm3,imh,i
1249  real :: pmax,pmin,tmp
1250  !
1251  IMH = IMR / 2
1252  JMR = JNP - 1
1253  IJM3 = IMR*(JMR-3)
1254  !
1255  IF(ID.EQ.2) THEN
1256  do i=1,IMR*(JMR-1)
1257  tmp = 0.25*(p(i,3) - p(i,1))
1258  Pmax = max(p(i,1),p(i,2),p(i,3)) - p(i,2)
1259  Pmin = p(i,2) - min(p(i,1),p(i,2),p(i,3))
1260  DC(i,2) = sign(min(abs(tmp),Pmin,Pmax),tmp)
1261  end do
1262  ELSE
1263  do i=1,IMH
1264  ! J=2
1265  tmp = (8.*(p(i,3) - p(i,1)) + p(i+IMH,2) - p(i,4))*R24
1266  Pmax = max(p(i,1),p(i,2),p(i,3)) - p(i,2)
1267  Pmin = p(i,2) - min(p(i,1),p(i,2),p(i,3))
1268  DC(i,2) = sign(min(abs(tmp),Pmin,Pmax),tmp)
1269  ! J=JMR
1270  tmp=(8.*(p(i,JNP)-p(i,JMR-1))+p(i,JMR-2)-p(i+IMH,JMR))*R24
1271  Pmax = max(p(i,JMR-1),p(i,JMR),p(i,JNP)) - p(i,JMR)
1272  Pmin = p(i,JMR) - min(p(i,JMR-1),p(i,JMR),p(i,JNP))
1273  DC(i,JMR) = sign(min(abs(tmp),Pmin,Pmax),tmp)
1274  end do
1275  do i=IMH+1,IMR
1276  ! J=2
1277  tmp = (8.*(p(i,3) - p(i,1)) + p(i-IMH,2) - p(i,4))*R24
1278  Pmax = max(p(i,1),p(i,2),p(i,3)) - p(i,2)
1279  Pmin = p(i,2) - min(p(i,1),p(i,2),p(i,3))
1280  DC(i,2) = sign(min(abs(tmp),Pmin,Pmax),tmp)
1281  ! J=JMR
1282  tmp=(8.*(p(i,JNP)-p(i,JMR-1))+p(i,JMR-2)-p(i-IMH,JMR))*R24
1283  Pmax = max(p(i,JMR-1),p(i,JMR),p(i,JNP)) - p(i,JMR)
1284  Pmin = p(i,JMR) - min(p(i,JMR-1),p(i,JMR),p(i,JNP))
1285  DC(i,JMR) = sign(min(abs(tmp),Pmin,Pmax),tmp)
1286  end do
1287  !
1288  do i=1,IJM3
1289  tmp = (8.*(p(i,4) - p(i,2)) + p(i,1) - p(i,5))*R24
1290  Pmax = max(p(i,2),p(i,3),p(i,4)) - p(i,3)
1291  Pmin = p(i,3) - min(p(i,2),p(i,3),p(i,4))
1292  DC(i,3) = sign(min(abs(tmp),Pmin,Pmax),tmp)
1293  end do
1294  ENDIF
1295  !
1296  if(j1.ne.2) then
1297  do i=1,IMR
1298  DC(i,1) = 0.
1299  DC(i,JNP) = 0.
1300  enddo
1301  else
1302  ! Determine slopes in polar caps for scalars!
1303  !
1304  do i=1,IMH
1305  ! South
1306  tmp = 0.25*(p(i,2) - p(i+imh,2))
1307  Pmax = max(p(i,2),p(i,1), p(i+imh,2)) - p(i,1)
1308  Pmin = p(i,1) - min(p(i,2),p(i,1), p(i+imh,2))
1309  DC(i,1)=sign(min(abs(tmp),Pmax,Pmin),tmp)
1310  ! North.
1311  tmp = 0.25*(p(i+imh,JMR) - p(i,JMR))
1312  Pmax = max(p(i+imh,JMR),p(i,jnp), p(i,JMR)) - p(i,JNP)
1313  Pmin = p(i,JNP) - min(p(i+imh,JMR),p(i,jnp), p(i,JMR))
1314  DC(i,JNP) = sign(min(abs(tmp),Pmax,pmin),tmp)
1315  end do
1316  !
1317  do i=imh+1,IMR
1318  DC(i,  1) =  - DC(i-imh,  1)
1319  DC(i,JNP) =  - DC(i-imh,JNP)
1320  end do
1321  endif
1322  return
1323end subroutine ymist
1324!
1325subroutine fyppm(VC,P,DC,flux,IMR,JNP,j1,j2,A6,AR,AL,JORD)
1326  implicit none
1327  integer :: IMR,JNP,j1,j2,JORD
1328  real,parameter :: R3 = 1./3., R23 = 2./3.
1329  real :: VC(IMR,*),flux(IMR,*),P(IMR,*),DC(IMR,*)
1330  ! Local work arrays.
1331  real :: AR(IMR,JNP),AL(IMR,JNP),A6(IMR,JNP)
1332  integer :: LMT,i
1333  integer :: IMH,JMR,j11,IMJM1,len
1334   ! logical first
1335   ! data first /.true./
1336   ! SAVE LMT
1337  !
1338  IMH = IMR / 2
1339  JMR = JNP - 1
1340  j11 = j1-1
1341  IMJM1 = IMR*(J2-J1+2)
1342  len   = IMR*(J2-J1+3)
1343   ! if(first) then
1344   ! IF(JORD.LE.0) then
1345   !       if(JMR.GE.90) then
1346   !             LMT = 0
1347   !       elseif(JMR.GE.45) then
1348   !             LMT = 1
1349   !       else
1350   !             LMT = 2
1351   !       endif
1352   ! else
1353   !       LMT = JORD - 3
1354   ! endif
1355  !
1356  !  first = .false.
1357  !  endif
1358  !
1359  ! modifs pour pouvoir choisir plusieurs schemas PPM
1360  LMT = JORD - 3
1361  !
1362  DO i=1,IMR*JMR
1363  AL(i,2) = 0.5*(p(i,1)+p(i,2)) + (DC(i,1) - DC(i,2))*R3
1364  AR(i,1) = AL(i,2)
1365  END DO
1366  !
1367  !Poles:
1368  !
1369  DO i=1,IMH
1370  AL(i,1) = AL(i+IMH,2)
1371  AL(i+IMH,1) = AL(i,2)
1372  !
1373  AR(i,JNP) = AR(i+IMH,JMR)
1374  AR(i+IMH,JNP) = AR(i,JMR)
1375  ENDDO
1376
1377  !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1378  !   Rajout pour LMDZ.3.3
1379  !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1380  AR(IMR,1)=AL(1,1)
1381  AR(IMR,JNP)=AL(1,JNP)
1382  !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1383
1384
1385  do i=1,len
1386    A6(i,j11) = 3.*(p(i,j11)+p(i,j11)  - (AL(i,j11)+AR(i,j11)))
1387  end do
1388  !
1389  if(LMT.le.2) call lmtppm(DC(1,j11),A6(1,j11),AR(1,j11) &
1390        ,AL(1,j11),P(1,j11),len,LMT)
1391  !
1392
1393  DO i=1,IMJM1
1394  IF(VC(i,j1).GT.0.) then
1395  flux(i,j1) = AR(i,j11) + 0.5*VC(i,j1)*(AL(i,j11) - AR(i,j11) + &
1396        A6(i,j11)*(1.-R23*VC(i,j1)) )
1397  else
1398  flux(i,j1) = AL(i,j1) - 0.5*VC(i,j1)*(AR(i,j1) - AL(i,j1) + &
1399        A6(i,j1)*(1.+R23*VC(i,j1)))
1400  endif
1401  END DO
1402  return
1403end subroutine fyppm
1404!
1405  subroutine yadv(IMR,JNP,j1,j2,p,VA,ady,wk,IAD)
1406    implicit none
1407    integer :: IMR,JNP,j1,j2,IAD
1408    REAL :: p(IMR,JNP),ady(IMR,JNP),VA(IMR,JNP)
1409    REAL :: WK(IMR,-1:JNP+2)
1410    INTEGER :: JMR,IMH,i,j,jp
1411    REAL :: rv,a1,b1,sum1,sum2
1412  !
1413    JMR = JNP-1
1414    IMH = IMR/2
1415    do j=1,JNP
1416    do i=1,IMR
1417    wk(i,j) = p(i,j)
1418    enddo
1419    enddo
1420  ! Poles:
1421    do i=1,IMH
1422    wk(i,   -1) = p(i+IMH,3)
1423    wk(i+IMH,-1) = p(i,3)
1424    wk(i,    0) = p(i+IMH,2)
1425    wk(i+IMH,0) = p(i,2)
1426    wk(i,JNP+1) = p(i+IMH,JMR)
1427    wk(i+IMH,JNP+1) = p(i,JMR)
1428    wk(i,JNP+2) = p(i+IMH,JNP-2)
1429    wk(i+IMH,JNP+2) = p(i,JNP-2)
1430    enddo
1431     ! write(*,*) 'toto 1'
1432  ! --------------------------------
1433  IF(IAD.eq.2) then
1434  do j=j1-1,j2+1
1435  do i=1,IMR
1436   ! write(*,*) 'avt NINT','i=',i,'j=',j
1437  JP = NINT(VA(i,j))
1438  rv = JP - VA(i,j)
1439   ! write(*,*) 'VA=',VA(i,j), 'JP1=',JP,'rv=',rv
1440  JP = j - JP
1441   ! write(*,*) 'JP2=',JP
1442  a1 = 0.5*(wk(i,jp+1)+wk(i,jp-1)) - wk(i,jp)
1443  b1 = 0.5*(wk(i,jp+1)-wk(i,jp-1))
1444   ! write(*,*) 'a1=',a1,'b1=',b1
1445  ady(i,j) = wk(i,jp) + rv*(a1*rv + b1) - wk(i,j)
1446  enddo
1447  enddo
1448   ! write(*,*) 'toto 2'
1449  !
1450  ELSEIF(IAD.eq.1) then
1451    do j=j1-1,j2+1
1452  do i=1,imr
1453  JP = REAL(j)-VA(i,j)
1454  ady(i,j) = VA(i,j)*(wk(i,jp)-wk(i,jp+1))
1455  enddo
1456  enddo
1457  ENDIF
1458  !
1459    if(j1.ne.2) then
1460    sum1 = 0.
1461    sum2 = 0.
1462  do i=1,imr
1463  sum1 = sum1 + ady(i,2)
1464  sum2 = sum2 + ady(i,JMR)
1465  enddo
1466    sum1 = sum1 / IMR
1467    sum2 = sum2 / IMR
1468  !
1469  do i=1,imr
1470  ady(i,  2) =  sum1
1471  ady(i,JMR) =  sum2
1472  ady(i,  1) =  sum1
1473  ady(i,JNP) =  sum2
1474  enddo
1475    else
1476  ! Poles:
1477    sum1 = 0.
1478    sum2 = 0.
1479  do i=1,imr
1480  sum1 = sum1 + ady(i,1)
1481  sum2 = sum2 + ady(i,JNP)
1482  enddo
1483    sum1 = sum1 / IMR
1484    sum2 = sum2 / IMR
1485  !
1486  do i=1,imr
1487  ady(i,  1) =  sum1
1488  ady(i,JNP) =  sum2
1489  enddo
1490    endif
1491  !
1492    return
1493end subroutine yadv
1494!
1495  subroutine xadv(IMR,JNP,j1,j2,p,UA,JS,JN,IML,adx,IAD)
1496    implicit none
1497    INTEGER :: IMR,JNP,j1,j2,JS,JN,IML,IAD
1498    REAL :: p(IMR,JNP),adx(IMR,JNP),qtmp(-IMR:IMR+IMR),UA(IMR,JNP)
1499    INTEGER :: JMR,j,i,ip,iu,iiu
1500    REAL :: ru,a1,b1
1501  !
1502    JMR = JNP-1
1503  do j=j1,j2
1504  if(J.GT.JS  .and. J.LT.JN) GO TO 1309
1505  !
1506  do i=1,IMR
1507  qtmp(i) = p(i,j)
1508  enddo
1509  !
1510  do i=-IML,0
1511  qtmp(i)       = p(IMR+i,j)
1512  qtmp(IMR+1-i) = p(1-i,j)
1513  enddo
1514  !
1515  IF(IAD.eq.2) THEN
1516  DO i=1,IMR
1517  IP = NINT(UA(i,j))
1518  ru = IP - UA(i,j)
1519  IP = i - IP
1520  a1 = 0.5*(qtmp(ip+1)+qtmp(ip-1)) - qtmp(ip)
1521  b1 = 0.5*(qtmp(ip+1)-qtmp(ip-1))
1522  adx(i,j) = qtmp(ip) + ru*(a1*ru + b1)
1523  enddo
1524  ELSEIF(IAD.eq.1) then
1525  DO i=1,IMR
1526  iu = UA(i,j)
1527  ru = UA(i,j) - iu
1528  iiu = i-iu
1529  if(UA(i,j).GE.0.) then
1530  adx(i,j) = qtmp(iiu)+ru*(qtmp(iiu-1)-qtmp(iiu))
1531  else
1532  adx(i,j) = qtmp(iiu)+ru*(qtmp(iiu)-qtmp(iiu+1))
1533  endif
1534  enddo
1535  ENDIF
1536  !
1537  do i=1,IMR
1538  adx(i,j) = adx(i,j) - p(i,j)
1539  enddo
15401309 continue
1541  end do
1542  !
1543  ! Eulerian upwind
1544  !
1545  do j=JS+1,JN-1
1546  !
1547  do i=1,IMR
1548  qtmp(i) = p(i,j)
1549  enddo
1550  !
1551  qtmp(0)     = p(IMR,J)
1552  qtmp(IMR+1) = p(1,J)
1553  !
1554  IF(IAD.eq.2) THEN
1555  qtmp(-1)     = p(IMR-1,J)
1556  qtmp(IMR+2) = p(2,J)
1557  do i=1,imr
1558  IP = NINT(UA(i,j))
1559  ru = IP - UA(i,j)
1560  IP = i - IP
1561  a1 = 0.5*(qtmp(ip+1)+qtmp(ip-1)) - qtmp(ip)
1562  b1 = 0.5*(qtmp(ip+1)-qtmp(ip-1))
1563  adx(i,j) = qtmp(ip)- p(i,j) + ru*(a1*ru + b1)
1564  enddo
1565  ELSEIF(IAD.eq.1) then
1566  ! 1st order
1567  DO i=1,IMR
1568  IP = i - UA(i,j)
1569  adx(i,j) = UA(i,j)*(qtmp(ip)-qtmp(ip+1))
1570  enddo
1571  ENDIF
1572  enddo
1573  !
1574    if(j1.ne.2) then
1575  do i=1,IMR
1576  adx(i,  2) = 0.
1577  adx(i,JMR) = 0.
1578  enddo
1579    endif
1580  ! set cross term due to x-adv at the poles to zero.
1581  do i=1,IMR
1582  adx(i,  1) = 0.
1583  adx(i,JNP) = 0.
1584  enddo
1585    return
1586end subroutine xadv
1587!
1588subroutine lmtppm(DC,A6,AR,AL,P,IM,LMT)
1589  implicit none
1590  !
1591  ! A6 =  CURVATURE OF THE TEST PARABOLA
1592  ! AR =  RIGHT EDGE VALUE OF THE TEST PARABOLA
1593  ! AL =  LEFT  EDGE VALUE OF THE TEST PARABOLA
1594  ! DC =  0.5 * MISMATCH
1595  ! P  =  CELL-AVERAGED VALUE
1596  ! IM =  VECTOR LENGTH
1597  !
1598  ! OPTIONS:
1599  !
1600  ! LMT = 0: FULL MONOTONICITY
1601  ! LMT = 1: SEMI-MONOTONIC CONSTRAINT (NO UNDERSHOOTS)
1602  ! LMT = 2: POSITIVE-DEFINITE CONSTRAINT
1603  !
1604  real,parameter :: R12 = 1./12.
1605  real :: A6(IM),AR(IM),AL(IM),P(IM),DC(IM)
1606  integer :: IM,LMT
1607  INTEGER :: i
1608  REAL :: da1,da2,a6da,fmin
1609  !
1610  if(LMT.eq.0) then
1611  ! Full constraint
1612  do i=1,IM
1613  if(DC(i).eq.0.) then
1614        AR(i) = p(i)
1615        AL(i) = p(i)
1616        A6(i) = 0.
1617  else
1618  da1  = AR(i) - AL(i)
1619  da2  = da1**2
1620  A6DA = A6(i)*da1
1621  if(A6DA .lt. -da2) then
1622        A6(i) = 3.*(AL(i)-p(i))
1623        AR(i) = AL(i) - A6(i)
1624  elseif(A6DA .gt. da2) then
1625        A6(i) = 3.*(AR(i)-p(i))
1626        AL(i) = AR(i) - A6(i)
1627  endif
1628  endif
1629  end do
1630  elseif(LMT.eq.1) then
1631  ! Semi-monotonic constraint
1632  do i=1,IM
1633  if(abs(AR(i)-AL(i)) .GE. -A6(i)) go to 150
1634  if(p(i).lt.AR(i) .and. p(i).lt.AL(i)) then
1635        AR(i) = p(i)
1636        AL(i) = p(i)
1637        A6(i) = 0.
1638  elseif(AR(i) .gt. AL(i)) then
1639        A6(i) = 3.*(AL(i)-p(i))
1640        AR(i) = AL(i) - A6(i)
1641  else
1642        A6(i) = 3.*(AR(i)-p(i))
1643        AL(i) = AR(i) - A6(i)
1644  endif
1645150 continue
1646  end do
1647  elseif(LMT.eq.2) then
1648  do i=1,IM
1649  if(abs(AR(i)-AL(i)) .GE. -A6(i)) go to 250
1650  fmin = p(i) + 0.25*(AR(i)-AL(i))**2/A6(i) + A6(i)*R12
1651  if(fmin.ge.0.) go to 250
1652  if(p(i).lt.AR(i) .and. p(i).lt.AL(i)) then
1653        AR(i) = p(i)
1654        AL(i) = p(i)
1655        A6(i) = 0.
1656  elseif(AR(i) .gt. AL(i)) then
1657        A6(i) = 3.*(AL(i)-p(i))
1658        AR(i) = AL(i) - A6(i)
1659  else
1660        A6(i) = 3.*(AR(i)-p(i))
1661        AL(i) = AR(i) - A6(i)
1662  endif
1663250 continue
1664  end do
1665  endif
1666  return
1667end subroutine lmtppm
1668!
1669subroutine A2C(U,V,IMR,JMR,j1,j2,CRX,CRY,dtdx5,DTDY5)
1670  implicit none
1671  integer :: IMR,JMR,j1,j2
1672  real :: U(IMR,*),V(IMR,*),CRX(IMR,*),CRY(IMR,*),DTDX5(*),DTDY5
1673  integer :: i,j
1674  !
1675  do j=j1,j2
1676  do i=2,IMR
1677    CRX(i,J) = dtdx5(j)*(U(i,j)+U(i-1,j))
1678  end do
1679  end do
1680  !
1681  do j=j1,j2
1682    CRX(1,J) = dtdx5(j)*(U(1,j)+U(IMR,j))
1683  end do
1684  !
1685  do i=1,IMR*JMR
1686    CRY(i,2) = DTDY5*(V(i,2)+V(i,1))
1687  end do
1688  return
1689end subroutine a2c
1690!
1691subroutine cosa(cosp,cose,JNP,PI,DP)
1692  implicit none
1693  integer :: JNP
1694  real :: cosp(*),cose(*),PI,DP
1695  integer :: JMR,j,jeq
1696  real :: ph5
1697  JMR = JNP-1
1698  do j=2,JNP
1699    ph5  =  -0.5*PI + (REAL(J-1)-0.5)*DP
1700    cose(j) = cos(ph5)
1701  end do
1702  !
1703  JEQ = (JNP+1) / 2
1704  if(JMR .eq. 2*(JMR/2) ) then
1705  do j=JNP, JEQ+1, -1
1706   cose(j) =  cose(JNP+2-j)
1707  enddo
1708  else
1709  ! cell edge at equator.
1710   cose(JEQ+1) =  1.
1711  do j=JNP, JEQ+2, -1
1712   cose(j) =  cose(JNP+2-j)
1713   enddo
1714  endif
1715  !
1716  do j=2,JMR
1717    cosp(j) = 0.5*(cose(j)+cose(j+1))
1718  end do
1719  cosp(1) = 0.
1720  cosp(JNP) = 0.
1721  return
1722end subroutine cosa
1723!
1724subroutine cosc(cosp,cose,JNP,PI,DP)
1725  implicit none
1726  integer :: JNP
1727  real :: cosp(*),cose(*),PI,DP
1728  real :: phi
1729  integer :: j
1730  !
1731  phi = -0.5*PI
1732  do j=2,JNP-1
1733  phi  =  phi + DP
1734    cosp(j) = cos(phi)
1735  end do
1736    cosp(  1) = 0.
1737    cosp(JNP) = 0.
1738  !
1739  do j=2,JNP
1740    cose(j) = 0.5*(cosp(j)+cosp(j-1))
1741  end do
1742  !
1743  do j=2,JNP-1
1744   cosp(j) = 0.5*(cose(j)+cose(j+1))
1745  end do
1746  return
1747end subroutine cosc
1748!
1749SUBROUTINE qckxyz (Q,qtmp,IMR,JNP,NLAY,j1,j2,cosp,acosp, &
1750        cross,IC,NSTEP)
1751  !
1752  real,parameter :: tiny = 1.E-60
1753  INTEGER :: IMR,JNP,NLAY,j1,j2,IC,NSTEP
1754  REAL :: Q(IMR,JNP,NLAY),qtmp(IMR,JNP),cosp(*),acosp(*)
1755  logical :: cross
1756  INTEGER :: NLAYM1,len,ip,L,icr,ipy,ipx,i
1757  real :: qup,qly,dup,sum
1758  !
1759  NLAYM1 = NLAY-1
1760  len = IMR*(j2-j1+1)
1761  ip = 0
1762  !
1763  ! Top layer
1764  L = 1
1765    icr = 1
1766  call filns(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,ipy,tiny)
1767  if(ipy.eq.0) goto 50
1768  call filew(q(1,1,L),qtmp,IMR,JNP,j1,j2,ipx,tiny)
1769  if(ipx.eq.0) goto 50
1770  !
1771  if(cross) then
1772  call filcr(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,icr,tiny)
1773  endif
1774  if(icr.eq.0) goto 50
1775  !
1776  ! Vertical filling...
1777  do i=1,len
1778  IF( Q(i,j1,1).LT.0.) THEN
1779  ip = ip + 1
1780      Q(i,j1,2) = Q(i,j1,2) + Q(i,j1,1)
1781      Q(i,j1,1) = 0.
1782  endif
1783  enddo
1784  !
178550   continue
1786  DO L = 2,NLAYM1
1787  icr = 1
1788  !
1789  call filns(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,ipy,tiny)
1790  if(ipy.eq.0) goto 225
1791  call filew(q(1,1,L),qtmp,IMR,JNP,j1,j2,ipx,tiny)
1792  if(ipx.eq.0) go to 225
1793  if(cross) then
1794  call filcr(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,icr,tiny)
1795  endif
1796  if(icr.eq.0) goto 225
1797  !
1798  do i=1,len
1799  IF( Q(I,j1,L).LT.0.) THEN
1800  !
1801  ip = ip + 1
1802  ! From above
1803      qup =  Q(I,j1,L-1)
1804      qly = -Q(I,j1,L)
1805      dup  = min(qly,qup)
1806      Q(I,j1,L-1) = qup - dup
1807      Q(I,j1,L  ) = dup-qly
1808  ! Below
1809      Q(I,j1,L+1) = Q(I,j1,L+1) + Q(I,j1,L)
1810      Q(I,j1,L)   = 0.
1811  ENDIF
1812  ENDDO
1813225 CONTINUE
1814  END DO
1815  !
1816  ! BOTTOM LAYER
1817  sum = 0.
1818  L = NLAY
1819  !
1820  call filns(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,ipy,tiny)
1821  if(ipy.eq.0) goto 911
1822  call filew(q(1,1,L),qtmp,IMR,JNP,j1,j2,ipx,tiny)
1823  if(ipx.eq.0) goto 911
1824  !
1825  call filcr(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,icr,tiny)
1826  if(icr.eq.0) goto 911
1827  !
1828  DO  I=1,len
1829  IF( Q(I,j1,L).LT.0.) THEN
1830  ip = ip + 1
1831  !
1832  ! From above
1833  !
1834      qup = Q(I,j1,NLAYM1)
1835      qly = -Q(I,j1,L)
1836      dup = min(qly,qup)
1837      Q(I,j1,NLAYM1) = qup - dup
1838  ! From "below" the surface.
1839      sum = sum + qly-dup
1840      Q(I,j1,L) = 0.
1841   ENDIF
1842  ENDDO
1843  !
1844911   continue
1845  !
1846  if(ip.gt.IMR) then
1847  write(6,*) 'IC=',IC,' STEP=',NSTEP, &
1848        ' Vertical filling pts=',ip
1849  endif
1850  !
1851  if(sum.gt.1.e-25) then
1852  write(6,*) IC,NSTEP,' Mass source from the ground=',sum
1853  endif
1854  RETURN
1855END SUBROUTINE qckxyz
1856!
1857subroutine filcr(q,IMR,JNP,j1,j2,cosp,acosp,icr,tiny)
1858  implicit none
1859  integer :: IMR,JNP,j1,j2,icr
1860  real :: q(IMR,*),cosp(*),acosp(*),tiny
1861  integer :: i,j
1862  real :: dq,dn,d0,d1,ds,d2
1863  icr = 0
1864  do j=j1+1,j2-1
1865  DO i=1,IMR-1
1866  IF(q(i,j).LT.0.) THEN
1867  icr =  1
1868  dq  = - q(i,j)*cosp(j)
1869  ! N-E
1870  dn = q(i+1,j+1)*cosp(j+1)
1871  d0 = max(0.,dn)
1872  d1 = min(dq,d0)
1873  q(i+1,j+1) = (dn - d1)*acosp(j+1)
1874  dq = dq - d1
1875  ! S-E
1876  ds = q(i+1,j-1)*cosp(j-1)
1877  d0 = max(0.,ds)
1878  d2 = min(dq,d0)
1879  q(i+1,j-1) = (ds - d2)*acosp(j-1)
1880  q(i,j) = (d2 - dq)*acosp(j) + tiny
1881  endif
188250  CONTINUE
1883  END DO
1884  if(icr.eq.0 .and. q(IMR,j).ge.0.) goto 65
1885  DO i=2,IMR
1886  IF(q(i,j).LT.0.) THEN
1887  icr =  1
1888  dq  = - q(i,j)*cosp(j)
1889  ! N-W
1890  dn = q(i-1,j+1)*cosp(j+1)
1891  d0 = max(0.,dn)
1892  d1 = min(dq,d0)
1893  q(i-1,j+1) = (dn - d1)*acosp(j+1)
1894  dq = dq - d1
1895  ! S-W
1896  ds = q(i-1,j-1)*cosp(j-1)
1897  d0 = max(0.,ds)
1898  d2 = min(dq,d0)
1899  q(i-1,j-1) = (ds - d2)*acosp(j-1)
1900  q(i,j) = (d2 - dq)*acosp(j) + tiny
1901  endif
1902  END DO
1903  ! *****************************************
1904  ! i=1
1905  i=1
1906  IF(q(i,j).LT.0.) THEN
1907  icr =  1
1908  dq  = - q(i,j)*cosp(j)
1909  ! N-W
1910  dn = q(IMR,j+1)*cosp(j+1)
1911  d0 = max(0.,dn)
1912  d1 = min(dq,d0)
1913  q(IMR,j+1) = (dn - d1)*acosp(j+1)
1914  dq = dq - d1
1915  ! S-W
1916  ds = q(IMR,j-1)*cosp(j-1)
1917  d0 = max(0.,ds)
1918  d2 = min(dq,d0)
1919  q(IMR,j-1) = (ds - d2)*acosp(j-1)
1920  q(i,j) = (d2 - dq)*acosp(j) + tiny
1921  endif
1922  ! *****************************************
1923  ! i=IMR
1924  i=IMR
1925  IF(q(i,j).LT.0.) THEN
1926  icr =  1
1927  dq  = - q(i,j)*cosp(j)
1928  ! N-E
1929  dn = q(1,j+1)*cosp(j+1)
1930  d0 = max(0.,dn)
1931  d1 = min(dq,d0)
1932  q(1,j+1) = (dn - d1)*acosp(j+1)
1933  dq = dq - d1
1934  ! S-E
1935  ds = q(1,j-1)*cosp(j-1)
1936  d0 = max(0.,ds)
1937  d2 = min(dq,d0)
1938  q(1,j-1) = (ds - d2)*acosp(j-1)
1939  q(i,j) = (d2 - dq)*acosp(j) + tiny
1940  endif
1941  ! *****************************************
194265  continue
1943  end do
1944  !
1945  do i=1,IMR
1946  if(q(i,j1).lt.0. .or. q(i,j2).lt.0.) then
1947  icr = 1
1948  goto 80
1949  endif
1950  enddo
1951  !
195280   continue
1953  !
1954  if(q(1,1).lt.0. .or. q(1,jnp).lt.0.) then
1955  icr = 1
1956  endif
1957  !
1958  return
1959end subroutine filcr
1960!
1961subroutine filns(q,IMR,JNP,j1,j2,cosp,acosp,ipy,tiny)
1962  implicit none
1963  integer :: IMR,JNP,j1,j2,ipy
1964  real :: q(IMR,*),cosp(*),acosp(*),tiny
1965  real :: DP,CAP1,dq,dn,d0,d1,ds,d2
1966  INTEGER :: i,j
1967   ! logical first
1968   ! data first /.true./
1969   ! save cap1
1970  !
1971  !  if(first) then
1972  DP = 4.*ATAN(1.)/REAL(JNP-1)
1973  CAP1 = IMR*(1.-COS((j1-1.5)*DP))/DP
1974   ! first = .false.
1975   ! endif
1976  !
1977  ipy = 0
1978  do j=j1+1,j2-1
1979  DO i=1,IMR
1980  IF(q(i,j).LT.0.) THEN
1981  ipy =  1
1982  dq  = - q(i,j)*cosp(j)
1983  ! North
1984  dn = q(i,j+1)*cosp(j+1)
1985  d0 = max(0.,dn)
1986  d1 = min(dq,d0)
1987  q(i,j+1) = (dn - d1)*acosp(j+1)
1988  dq = dq - d1
1989  ! South
1990  ds = q(i,j-1)*cosp(j-1)
1991  d0 = max(0.,ds)
1992  d2 = min(dq,d0)
1993  q(i,j-1) = (ds - d2)*acosp(j-1)
1994  q(i,j) = (d2 - dq)*acosp(j) + tiny
1995  endif
1996  END DO
1997  end do
1998  !
1999  do i=1,imr
2000  IF(q(i,j1).LT.0.) THEN
2001  ipy =  1
2002  dq  = - q(i,j1)*cosp(j1)
2003  ! North
2004  dn = q(i,j1+1)*cosp(j1+1)
2005  d0 = max(0.,dn)
2006  d1 = min(dq,d0)
2007  q(i,j1+1) = (dn - d1)*acosp(j1+1)
2008  q(i,j1) = (d1 - dq)*acosp(j1) + tiny
2009  endif
2010  enddo
2011  !
2012  j = j2
2013  do i=1,imr
2014  IF(q(i,j).LT.0.) THEN
2015  ipy =  1
2016  dq  = - q(i,j)*cosp(j)
2017  ! South
2018  ds = q(i,j-1)*cosp(j-1)
2019  d0 = max(0.,ds)
2020  d2 = min(dq,d0)
2021  q(i,j-1) = (ds - d2)*acosp(j-1)
2022  q(i,j) = (d2 - dq)*acosp(j) + tiny
2023  endif
2024  enddo
2025  !
2026  ! Check Poles.
2027  if(q(1,1).lt.0.) then
2028  dq = q(1,1)*cap1/REAL(IMR)*acosp(j1)
2029  do i=1,imr
2030  q(i,1) = 0.
2031  q(i,j1) = q(i,j1) + dq
2032  if(q(i,j1).lt.0.) ipy = 1
2033  enddo
2034  endif
2035  !
2036  if(q(1,JNP).lt.0.) then
2037  dq = q(1,JNP)*cap1/REAL(IMR)*acosp(j2)
2038  do i=1,imr
2039  q(i,JNP) = 0.
2040  q(i,j2) = q(i,j2) + dq
2041  if(q(i,j2).lt.0.) ipy = 1
2042  enddo
2043  endif
2044  !
2045  return
2046end subroutine filns
2047!
2048subroutine filew(q,qtmp,IMR,JNP,j1,j2,ipx,tiny)
2049  implicit none
2050  integer :: IMR,JNP,j1,j2,ipx
2051  real :: q(IMR,*),qtmp(JNP,IMR),tiny
2052  integer :: i,j
2053  real :: d0,d1,d2
2054  !
2055  ipx = 0
2056  ! Copy & swap direction for vectorization.
2057  do i=1,imr
2058  do j=j1,j2
2059    qtmp(j,i) = q(i,j)
2060  end do
2061  end do
2062  !
2063  do i=2,imr-1
2064  do j=j1,j2
2065  if(qtmp(j,i).lt.0.) then
2066  ipx =  1
2067  ! west
2068  d0 = max(0.,qtmp(j,i-1))
2069  d1 = min(-qtmp(j,i),d0)
2070  qtmp(j,i-1) = qtmp(j,i-1) - d1
2071  qtmp(j,i) = qtmp(j,i) + d1
2072  ! east
2073  d0 = max(0.,qtmp(j,i+1))
2074  d2 = min(-qtmp(j,i),d0)
2075  qtmp(j,i+1) = qtmp(j,i+1) - d2
2076  qtmp(j,i) = qtmp(j,i) + d2 + tiny
2077  endif
2078  end do
2079  end do
2080  !
2081  i=1
2082  do j=j1,j2
2083  if(qtmp(j,i).lt.0.) then
2084  ipx =  1
2085  ! west
2086  d0 = max(0.,qtmp(j,imr))
2087  d1 = min(-qtmp(j,i),d0)
2088  qtmp(j,imr) = qtmp(j,imr) - d1
2089  qtmp(j,i) = qtmp(j,i) + d1
2090  ! east
2091  d0 = max(0.,qtmp(j,i+1))
2092  d2 = min(-qtmp(j,i),d0)
2093  qtmp(j,i+1) = qtmp(j,i+1) - d2
2094  !
2095  qtmp(j,i) = qtmp(j,i) + d2 + tiny
2096  endif
209765  continue
2098  end do
2099  i=IMR
2100  do j=j1,j2
2101  if(qtmp(j,i).lt.0.) then
2102  ipx =  1
2103  ! west
2104  d0 = max(0.,qtmp(j,i-1))
2105  d1 = min(-qtmp(j,i),d0)
2106  qtmp(j,i-1) = qtmp(j,i-1) - d1
2107  qtmp(j,i) = qtmp(j,i) + d1
2108  ! east
2109  d0 = max(0.,qtmp(j,1))
2110  d2 = min(-qtmp(j,i),d0)
2111  qtmp(j,1) = qtmp(j,1) - d2
2112  !
2113  qtmp(j,i) = qtmp(j,i) + d2 + tiny
2114  endif
2115  end do
2116  !
2117  if(ipx.ne.0) then
2118  do j=j1,j2
2119  do i=1,imr
2120    q(i,j) = qtmp(j,i)
2121  end do
2122  end do
2123  else
2124  !
2125  ! Poles.
2126  if(q(1,1).lt.0.OR. q(1,JNP).lt.0.) ipx = 1
2127  endif
2128  return
2129end subroutine filew
2130!
2131subroutine zflip(q,im,km,nc)
2132  implicit none
2133  ! This routine flip the array q (in the vertical).
2134  integer :: im,km,nc
2135  real :: q(im,km,nc)
2136  ! local dynamic array
2137  real :: qtmp(im,km)
2138  integer :: IC,k,i
2139  !
2140  do IC = 1, nc
2141  !
2142  do k=1,km
2143  do i=1,im
2144  qtmp(i,k) = q(i,km+1-k,IC)
2145  end do
2146  end do
2147  !
2148  do i=1,im*km
2149    q(i,1,IC) = qtmp(i,1)
21502000 continue
2151  end do
2152  end do
2153  return
2154end subroutine zflip
Note: See TracBrowser for help on using the repository browser.