Ignore:
Timestamp:
Sep 23, 2013, 9:56:47 AM (11 years ago)
Author:
emillour
Message:

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/vlz_fi.F

    r38 r1047  
    1       SUBROUTINE vlz_fi(ngrid,q,pente_max,masse,w,wq)
     1      SUBROUTINE vlz_fi(ngrid,nlay,q,pente_max,masse,w,wq)
    22c
    33c     Auteurs:   P.Le Van, F.Hourdin, F.Forget
     
    1616      IMPLICIT NONE
    1717c
    18 #include "dimensions.h"
    19 #include "dimphys.h"
     18!#include "dimensions.h"
     19!#include "dimphys.h"
    2020
    2121c
     
    2323c   Arguments:
    2424c   ----------
    25       integer ngrid
    26       real masse(ngrid,llm),pente_max
    27       REAL q(ngrid,llm)
    28       REAL w(ngrid,llm)
    29       REAL wq(ngrid,llm+1)
     25      integer,intent(in) :: ngrid ! number of atmospheric columns
     26      integer,intent(in) :: nlay ! number of atmospheric layers
     27      real masse(ngrid,nlay),pente_max
     28      REAL q(ngrid,nlay)
     29      REAL w(ngrid,nlay)
     30      REAL wq(ngrid,nlay+1)
    3031c
    3132c      Local
     
    3536c
    3637
    37       real dzq(ngridmx,llm),dzqw(ngridmx,llm),adzqw(ngridmx,llm),dzqmax
     38      real dzq(ngrid,nlay),dzqw(ngrid,nlay),adzqw(ngrid,nlay),dzqmax
    3839      real newmasse
    3940      real sigw, Mtot, MQtot
     
    4748c    sens de W
    4849
    49       do l=2,llm
     50      do l=2,nlay
    5051         do ij=1,ngrid
    5152            dzqw(ij,l)=q(ij,l-1)-q(ij,l)
     
    5455      enddo
    5556
    56       do l=2,llm-1
     57      do l=2,nlay-1
    5758         do ij=1,ngrid
    5859#ifdef CRAY
     
    7374      do ij=1,ngrid
    7475         dzq(ij,1)=0.
    75          dzq(ij,llm)=0.
     76         dzq(ij,nlay)=0.
    7677      enddo
    7778c ---------------------------------------------------------------
     
    8384c      No flux at the model top:
    8485       do ij=1,ngrid
    85           wq(ij,llm+1)=0.
     86          wq(ij,nlay+1)=0.
    8687       enddo
    8788
     
    8990c      ===============================
    9091
    91        do l = 1,llm          ! loop different than when w<0
     92       do l = 1,nlay          ! loop different than when w<0
    9293        do  ij = 1,ngrid
    9394
     
    107108            Mtot = masse(ij,m)
    108109            MQtot = masse(ij,m)*q(ij,m)
    109             if(m.ge.llm)goto 88
     110            if(m.ge.nlay)goto 88
    110111            do while(w(ij,l).gt.(Mtot+masse(ij,m+1)))
    111112                m=m+1
    112113                Mtot = Mtot + masse(ij,m)
    113114                MQtot = MQtot + masse(ij,m)*q(ij,m)
    114                 if(m.ge.llm)goto 88
     115                if(m.ge.nlay)goto 88
    115116            end do
    116117 88         continue
    117             if (m.lt.llm) then
     118            if (m.lt.nlay) then
    118119                sigw=(w(ij,l)-Mtot)/masse(ij,m+1)
    119120                wq(ij,l)=(MQtot + (w(ij,l)-Mtot)*
     
    137138       end do
    138139
    139        do l = 1,llm-1  ! loop different than when w>0
     140       do l = 1,nlay-1  ! loop different than when w>0
    140141        do  ij = 1,ngrid
    141142         if(w(ij,l+1).le.0)then
     
    176177 99    continue
    177178
    178       do l=1,llm
     179      do l=1,nlay
    179180         do ij=1,ngrid
    180181
Note: See TracChangeset for help on using the changeset viewer.