Ignore:
Timestamp:
Jul 24, 2024, 6:46:45 PM (4 months ago)
Author:
abarral
Message:

enforce PRIVATE by default in several modules, expose PUBLIC as needed
move eigen.f90 to obsolete/
(lint) aslong the way

Location:
LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat/calfis.f90

    r5118 r5119  
    1 
    21! $Id$
    32
     
    3635  USE comvert_mod, ONLY: preff, presnivs
    3736  USE lmdz_iniprint, ONLY: lunout, prt_level
     37  USE lmdz_ssum_scopy, ONLY: scopy
    3838
    3939  IMPLICIT NONE
     
    9595
    9696  INTEGER :: ngridmx
    97   PARAMETER( ngridmx = 2+(jjm-1)*iim - 1/jjm   )
     97  PARAMETER(ngridmx = 2 + (jjm - 1) * iim - 1 / jjm)
    9898
    9999  include "comgeom2.h"
     
    101101  !    Arguments :
    102102  !    -----------
    103   LOGICAL,INTENT(IN) :: lafin ! .TRUE. for the very last CALL to physics
    104   REAL,INTENT(IN):: jD_cur, jH_cur
    105   REAL,INTENT(IN) :: pvcov(iip1,jjm,llm) ! covariant meridional velocity
    106   REAL,INTENT(IN) :: pucov(iip1,jjp1,llm) ! covariant zonal velocity
    107   REAL,INTENT(IN) :: pteta(iip1,jjp1,llm) ! potential temperature
    108   REAL,INTENT(IN) :: pmasse(iip1,jjp1,llm) ! mass in each cell ! not used
    109   REAL,INTENT(IN) :: pq(iip1,jjp1,llm,nqtot) ! tracers
    110   REAL,INTENT(IN) :: pphis(iip1,jjp1) ! surface geopotential
    111   REAL,INTENT(IN) :: pphi(iip1,jjp1,llm) ! geopotential
    112 
    113   REAL,INTENT(IN) :: pdvcov(iip1,jjm,llm) ! dynamical tendency on vcov
    114   REAL,INTENT(IN) :: pducov(iip1,jjp1,llm) ! dynamical tendency on ucov
    115   REAL,INTENT(IN) :: pdteta(iip1,jjp1,llm) ! dynamical tendency on teta
     103  LOGICAL, INTENT(IN) :: lafin ! .TRUE. for the very last CALL to physics
     104  REAL, INTENT(IN) :: jD_cur, jH_cur
     105  REAL, INTENT(IN) :: pvcov(iip1, jjm, llm) ! covariant meridional velocity
     106  REAL, INTENT(IN) :: pucov(iip1, jjp1, llm) ! covariant zonal velocity
     107  REAL, INTENT(IN) :: pteta(iip1, jjp1, llm) ! potential temperature
     108  REAL, INTENT(IN) :: pmasse(iip1, jjp1, llm) ! mass in each cell ! not used
     109  REAL, INTENT(IN) :: pq(iip1, jjp1, llm, nqtot) ! tracers
     110  REAL, INTENT(IN) :: pphis(iip1, jjp1) ! surface geopotential
     111  REAL, INTENT(IN) :: pphi(iip1, jjp1, llm) ! geopotential
     112
     113  REAL, INTENT(IN) :: pdvcov(iip1, jjm, llm) ! dynamical tendency on vcov
     114  REAL, INTENT(IN) :: pducov(iip1, jjp1, llm) ! dynamical tendency on ucov
     115  REAL, INTENT(IN) :: pdteta(iip1, jjp1, llm) ! dynamical tendency on teta
    116116  ! NB: pdteta is used only to compute pcvgt which is in fact not used...
    117   REAL,INTENT(IN) :: pdq(iip1,jjp1,llm,nqtot) ! dynamical tendency on tracers
     117  REAL, INTENT(IN) :: pdq(iip1, jjp1, llm, nqtot) ! dynamical tendency on tracers
    118118  ! NB: pdq is only used to compute pcvgq which is in fact not used...
    119119
    120   REAL,INTENT(IN) :: pps(iip1,jjp1) ! surface pressure (Pa)
    121   REAL,INTENT(IN) :: pp(iip1,jjp1,llmp1) ! pressure at mesh interfaces (Pa)
    122   REAL,INTENT(IN) :: ppk(iip1,jjp1,llm) ! Exner at mid-layer
    123   REAL,INTENT(IN) :: flxw(iip1,jjp1,llm) ! Vertical mass flux on lower mesh interfaces (kg/s) (on llm because flxw(:,:,llm+1)=0)
     120  REAL, INTENT(IN) :: pps(iip1, jjp1) ! surface pressure (Pa)
     121  REAL, INTENT(IN) :: pp(iip1, jjp1, llmp1) ! pressure at mesh interfaces (Pa)
     122  REAL, INTENT(IN) :: ppk(iip1, jjp1, llm) ! Exner at mid-layer
     123  REAL, INTENT(IN) :: flxw(iip1, jjp1, llm) ! Vertical mass flux on lower mesh interfaces (kg/s) (on llm because flxw(:,:,llm+1)=0)
    124124
    125125  ! tendencies (in */s) from the physics
    126   REAL,INTENT(OUT) :: pdvfi(iip1,jjm,llm) ! tendency on covariant meridional wind
    127   REAL,INTENT(OUT) :: pdufi(iip1,jjp1,llm) ! tendency on covariant zonal wind
    128   REAL,INTENT(OUT) :: pdhfi(iip1,jjp1,llm) ! tendency on potential temperature (K/s)
    129   REAL,INTENT(OUT) :: pdqfi(iip1,jjp1,llm,nqtot) ! tendency on tracers
    130   REAL,INTENT(OUT) :: pdpsfi(iip1,jjp1) ! tendency on surface pressure (Pa/s)
     126  REAL, INTENT(OUT) :: pdvfi(iip1, jjm, llm) ! tendency on covariant meridional wind
     127  REAL, INTENT(OUT) :: pdufi(iip1, jjp1, llm) ! tendency on covariant zonal wind
     128  REAL, INTENT(OUT) :: pdhfi(iip1, jjp1, llm) ! tendency on potential temperature (K/s)
     129  REAL, INTENT(OUT) :: pdqfi(iip1, jjp1, llm, nqtot) ! tendency on tracers
     130  REAL, INTENT(OUT) :: pdpsfi(iip1, jjp1) ! tendency on surface pressure (Pa/s)
    131131
    132132
     
    134134  !    -----------------
    135135
    136   INTEGER :: i,j,l,ig0,ig,iq,itr
     136  INTEGER :: i, j, l, ig0, ig, iq, itr
    137137  REAL :: zpsrf(ngridmx)
    138   REAL :: zplev(ngridmx,llm+1),zplay(ngridmx,llm)
    139   REAL :: zphi(ngridmx,llm),zphis(ngridmx)
    140   !
    141   REAL :: zrot(iip1,jjm,llm) ! AdlC May 2014
    142   REAL :: zufi(ngridmx,llm), zvfi(ngridmx,llm)
    143   REAL :: zrfi(ngridmx,llm) ! relative wind vorticity
    144   REAL :: ztfi(ngridmx,llm),zqfi(ngridmx,llm,nqtot)
    145   REAL :: zpk(ngridmx,llm)
    146   !
    147   REAL :: pcvgu(ngridmx,llm), pcvgv(ngridmx,llm)
    148   REAL :: pcvgt(ngridmx,llm), pcvgq(ngridmx,llm,2)
    149   !
    150   REAL :: zdufi(ngridmx,llm),zdvfi(ngridmx,llm)
    151   REAL :: zdtfi(ngridmx,llm),zdqfi(ngridmx,llm,nqtot)
     138  REAL :: zplev(ngridmx, llm + 1), zplay(ngridmx, llm)
     139  REAL :: zphi(ngridmx, llm), zphis(ngridmx)
     140  !
     141  REAL :: zrot(iip1, jjm, llm) ! AdlC May 2014
     142  REAL :: zufi(ngridmx, llm), zvfi(ngridmx, llm)
     143  REAL :: zrfi(ngridmx, llm) ! relative wind vorticity
     144  REAL :: ztfi(ngridmx, llm), zqfi(ngridmx, llm, nqtot)
     145  REAL :: zpk(ngridmx, llm)
     146  !
     147  REAL :: pcvgu(ngridmx, llm), pcvgv(ngridmx, llm)
     148  REAL :: pcvgt(ngridmx, llm), pcvgq(ngridmx, llm, 2)
     149  !
     150  REAL :: zdufi(ngridmx, llm), zdvfi(ngridmx, llm)
     151  REAL :: zdtfi(ngridmx, llm), zdqfi(ngridmx, llm, nqtot)
    152152  REAL :: zdpsrf(ngridmx)
    153153  !
    154   REAL :: zdufic(ngridmx,llm),zdvfic(ngridmx,llm)
    155   REAL :: zdtfic(ngridmx,llm),zdqfic(ngridmx,llm,nqtot)
    156   REAL :: jH_cur_split,zdt_split
    157   LOGICAL :: debut_split,lafin_split
     154  REAL :: zdufic(ngridmx, llm), zdvfic(ngridmx, llm)
     155  REAL :: zdtfic(ngridmx, llm), zdqfic(ngridmx, llm, nqtot)
     156  REAL :: jH_cur_split, zdt_split
     157  LOGICAL :: debut_split, lafin_split
    158158  INTEGER :: isplit
    159159
    160   REAL :: zsin(iim),zcos(iim),z1(iim)
    161   REAL :: zsinbis(iim),zcosbis(iim),z1bis(iim)
     160  REAL :: zsin(iim), zcos(iim), z1(iim)
     161  REAL :: zsinbis(iim), zcosbis(iim), z1bis(iim)
    162162  REAL :: unskap, pksurcp
    163163  !
    164   REAL :: flxwfi(ngridmx,llm)  ! Flux de masse verticale sur la grille physiq
     164  REAL :: flxwfi(ngridmx, llm)  ! Flux de masse verticale sur la grille physiq
    165165  !
    166166
    167167  REAL :: SSUM
    168168
    169   LOGICAL,SAVE :: firstcal=.TRUE., debut=.TRUE.
    170    ! REAL rdayvrai
     169  LOGICAL, SAVE :: firstcal = .TRUE., debut = .TRUE.
     170  ! REAL rdayvrai
    171171
    172172  !
     
    177177  !
    178178  !
    179   IF ( firstcal )  THEN
     179  IF (firstcal)  THEN
    180180    debut = .TRUE.
    181     IF (ngridmx/=2+(jjm-1)*iim) THEN
    182      WRITE(lunout,*) 'STOP dans calfis'
    183      WRITE(lunout,*) &
    184            'La dimension ngridmx doit etre egale a 2 + (jjm-1)*iim'
    185      WRITE(lunout,*) '  ngridmx  jjm   iim   '
    186      WRITE(lunout,*) ngridmx,jjm,iim
    187      CALL abort_gcm("calfis", "", 1)
     181    IF (ngridmx/=2 + (jjm - 1) * iim) THEN
     182      WRITE(lunout, *) 'STOP dans calfis'
     183      WRITE(lunout, *) &
     184              'La dimension ngridmx doit etre egale a 2 + (jjm-1)*iim'
     185      WRITE(lunout, *) '  ngridmx  jjm   iim   '
     186      WRITE(lunout, *) ngridmx, jjm, iim
     187      CALL abort_gcm("calfis", "", 1)
    188188    ENDIF
    189189  ELSE
     
    200200  !   ----------------------------------
    201201
    202 
    203   zpsrf(1) = pps(1,1)
    204 
    205   ig0  = 2
    206   DO j = 2,jjm
    207      CALL SCOPY( iim,pps(1,j),1,zpsrf(ig0), 1 )
    208      ig0 = ig0+iim
    209   ENDDO
    210 
    211   zpsrf(ngridmx) = pps(1,jjp1)
     202  zpsrf(1) = pps(1, 1)
     203
     204  ig0 = 2
     205  DO j = 2, jjm
     206    CALL SCOPY(iim, pps(1, j), 1, zpsrf(ig0), 1)
     207    ig0 = ig0 + iim
     208  ENDDO
     209
     210  zpsrf(ngridmx) = pps(1, jjp1)
    212211
    213212
     
    221220  !    ...    Exner = cp * ( p(l) / preff ) ** kappa     ....
    222221  !
    223    unskap   = 1./ kappa
    224   !
    225   DO l = 1, llm
    226     zpk(   1,l ) = ppk(1,1,l)
    227     zplev( 1,l ) = pp(1,1,l)
     222  unskap = 1. / kappa
     223  !
     224  DO l = 1, llm
     225    zpk(1, l) = ppk(1, 1, l)
     226    zplev(1, l) = pp(1, 1, l)
    228227    ig0 = 2
    229       DO j = 2, jjm
    230          DO i =1, iim
    231           zpk(   ig0,l ) = ppk(i,j,l)
    232           zplev( ig0,l ) = pp(i,j,l)
    233           ig0 = ig0 +1
    234          ENDDO
    235       ENDDO
    236     zpk(   ngridmx,l ) = ppk(1,jjp1,l)
    237     zplev( ngridmx,l ) = pp(1,jjp1,l)
    238   ENDDO
    239     zplev( 1,llmp1 ) = pp(1,1,llmp1)
    240     ig0 = 2
    241       DO j = 2, jjm
    242          DO i =1, iim
    243           zplev( ig0,llmp1 ) = pp(i,j,llmp1)
    244           ig0 = ig0 +1
    245          ENDDO
    246       ENDDO
    247     zplev( ngridmx,llmp1 ) = pp(1,jjp1,llmp1)
     228    DO j = 2, jjm
     229      DO i = 1, iim
     230        zpk(ig0, l) = ppk(i, j, l)
     231        zplev(ig0, l) = pp(i, j, l)
     232        ig0 = ig0 + 1
     233      ENDDO
     234    ENDDO
     235    zpk(ngridmx, l) = ppk(1, jjp1, l)
     236    zplev(ngridmx, l) = pp(1, jjp1, l)
     237  ENDDO
     238  zplev(1, llmp1) = pp(1, 1, llmp1)
     239  ig0 = 2
     240  DO j = 2, jjm
     241    DO i = 1, iim
     242      zplev(ig0, llmp1) = pp(i, j, llmp1)
     243      ig0 = ig0 + 1
     244    ENDDO
     245  ENDDO
     246  zplev(ngridmx, llmp1) = pp(1, jjp1, llmp1)
    248247  !
    249248  !
     
    252251  !   ---------------------------------------------------------------
    253252
    254   DO l=1,llm
    255 
    256      pksurcp     =  ppk(1,1,l) / cpp
    257      zplay(1,l)  = preff * pksurcp ** unskap
    258      ztfi(1,l)   =  pteta(1,1,l) * pksurcp
    259      pcvgt(1,l)  =  pdteta(1,1,l) * pksurcp / pmasse(1,1,l)
    260      ig0        = 2
    261 
    262      DO j = 2, jjm
    263         DO i = 1, iim
    264           pksurcp        = ppk(i,j,l) / cpp
    265           zplay(ig0,l)  = preff * pksurcp ** unskap
    266           ztfi(ig0,l)    = pteta(i,j,l) * pksurcp
    267           pcvgt(ig0,l)   = pdteta(i,j,l) * pksurcp / pmasse(i,j,l)
    268           ig0            = ig0 + 1
    269         ENDDO
    270      ENDDO
    271 
    272      pksurcp       = ppk(1,jjp1,l) / cpp
    273      zplay(ig0,l) = preff * pksurcp ** unskap
    274      ztfi (ig0,l)  = pteta(1,jjp1,l) * pksurcp
    275      pcvgt(ig0,l)  = pdteta(1,jjp1,l) * pksurcp/ pmasse(1,jjp1,l)
     253  DO l = 1, llm
     254
     255    pksurcp = ppk(1, 1, l) / cpp
     256    zplay(1, l) = preff * pksurcp ** unskap
     257    ztfi(1, l) = pteta(1, 1, l) * pksurcp
     258    pcvgt(1, l) = pdteta(1, 1, l) * pksurcp / pmasse(1, 1, l)
     259    ig0 = 2
     260
     261    DO j = 2, jjm
     262      DO i = 1, iim
     263        pksurcp = ppk(i, j, l) / cpp
     264        zplay(ig0, l) = preff * pksurcp ** unskap
     265        ztfi(ig0, l) = pteta(i, j, l) * pksurcp
     266        pcvgt(ig0, l) = pdteta(i, j, l) * pksurcp / pmasse(i, j, l)
     267        ig0 = ig0 + 1
     268      ENDDO
     269    ENDDO
     270
     271    pksurcp = ppk(1, jjp1, l) / cpp
     272    zplay(ig0, l) = preff * pksurcp ** unskap
     273    ztfi (ig0, l) = pteta(1, jjp1, l) * pksurcp
     274    pcvgt(ig0, l) = pdteta(1, jjp1, l) * pksurcp / pmasse(1, jjp1, l)
    276275
    277276  ENDDO
     
    280279  !   ---------------
    281280  !
    282   itr=0
    283   DO iq=1,nqtot
    284      IF(.NOT.tracers(iq)%isAdvected) CYCLE
    285      itr = itr + 1
    286      DO l=1,llm
    287         zqfi(1,l,itr) = pq(1,1,l,iq)
    288         ig0           = 2
    289         DO j=2,jjm
    290            DO i = 1, iim
    291               zqfi(ig0,l,itr) = pq(i,j,l,iq)
    292               ig0             = ig0 + 1
    293            ENDDO
     281  itr = 0
     282  DO iq = 1, nqtot
     283    IF(.NOT.tracers(iq)%isAdvected) CYCLE
     284    itr = itr + 1
     285    DO l = 1, llm
     286      zqfi(1, l, itr) = pq(1, 1, l, iq)
     287      ig0 = 2
     288      DO j = 2, jjm
     289        DO i = 1, iim
     290          zqfi(ig0, l, itr) = pq(i, j, l, iq)
     291          ig0 = ig0 + 1
    294292        ENDDO
    295         zqfi(ig0,l,itr) = pq(1,jjp1,l,iq)
    296      ENDDO
     293      ENDDO
     294      zqfi(ig0, l, itr) = pq(1, jjp1, l, iq)
     295    ENDDO
    297296  ENDDO
    298297
    299298  !   convergence dynamique pour les traceurs "EAU"
    300299  ! Earth-specific treatment of first 2 tracers (water)
    301    IF (planet_type=="earth") THEN
    302     DO iq=1,2
    303      DO l=1,llm
    304         pcvgq(1,l,iq)= pdq(1,1,l,iq) / pmasse(1,1,l)
    305         ig0          = 2
    306         DO j=2,jjm
    307            DO i = 1, iim
    308               pcvgq(ig0,l,iq) = pdq(i,j,l,iq) / pmasse(i,j,l)
    309               ig0            = ig0 + 1
    310            ENDDO
     300  IF (planet_type=="earth") THEN
     301    DO iq = 1, 2
     302      DO l = 1, llm
     303        pcvgq(1, l, iq) = pdq(1, 1, l, iq) / pmasse(1, 1, l)
     304        ig0 = 2
     305        DO j = 2, jjm
     306          DO i = 1, iim
     307            pcvgq(ig0, l, iq) = pdq(i, j, l, iq) / pmasse(i, j, l)
     308            ig0 = ig0 + 1
     309          ENDDO
    311310        ENDDO
    312         pcvgq(ig0,l,iq)= pdq(1,jjp1,l,iq) / pmasse(1,jjp1,l)
    313      ENDDO
    314     ENDDO
    315    endif ! of if (planet_type=="earth")
     311        pcvgq(ig0, l, iq) = pdq(1, jjp1, l, iq) / pmasse(1, jjp1, l)
     312      ENDDO
     313    ENDDO
     314  endif ! of if (planet_type=="earth")
    316315
    317316
     
    319318  !   -----------------------------------------------------
    320319
    321   CALL gr_dyn_fi(llm,iip1,jjp1,ngridmx,pphi,zphi)
    322   CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,pphis,zphis)
    323   DO l=1,llm
    324      DO ig=1,ngridmx
    325        zphi(ig,l)=zphi(ig,l)-zphis(ig)
    326      ENDDO
     320  CALL gr_dyn_fi(llm, iip1, jjp1, ngridmx, pphi, zphi)
     321  CALL gr_dyn_fi(1, iip1, jjp1, ngridmx, pphis, zphis)
     322  DO l = 1, llm
     323    DO ig = 1, ngridmx
     324      zphi(ig, l) = zphi(ig, l) - zphis(ig)
     325    ENDDO
    327326  ENDDO
    328327
     
    330329  ! JG : ancien calcule de omega utilise dans physiq.F. Maintenant le flux
    331330  !    de masse est calclue dans advtrac.F
    332    ! DO l=1,llm
    333    !   pvervel(1,l)=pw(1,1,l) * g /apoln
    334    !   ig0=2
    335    !  DO j=2,jjm
    336    !      DO i = 1, iim
    337    !         pvervel(ig0,l) = pw(i,j,l) * g * unsaire(i,j)
    338    !         ig0 = ig0 + 1
    339    !      ENDDO
     331  ! DO l=1,llm
     332  !   pvervel(1,l)=pw(1,1,l) * g /apoln
     333  !   ig0=2
     334  !  DO j=2,jjm
     335  !      DO i = 1, iim
     336  !         pvervel(ig0,l) = pw(i,j,l) * g * unsaire(i,j)
     337  !         ig0 = ig0 + 1
     338  !      ENDDO
    340339  !       ENDDO
    341    !   pvervel(ig0,l)=pw(1,jjp1,l) * g /apols
    342    ! ENDDO
     340  !   pvervel(ig0,l)=pw(1,jjp1,l) * g /apols
     341  ! ENDDO
    343342
    344343  !
     
    346345  !   ------------
    347346
    348   DO l=1,llm
    349 
    350      DO j=2,jjm
    351         ig0 = 1+(j-2)*iim
    352         zufi(ig0+1,l)= 0.5 * &
    353               ( pucov(iim,j,l)/cu(iim,j) + pucov(1,j,l)/cu(1,j) )
    354         pcvgu(ig0+1,l)= 0.5 * &
    355               ( pducov(iim,j,l)/cu(iim,j) + pducov(1,j,l)/cu(1,j) )
    356         DO i=2,iim
    357            zufi(ig0+i,l)= 0.5 * &
    358                  ( pucov(i-1,j,l)/cu(i-1,j) + pucov(i,j,l)/cu(i,j) )
    359            pcvgu(ig0+i,l)= 0.5 * &
    360                  ( pducov(i-1,j,l)/cu(i-1,j) + pducov(i,j,l)/cu(i,j) )
    361   END DO
    362   END DO
     347  DO l = 1, llm
     348
     349    DO j = 2, jjm
     350      ig0 = 1 + (j - 2) * iim
     351      zufi(ig0 + 1, l) = 0.5 * &
     352              (pucov(iim, j, l) / cu(iim, j) + pucov(1, j, l) / cu(1, j))
     353      pcvgu(ig0 + 1, l) = 0.5 * &
     354              (pducov(iim, j, l) / cu(iim, j) + pducov(1, j, l) / cu(1, j))
     355      DO i = 2, iim
     356        zufi(ig0 + i, l) = 0.5 * &
     357                (pucov(i - 1, j, l) / cu(i - 1, j) + pucov(i, j, l) / cu(i, j))
     358        pcvgu(ig0 + i, l) = 0.5 * &
     359                (pducov(i - 1, j, l) / cu(i - 1, j) + pducov(i, j, l) / cu(i, j))
     360      END DO
     361    END DO
    363362
    364363  END DO
     
    368367  !  46.1 Calcul de la vorticite et passage sur la grille physique
    369368  !  --------------------------------------------------------------
    370   DO l=1,llm
    371     do i=1,iim
    372       do j=1,jjm
    373         zrot(i,j,l) = (pvcov(i+1,j,l) - pvcov(i,j,l) &
    374               + pucov(i,j+1,l) - pucov(i,j,l)) &
    375               / (cu(i,j)+cu(i,j+1)) &
    376               / (cv(i+1,j)+cv(i,j)) *4
     369  DO l = 1, llm
     370    do i = 1, iim
     371      do j = 1, jjm
     372        zrot(i, j, l) = (pvcov(i + 1, j, l) - pvcov(i, j, l) &
     373                + pucov(i, j + 1, l) - pucov(i, j, l)) &
     374                / (cu(i, j) + cu(i, j + 1)) &
     375                / (cv(i + 1, j) + cv(i, j)) * 4
    377376      enddo
    378377    enddo
     
    382381  !   -----------
    383382
    384   DO l=1,llm
    385      DO j=2,jjm
    386         ig0=1+(j-2)*iim
    387         DO i=1,iim
    388            zvfi(ig0+i,l)= 0.5 * &
    389                  ( pvcov(i,j-1,l)/cv(i,j-1) + pvcov(i,j,l)/cv(i,j) )
    390            pcvgv(ig0+i,l)= 0.5 * &
    391                  ( pdvcov(i,j-1,l)/cv(i,j-1) + pdvcov(i,j,l)/cv(i,j) )
    392         ENDDO
    393            zrfi(ig0 + 1,l)= 0.25 *(zrot(iim,j-1,l)+zrot(iim,j,l) &
    394                  +zrot(1,j-1,l)+zrot(1,j,l))
    395         DO i=2,iim
    396            zrfi(ig0 + i,l)= 0.25 *(zrot(i-1,j-1,l)+zrot(i-1,j,l) &
    397                  +zrot(i,j-1,l)+zrot(i,j,l))   !  AdlC MAY 2014
    398         ENDDO
    399      ENDDO
     383  DO l = 1, llm
     384    DO j = 2, jjm
     385      ig0 = 1 + (j - 2) * iim
     386      DO i = 1, iim
     387        zvfi(ig0 + i, l) = 0.5 * &
     388                (pvcov(i, j - 1, l) / cv(i, j - 1) + pvcov(i, j, l) / cv(i, j))
     389        pcvgv(ig0 + i, l) = 0.5 * &
     390                (pdvcov(i, j - 1, l) / cv(i, j - 1) + pdvcov(i, j, l) / cv(i, j))
     391      ENDDO
     392      zrfi(ig0 + 1, l) = 0.25 * (zrot(iim, j - 1, l) + zrot(iim, j, l) &
     393              + zrot(1, j - 1, l) + zrot(1, j, l))
     394      DO i = 2, iim
     395        zrfi(ig0 + i, l) = 0.25 * (zrot(i - 1, j - 1, l) + zrot(i - 1, j, l) &
     396                + zrot(i, j - 1, l) + zrot(i, j, l))   !  AdlC MAY 2014
     397      ENDDO
     398    ENDDO
    400399  ENDDO
    401400
     
    403402  !   47. champs de vents aux pole nord
    404403  !   ------------------------------
    405      ! U = 1 / pi  *  integrale [ v * cos(long) * d long ]
    406      ! V = 1 / pi  *  integrale [ v * sin(long) * d long ]
    407 
    408   DO l=1,llm
    409 
    410      z1(1)   =(rlonu(1)-rlonu(iim)+2.*pi)*pvcov(1,1,l)/cv(1,1)
    411      z1bis(1)=(rlonu(1)-rlonu(iim)+2.*pi)*pdvcov(1,1,l)/cv(1,1)
    412      DO i=2,iim
    413         z1(i)   =(rlonu(i)-rlonu(i-1))*pvcov(i,1,l)/cv(i,1)
    414         z1bis(i)=(rlonu(i)-rlonu(i-1))*pdvcov(i,1,l)/cv(i,1)
    415      ENDDO
    416 
    417      DO i=1,iim
    418         zcos(i)   = COS(rlonv(i))*z1(i)
    419         zcosbis(i)= COS(rlonv(i))*z1bis(i)
    420         zsin(i)   = SIN(rlonv(i))*z1(i)
    421         zsinbis(i)= SIN(rlonv(i))*z1bis(i)
    422      ENDDO
    423 
    424      zufi(1,l)  = SSUM(iim,zcos,1)/pi
    425      pcvgu(1,l) = SSUM(iim,zcosbis,1)/pi
    426      zvfi(1,l)  = SSUM(iim,zsin,1)/pi
    427      pcvgv(1,l) = SSUM(iim,zsinbis,1)/pi
    428      zrfi(1, l) = 0.
     404  ! U = 1 / pi  *  integrale [ v * cos(long) * d long ]
     405  ! V = 1 / pi  *  integrale [ v * sin(long) * d long ]
     406
     407  DO l = 1, llm
     408
     409    z1(1) = (rlonu(1) - rlonu(iim) + 2. * pi) * pvcov(1, 1, l) / cv(1, 1)
     410    z1bis(1) = (rlonu(1) - rlonu(iim) + 2. * pi) * pdvcov(1, 1, l) / cv(1, 1)
     411    DO i = 2, iim
     412      z1(i) = (rlonu(i) - rlonu(i - 1)) * pvcov(i, 1, l) / cv(i, 1)
     413      z1bis(i) = (rlonu(i) - rlonu(i - 1)) * pdvcov(i, 1, l) / cv(i, 1)
     414    ENDDO
     415
     416    DO i = 1, iim
     417      zcos(i) = COS(rlonv(i)) * z1(i)
     418      zcosbis(i) = COS(rlonv(i)) * z1bis(i)
     419      zsin(i) = SIN(rlonv(i)) * z1(i)
     420      zsinbis(i) = SIN(rlonv(i)) * z1bis(i)
     421    ENDDO
     422
     423    zufi(1, l) = SSUM(iim, zcos, 1) / pi
     424    pcvgu(1, l) = SSUM(iim, zcosbis, 1) / pi
     425    zvfi(1, l) = SSUM(iim, zsin, 1) / pi
     426    pcvgv(1, l) = SSUM(iim, zsinbis, 1) / pi
     427    zrfi(1, l) = 0.
    429428  ENDDO
    430429
     
    432431  !   48. champs de vents aux pole sud:
    433432  !   ---------------------------------
    434      ! U = 1 / pi  *  integrale [ v * cos(long) * d long ]
    435      ! V = 1 / pi  *  integrale [ v * sin(long) * d long ]
    436 
    437   DO l=1,llm
    438 
    439      z1(1)   =(rlonu(1)-rlonu(iim)+2.*pi)*pvcov(1,jjm,l)/cv(1,jjm)
    440      z1bis(1)=(rlonu(1)-rlonu(iim)+2.*pi)*pdvcov(1,jjm,l)/cv(1,jjm)
    441      DO i=2,iim
    442         z1(i)   =(rlonu(i)-rlonu(i-1))*pvcov(i,jjm,l)/cv(i,jjm)
    443         z1bis(i)=(rlonu(i)-rlonu(i-1))*pdvcov(i,jjm,l)/cv(i,jjm)
    444      ENDDO
    445 
    446      DO i=1,iim
    447         zcos(i)    = COS(rlonv(i))*z1(i)
    448         zcosbis(i) = COS(rlonv(i))*z1bis(i)
    449         zsin(i)    = SIN(rlonv(i))*z1(i)
    450         zsinbis(i) = SIN(rlonv(i))*z1bis(i)
    451      ENDDO
    452 
    453      zufi(ngridmx,l)  = SSUM(iim,zcos,1)/pi
    454      pcvgu(ngridmx,l) = SSUM(iim,zcosbis,1)/pi
    455      zvfi(ngridmx,l)  = SSUM(iim,zsin,1)/pi
    456      pcvgv(ngridmx,l) = SSUM(iim,zsinbis,1)/pi
    457      zrfi(ngridmx, l) = 0.
     433  ! U = 1 / pi  *  integrale [ v * cos(long) * d long ]
     434  ! V = 1 / pi  *  integrale [ v * sin(long) * d long ]
     435
     436  DO l = 1, llm
     437
     438    z1(1) = (rlonu(1) - rlonu(iim) + 2. * pi) * pvcov(1, jjm, l) / cv(1, jjm)
     439    z1bis(1) = (rlonu(1) - rlonu(iim) + 2. * pi) * pdvcov(1, jjm, l) / cv(1, jjm)
     440    DO i = 2, iim
     441      z1(i) = (rlonu(i) - rlonu(i - 1)) * pvcov(i, jjm, l) / cv(i, jjm)
     442      z1bis(i) = (rlonu(i) - rlonu(i - 1)) * pdvcov(i, jjm, l) / cv(i, jjm)
     443    ENDDO
     444
     445    DO i = 1, iim
     446      zcos(i) = COS(rlonv(i)) * z1(i)
     447      zcosbis(i) = COS(rlonv(i)) * z1bis(i)
     448      zsin(i) = SIN(rlonv(i)) * z1(i)
     449      zsinbis(i) = SIN(rlonv(i)) * z1bis(i)
     450    ENDDO
     451
     452    zufi(ngridmx, l) = SSUM(iim, zcos, 1) / pi
     453    pcvgu(ngridmx, l) = SSUM(iim, zcosbis, 1) / pi
     454    zvfi(ngridmx, l) = SSUM(iim, zsin, 1) / pi
     455    pcvgv(ngridmx, l) = SSUM(iim, zsinbis, 1) / pi
     456    zrfi(ngridmx, l) = 0.
    458457  ENDDO
    459458  !
    460459  ! On change de grille, dynamique vers physiq, pour le flux de masse verticale
    461   CALL gr_dyn_fi(llm,iip1,jjp1,ngridmx,flxw,flxwfi)
     460  CALL gr_dyn_fi(llm, iip1, jjp1, ngridmx, flxw, flxwfi)
    462461
    463462  !-----------------------------------------------------------------------
     
    467466
    468467
    469    ! WRITE(lunout,*) 'PHYSIQUE AVEC NSPLIT_PHYS=',nsplit_phys
    470   zdt_split=dtphys/nsplit_phys
    471   zdufic(:,:)=0.
    472   zdvfic(:,:)=0.
    473   zdtfic(:,:)=0.
    474   zdqfic(:,:,:)=0.
    475 
    476    IF (CPPKEY_PHYS) THEN
    477 
    478    do isplit=1,nsplit_phys
    479 
    480      jH_cur_split=jH_cur+(isplit-1) * dtvr / (daysec *nsplit_phys)
    481      debut_split=debut.AND.isplit==1
    482      lafin_split=lafin.AND.isplit==nsplit_phys
    483 
    484    ! if (planet_type=="earth") THEN
    485     CALL call_physiq(ngridmx,llm,nqtot,tracers(:)%name, &
    486           debut_split,lafin_split, &
    487           jD_cur,jH_cur_split,zdt_split, &
    488           zplev,zplay, &
    489           zpk,zphi,zphis, &
    490           presnivs, &
    491           zufi,zvfi,zrfi,ztfi,zqfi, &
    492           flxwfi,pducov, &
    493           zdufi,zdvfi,zdtfi,zdqfi,zdpsrf)
    494 
    495    ! ELSE IF ( planet_type=="generic" ) THEN
    496    !    CALL physiq (ngridmx,     !! ngrid
    497   ! .             llm,            !! nlayer
    498   ! .             nqtot,          !! nq
    499   ! .             tracers(:)%name,!! tracer names from dynamical core (given in infotrac)
    500   ! .             debut_split,    !! firstcall
    501   ! .             lafin_split,    !! lastcall
    502   ! .             jD_cur,         !! pday. see leapfrog
    503   ! .             jH_cur_split,   !! ptime "fraction of day"
    504   ! .             zdt_split,      !! ptimestep
    505   ! .             zplev,          !! pplev
    506   ! .             zplay,          !! pplay
    507   ! .             zphi,           !! pphi
    508   ! .             zufi,           !! pu
    509   ! .             zvfi,           !! pv
    510   ! .             ztfi,           !! pt
    511   ! .             zqfi,           !! pq
    512   ! .             flxwfi,         !! pw !! or 0. anyway this is for diagnostic. not used in physiq.
    513   ! .             zdufi,          !! pdu
    514   ! .             zdvfi,          !! pdv
    515   ! .             zdtfi,          !! pdt
    516   ! .             zdqfi,          !! pdq
    517   ! .             zdpsrf,         !! pdpsrf
    518   ! .             tracerdyn)      !! tracerdyn <-- utilite ???
    519 
    520   !  ENDIF ! of if (planet_type=="earth")
    521 
    522      zufi(:,:)=zufi(:,:)+zdufi(:,:)*zdt_split
    523      zvfi(:,:)=zvfi(:,:)+zdvfi(:,:)*zdt_split
    524      ztfi(:,:)=ztfi(:,:)+zdtfi(:,:)*zdt_split
    525      zqfi(:,:,:)=zqfi(:,:,:)+zdqfi(:,:,:)*zdt_split
    526 
    527      zdufic(:,:)=zdufic(:,:)+zdufi(:,:)
    528      zdvfic(:,:)=zdvfic(:,:)+zdvfi(:,:)
    529      zdtfic(:,:)=zdtfic(:,:)+zdtfi(:,:)
    530      zdqfic(:,:,:)=zdqfic(:,:,:)+zdqfi(:,:,:)
    531 
    532    enddo ! of do isplit=1,nsplit_phys
    533 
    534    END IF
    535 
    536   zdufi(:,:)=zdufic(:,:)/nsplit_phys
    537   zdvfi(:,:)=zdvfic(:,:)/nsplit_phys
    538   zdtfi(:,:)=zdtfic(:,:)/nsplit_phys
    539   zdqfi(:,:,:)=zdqfic(:,:,:)/nsplit_phys
     468  ! WRITE(lunout,*) 'PHYSIQUE AVEC NSPLIT_PHYS=',nsplit_phys
     469  zdt_split = dtphys / nsplit_phys
     470  zdufic(:, :) = 0.
     471  zdvfic(:, :) = 0.
     472  zdtfic(:, :) = 0.
     473  zdqfic(:, :, :) = 0.
     474
     475  IF (CPPKEY_PHYS) THEN
     476
     477    do isplit = 1, nsplit_phys
     478
     479      jH_cur_split = jH_cur + (isplit - 1) * dtvr / (daysec * nsplit_phys)
     480      debut_split = debut.AND.isplit==1
     481      lafin_split = lafin.AND.isplit==nsplit_phys
     482
     483      ! if (planet_type=="earth") THEN
     484      CALL call_physiq(ngridmx, llm, nqtot, tracers(:)%name, &
     485              debut_split, lafin_split, &
     486              jD_cur, jH_cur_split, zdt_split, &
     487              zplev, zplay, &
     488              zpk, zphi, zphis, &
     489              presnivs, &
     490              zufi, zvfi, zrfi, ztfi, zqfi, &
     491              flxwfi, pducov, &
     492              zdufi, zdvfi, zdtfi, zdqfi, zdpsrf)
     493
     494      ! ELSE IF ( planet_type=="generic" ) THEN
     495      !    CALL physiq (ngridmx,     !! ngrid
     496      ! .             llm,            !! nlayer
     497      ! .             nqtot,          !! nq
     498      ! .             tracers(:)%name,!! tracer names from dynamical core (given in infotrac)
     499      ! .             debut_split,    !! firstcall
     500      ! .             lafin_split,    !! lastcall
     501      ! .             jD_cur,         !! pday. see leapfrog
     502      ! .             jH_cur_split,   !! ptime "fraction of day"
     503      ! .             zdt_split,      !! ptimestep
     504      ! .             zplev,          !! pplev
     505      ! .             zplay,          !! pplay
     506      ! .             zphi,           !! pphi
     507      ! .             zufi,           !! pu
     508      ! .             zvfi,           !! pv
     509      ! .             ztfi,           !! pt
     510      ! .             zqfi,           !! pq
     511      ! .             flxwfi,         !! pw !! or 0. anyway this is for diagnostic. not used in physiq.
     512      ! .             zdufi,          !! pdu
     513      ! .             zdvfi,          !! pdv
     514      ! .             zdtfi,          !! pdt
     515      ! .             zdqfi,          !! pdq
     516      ! .             zdpsrf,         !! pdpsrf
     517      ! .             tracerdyn)      !! tracerdyn <-- utilite ???
     518
     519      !  ENDIF ! of if (planet_type=="earth")
     520
     521      zufi(:, :) = zufi(:, :) + zdufi(:, :) * zdt_split
     522      zvfi(:, :) = zvfi(:, :) + zdvfi(:, :) * zdt_split
     523      ztfi(:, :) = ztfi(:, :) + zdtfi(:, :) * zdt_split
     524      zqfi(:, :, :) = zqfi(:, :, :) + zdqfi(:, :, :) * zdt_split
     525
     526      zdufic(:, :) = zdufic(:, :) + zdufi(:, :)
     527      zdvfic(:, :) = zdvfic(:, :) + zdvfi(:, :)
     528      zdtfic(:, :) = zdtfic(:, :) + zdtfi(:, :)
     529      zdqfic(:, :, :) = zdqfic(:, :, :) + zdqfi(:, :, :)
     530
     531    enddo ! of do isplit=1,nsplit_phys
     532
     533  END IF
     534
     535  zdufi(:, :) = zdufic(:, :) / nsplit_phys
     536  zdvfi(:, :) = zdvfic(:, :) / nsplit_phys
     537  zdtfi(:, :) = zdtfic(:, :) / nsplit_phys
     538  zdqfi(:, :, :) = zdqfic(:, :, :) / nsplit_phys
    540539
    541540  !-----------------------------------------------------------------------
     
    546545  !  -----------------------------------
    547546
    548   CALL gr_fi_dyn(1,ngridmx,iip1,jjp1,zdpsrf,pdpsfi)
     547  CALL gr_fi_dyn(1, ngridmx, iip1, jjp1, zdpsrf, pdpsfi)
    549548  !
    550549  !   62. enthalpie potentielle
    551550  !   ---------------------
    552551
    553   DO l=1,llm
    554 
    555      DO i=1,iip1
    556       pdhfi(i,1,l)    = cpp *  zdtfi(1,l)      / ppk(i, 1  ,l)
    557       pdhfi(i,jjp1,l) = cpp *  zdtfi(ngridmx,l)/ ppk(i,jjp1,l)
    558      ENDDO
    559 
    560      DO j=2,jjm
    561         ig0=1+(j-2)*iim
    562         DO i=1,iim
    563            pdhfi(i,j,l) = cpp * zdtfi(ig0+i,l) / ppk(i,j,l)
    564         ENDDO
    565            pdhfi(iip1,j,l) =  pdhfi(1,j,l)
    566      ENDDO
     552  DO l = 1, llm
     553
     554    DO i = 1, iip1
     555      pdhfi(i, 1, l) = cpp * zdtfi(1, l) / ppk(i, 1, l)
     556      pdhfi(i, jjp1, l) = cpp * zdtfi(ngridmx, l) / ppk(i, jjp1, l)
     557    ENDDO
     558
     559    DO j = 2, jjm
     560      ig0 = 1 + (j - 2) * iim
     561      DO i = 1, iim
     562        pdhfi(i, j, l) = cpp * zdtfi(ig0 + i, l) / ppk(i, j, l)
     563      ENDDO
     564      pdhfi(iip1, j, l) = pdhfi(1, j, l)
     565    ENDDO
    567566
    568567  ENDDO
     
    572571  !   ---------------------
    573572  ! Ehouarn: removed this useless bit: was overwritten at step 63 anyways
    574    ! DO iq=1,nqtot
    575    !    DO l=1,llm
    576    !       DO i=1,iip1
    577    !          pdqfi(i,1,l,iq)    = zdqfi(1,l,iq)
    578    !          pdqfi(i,jjp1,l,iq) = zdqfi(ngridmx,l,iq)
    579    !       ENDDO
    580    !       DO j=2,jjm
    581    !          ig0=1+(j-2)*iim
    582    !          DO i=1,iim
    583    !             pdqfi(i,j,l,iq) = zdqfi(ig0+i,l,iq)
    584    !          ENDDO
    585    !          pdqfi(iip1,j,l,iq) = pdqfi(1,j,l,iq)
    586    !       ENDDO
    587    !    ENDDO
    588    ! ENDDO
     573  ! DO iq=1,nqtot
     574  !    DO l=1,llm
     575  !       DO i=1,iip1
     576  !          pdqfi(i,1,l,iq)    = zdqfi(1,l,iq)
     577  !          pdqfi(i,jjp1,l,iq) = zdqfi(ngridmx,l,iq)
     578  !       ENDDO
     579  !       DO j=2,jjm
     580  !          ig0=1+(j-2)*iim
     581  !          DO i=1,iim
     582  !             pdqfi(i,j,l,iq) = zdqfi(ig0+i,l,iq)
     583  !          ENDDO
     584  !          pdqfi(iip1,j,l,iq) = pdqfi(1,j,l,iq)
     585  !       ENDDO
     586  !    ENDDO
     587  ! ENDDO
    589588
    590589  !   63. traceurs
    591590  !   ------------
    592591  ! initialisation des tendances
    593   pdqfi(:,:,:,:)=0.
     592  pdqfi(:, :, :, :) = 0.
    594593  !
    595594  itr = 0
    596   DO iq=1,nqtot
    597      IF(.NOT.tracers(iq)%isAdvected) CYCLE
    598      itr = itr + 1
    599      DO l=1,llm
    600         DO i=1,iip1
    601            pdqfi(i,1,l,iq)    = zdqfi(1,l,itr)
    602            pdqfi(i,jjp1,l,iq) = zdqfi(ngridmx,l,itr)
     595  DO iq = 1, nqtot
     596    IF(.NOT.tracers(iq)%isAdvected) CYCLE
     597    itr = itr + 1
     598    DO l = 1, llm
     599      DO i = 1, iip1
     600        pdqfi(i, 1, l, iq) = zdqfi(1, l, itr)
     601        pdqfi(i, jjp1, l, iq) = zdqfi(ngridmx, l, itr)
     602      ENDDO
     603      DO j = 2, jjm
     604        ig0 = 1 + (j - 2) * iim
     605        DO i = 1, iim
     606          pdqfi(i, j, l, iq) = zdqfi(ig0 + i, l, itr)
    603607        ENDDO
    604         DO j=2,jjm
    605            ig0=1+(j-2)*iim
    606            DO i=1,iim
    607               pdqfi(i,j,l,iq) = zdqfi(ig0+i,l,itr)
    608            ENDDO
    609            pdqfi(iip1,j,l,iq) = pdqfi(1,j,l,itr)
    610         ENDDO
    611      ENDDO
     608        pdqfi(iip1, j, l, iq) = pdqfi(1, j, l, itr)
     609      ENDDO
     610    ENDDO
    612611  ENDDO
    613612
     
    615614  !   ------------
    616615
    617   DO l=1,llm
    618 
    619      DO i=1,iip1
    620         pdufi(i,1,l)    = 0.
    621         pdufi(i,jjp1,l) = 0.
    622      ENDDO
    623 
    624      DO j=2,jjm
    625         ig0=1+(j-2)*iim
    626         DO i=1,iim-1
    627            pdufi(i,j,l)= &
    628                  0.5*(zdufi(ig0+i,l)+zdufi(ig0+i+1,l))*cu(i,j)
    629         ENDDO
    630         pdufi(iim,j,l)= &
    631               0.5*(zdufi(ig0+1,l)+zdufi(ig0+iim,l))*cu(iim,j)
    632         pdufi(iip1,j,l)=pdufi(1,j,l)
    633      ENDDO
     616  DO l = 1, llm
     617
     618    DO i = 1, iip1
     619      pdufi(i, 1, l) = 0.
     620      pdufi(i, jjp1, l) = 0.
     621    ENDDO
     622
     623    DO j = 2, jjm
     624      ig0 = 1 + (j - 2) * iim
     625      DO i = 1, iim - 1
     626        pdufi(i, j, l) = &
     627                0.5 * (zdufi(ig0 + i, l) + zdufi(ig0 + i + 1, l)) * cu(i, j)
     628      ENDDO
     629      pdufi(iim, j, l) = &
     630              0.5 * (zdufi(ig0 + 1, l) + zdufi(ig0 + iim, l)) * cu(iim, j)
     631      pdufi(iip1, j, l) = pdufi(1, j, l)
     632    ENDDO
    634633
    635634  ENDDO
     
    639638  !   ------------
    640639
    641   DO l=1,llm
    642 
    643      DO j=2,jjm-1
    644         ig0=1+(j-2)*iim
    645         DO i=1,iim
    646            pdvfi(i,j,l)= &
    647                  0.5*(zdvfi(ig0+i,l)+zdvfi(ig0+i+iim,l))*cv(i,j)
    648         ENDDO
    649         pdvfi(iip1,j,l) = pdvfi(1,j,l)
    650      ENDDO
     640  DO l = 1, llm
     641
     642    DO j = 2, jjm - 1
     643      ig0 = 1 + (j - 2) * iim
     644      DO i = 1, iim
     645        pdvfi(i, j, l) = &
     646                0.5 * (zdvfi(ig0 + i, l) + zdvfi(ig0 + i + iim, l)) * cv(i, j)
     647      ENDDO
     648      pdvfi(iip1, j, l) = pdvfi(1, j, l)
     649    ENDDO
    651650  ENDDO
    652651
     
    654653  !   68. champ v pres des poles:
    655654  !   ---------------------------
    656    ! v = U * cos(long) + V * SIN(long)
    657 
    658   DO l=1,llm
    659 
    660      DO i=1,iim
    661         pdvfi(i,1,l)= &
    662               zdufi(1,l)*COS(rlonv(i))+zdvfi(1,l)*SIN(rlonv(i))
    663         pdvfi(i,jjm,l)=zdufi(ngridmx,l)*COS(rlonv(i)) &
    664               +zdvfi(ngridmx,l)*SIN(rlonv(i))
    665         pdvfi(i,1,l)= &
    666               0.5*(pdvfi(i,1,l)+zdvfi(i+1,l))*cv(i,1)
    667         pdvfi(i,jjm,l)= &
    668               0.5*(pdvfi(i,jjm,l)+zdvfi(ngridmx-iip1+i,l))*cv(i,jjm)
    669       ENDDO
    670 
    671      pdvfi(iip1,1,l)  = pdvfi(1,1,l)
    672      pdvfi(iip1,jjm,l)= pdvfi(1,jjm,l)
     655  ! v = U * cos(long) + V * SIN(long)
     656
     657  DO l = 1, llm
     658
     659    DO i = 1, iim
     660      pdvfi(i, 1, l) = &
     661              zdufi(1, l) * COS(rlonv(i)) + zdvfi(1, l) * SIN(rlonv(i))
     662      pdvfi(i, jjm, l) = zdufi(ngridmx, l) * COS(rlonv(i)) &
     663              + zdvfi(ngridmx, l) * SIN(rlonv(i))
     664      pdvfi(i, 1, l) = &
     665              0.5 * (pdvfi(i, 1, l) + zdvfi(i + 1, l)) * cv(i, 1)
     666      pdvfi(i, jjm, l) = &
     667              0.5 * (pdvfi(i, jjm, l) + zdvfi(ngridmx - iip1 + i, l)) * cv(i, jjm)
     668    ENDDO
     669
     670    pdvfi(iip1, 1, l) = pdvfi(1, 1, l)
     671    pdvfi(iip1, jjm, l) = pdvfi(1, jjm, l)
    673672
    674673  ENDDO
  • LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat/gr_dyn_fi.f90

    r5116 r5119  
    1 
    21! $Header$
    32
    4 SUBROUTINE gr_dyn_fi(nfield,im,jm,ngrid,pdyn,pfi)
     3SUBROUTINE gr_dyn_fi(nfield, im, jm, ngrid, pdyn, pfi)
     4  USE lmdz_ssum_scopy, ONLY: scopy
     5
    56  IMPLICIT NONE
    67  !=======================================================================
     
    1213  !   -------------
    1314
    14   INTEGER :: im,jm,ngrid,nfield
    15   REAL :: pdyn(im,jm,nfield)
    16   REAL :: pfi(ngrid,nfield)
     15  INTEGER :: im, jm, ngrid, nfield
     16  REAL :: pdyn(im, jm, nfield)
     17  REAL :: pfi(ngrid, nfield)
    1718
    18   INTEGER :: j,ifield,ig
     19  INTEGER :: j, ifield, ig
    1920
    2021  !-----------------------------------------------------------------------
     
    2223  !   -------
    2324
    24   IF (ngrid/=2+(jm-2)*(im-1)) THEN
    25      CALL abort_gcm("gr_dyn_fi", 'probleme de dim', 1)
     25  IF (ngrid/=2 + (jm - 2) * (im - 1)) THEN
     26    CALL abort_gcm("gr_dyn_fi", 'probleme de dim', 1)
    2627  end if
    2728  !   traitement des poles
    28   CALL SCOPY(nfield,pdyn,im*jm,pfi,ngrid)
    29   CALL SCOPY(nfield,pdyn(1,jm,1),im*jm,pfi(ngrid,1),ngrid)
     29  CALL SCOPY(nfield, pdyn, im * jm, pfi, ngrid)
     30  CALL SCOPY(nfield, pdyn(1, jm, 1), im * jm, pfi(ngrid, 1), ngrid)
    3031
    3132  !   traitement des point normaux
    32   DO ifield=1,nfield
    33      DO j=2,jm-1
    34         ig=2+(j-2)*(im-1)
    35         CALL SCOPY(im-1,pdyn(1,j,ifield),1,pfi(ig,ifield),1)
    36      ENDDO
     33  DO ifield = 1, nfield
     34    DO j = 2, jm - 1
     35      ig = 2 + (j - 2) * (im - 1)
     36      CALL SCOPY(im - 1, pdyn(1, j, ifield), 1, pfi(ig, ifield), 1)
     37    ENDDO
    3738  ENDDO
    3839
  • LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat/gr_fi_dyn.f90

    r5105 r5119  
    1 
    21! $Header$
    32
    4 SUBROUTINE gr_fi_dyn(nfield,ngrid,im,jm,pfi,pdyn)
     3SUBROUTINE gr_fi_dyn(nfield, ngrid, im, jm, pfi, pdyn)
     4  USE lmdz_ssum_scopy, ONLY: scopy
     5
    56  IMPLICIT NONE
    67  !=======================================================================
     
    1213  !   -------------
    1314
    14   INTEGER :: im,jm,ngrid,nfield
    15   REAL :: pdyn(im,jm,nfield)
    16   REAL :: pfi(ngrid,nfield)
     15  INTEGER :: im, jm, ngrid, nfield
     16  REAL :: pdyn(im, jm, nfield)
     17  REAL :: pfi(ngrid, nfield)
    1718
    18   INTEGER :: i,j,ifield,ig
     19  INTEGER :: i, j, ifield, ig
    1920
    2021  !-----------------------------------------------------------------------
     
    2223  !   -------
    2324
    24   DO ifield=1,nfield
    25   !   traitement des poles
    26      DO i=1,im
    27         pdyn(i,1,ifield)=pfi(1,ifield)
    28         pdyn(i,jm,ifield)=pfi(ngrid,ifield)
    29      ENDDO
     25  DO ifield = 1, nfield
     26    !   traitement des poles
     27    DO i = 1, im
     28      pdyn(i, 1, ifield) = pfi(1, ifield)
     29      pdyn(i, jm, ifield) = pfi(ngrid, ifield)
     30    ENDDO
    3031
    31   !   traitement des point normaux
    32      DO j=2,jm-1
    33         ig=2+(j-2)*(im-1)
    34         CALL SCOPY(im-1,pfi(ig,ifield),1,pdyn(1,j,ifield),1)
    35         pdyn(im,j,ifield)=pdyn(1,j,ifield)
    36      ENDDO
     32    !   traitement des point normaux
     33    DO j = 2, jm - 1
     34      ig = 2 + (j - 2) * (im - 1)
     35      CALL SCOPY(im - 1, pfi(ig, ifield), 1, pdyn(1, j, ifield), 1)
     36      pdyn(im, j, ifield) = pdyn(1, j, ifield)
     37    ENDDO
    3738  ENDDO
    3839
  • LMDZ6/branches/Amaury_dev/libf/dynphy_lonlat/phylmd/test_disvert_m.F90

    r5117 r5119  
    33  IMPLICIT NONE
    44
    5 contains
     5CONTAINS
    66
    77  SUBROUTINE test_disvert
     
    6565  END SUBROUTINE  test_disvert
    6666
    67 end module test_disvert_m
     67END MODULE test_disvert_m
Note: See TracChangeset for help on using the changeset viewer.