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/phyetat0.F

    r1036 r1047  
    1       SUBROUTINE phyetat0 (fichnom,tab0,Lmodif,nsoil,nq,
     1      SUBROUTINE phyetat0 (fichnom,tab0,Lmodif,nsoil,ngrid,nlay,nq,
    22     .           day_ini,time0,
    33     .           tsurf,tsoil,emis,q2,qsurf,co2ice)
     
    55      use netcdf
    66      use infotrac, only: nqtot, tnom
     7      use surfdat_h, only: phisfi, albedodat, z0, z0_default,
     8     &                     zmea, zstd, zsig, zgam, zthe
    79
    810      implicit none
     
    1719c======================================================================
    1820!#include "netcdf.inc"
    19 #include "dimensions.h"
    20 #include "dimphys.h"
     21!#include "dimensions.h"
     22!#include "dimphys.h"
    2123!#include "comgeomfi.h"
    22 #include "surfdat.h"
     24!#include "surfdat.h"
    2325#include "planete.h"
    24 #include "dimradmars.h"
    25 #include "yomaer.h"
     26!#include "dimradmars.h"
     27!#include "yomaer.h"
    2628#include "comcstfi.h"
    2729!#include "tracer.h"
     
    3537!  ---------
    3638!  inputs:
    37       character*(*) fichnom ! "startfi.nc" file
    38       integer tab0
    39       integer Lmodif
    40       integer nsoil ! # of soil layers
    41       integer nq
    42       integer day_ini
    43       real time0
     39      character*(*),intent(in) :: fichnom ! "startfi.nc" file
     40      integer,intent(in) :: tab0
     41      integer,intent(in) :: Lmodif
     42      integer,intent(in) :: nsoil ! # of soil layers
     43      integer,intent(in) :: ngrid ! # of atmospheric columns
     44      integer,intent(in) :: nlay ! # of atmospheric layers
     45      integer,intent(in) :: nq
     46      integer :: day_ini
     47      real :: time0
    4448
    4549!  outputs:
    46       real tsurf(ngridmx) ! surface temperature
    47       real tsoil(ngridmx,nsoil) ! soil temperature
    48       real emis(ngridmx) ! surface emissivity
    49       real q2(ngridmx, llm+1) !
    50       real qsurf(ngridmx,nq) ! tracers on surface
    51       real co2ice(ngridmx) ! co2 ice cover
     50      real,intent(out) :: tsurf(ngrid) ! surface temperature
     51      real,intent(out) :: tsoil(ngrid,nsoil) ! soil temperature
     52      real,intent(out) :: emis(ngrid) ! surface emissivity
     53      real,intent(out) :: q2(ngrid,nlay+1) !
     54      real,intent(out) :: qsurf(ngrid,nq) ! tracers on surface
     55      real,intent(out) :: co2ice(ngrid) ! co2 ice cover
    5256
    5357!======================================================================
    5458!  Local variables:
    5559
    56       real surffield(ngridmx) ! to temporarily store a surface field
     60      real surffield(ngrid) ! to temporarily store a surface field
    5761      real xmin,xmax ! to display min and max of a field
    5862c
     
    248252         CALL abort
    249253      ENDIF
    250       xmin = 1.0E+20
    251       xmax = -1.0E+20
    252       DO i = 1, ngridmx
    253          xmin = MIN(zmea(i),xmin)
    254          xmax = MAX(zmea(i),xmax)
    255       ENDDO
     254      xmin = MINVAL(zmea)
     255      xmax = MAXVAL(zmea)
    256256      PRINT*,'<zmea>:', xmin, xmax
    257257c
     
    270270         CALL abort
    271271      ENDIF
    272       xmin = 1.0E+20
    273       xmax = -1.0E+20
    274       DO i = 1, ngridmx
    275          xmin = MIN(zstd(i),xmin)
    276          xmax = MAX(zstd(i),xmax)
    277       ENDDO
     272      xmin = MINVAL(zstd)
     273      xmax = MAXVAL(zstd)
    278274      PRINT*,'<zstd>:', xmin, xmax
    279275c
     
    292288         CALL abort
    293289      ENDIF
    294       xmin = 1.0E+20
    295       xmax = -1.0E+20
    296       DO i = 1, ngridmx
    297          xmin = MIN(zsig(i),xmin)
    298          xmax = MAX(zsig(i),xmax)
    299       ENDDO
     290      xmin = MINVAL(zsig)
     291      xmax = MAXVAL(zsig)
    300292      PRINT*,'<zsig>:', xmin, xmax
    301293c
     
    314306         CALL abort
    315307      ENDIF
    316       xmin = 1.0E+20
    317       xmax = -1.0E+20
    318       DO i = 1, ngridmx
    319          xmin = MIN(zgam(i),xmin)
    320          xmax = MAX(zgam(i),xmax)
    321       ENDDO
     308      xmin = MINVAL(zgam)
     309      xmax = MAXVAL(zgam)
    322310      PRINT*,'<zgam>:', xmin, xmax
    323311c
     
    336324         CALL abort
    337325      ENDIF
    338       xmin = 1.0E+20
    339       xmax = -1.0E+20
    340       DO i = 1, ngridmx
    341          xmin = MIN(zthe(i),xmin)
    342          xmax = MAX(zthe(i),xmax)
    343       ENDDO
     326      xmin = MINVAL(zthe)
     327      xmax = MAXVAL(zthe)
    344328      PRINT*,'<zthe>:', xmin, xmax
    345329     
     
    417401      corner(1)=1
    418402      corner(2)=indextime
    419       edges(1)=ngridmx
     403      edges(1)=ngrid
    420404      edges(2)=1
    421405      ierr=nf90_inq_varid(nid,"co2ice",nvarid)
     
    486470!         IF (nbsrf >= 2) THEN
    487471!            DO nsrf = 2, nbsrf
    488 !               DO i = 1, ngridmx
     472!               DO i = 1, ngrid
    489473!                  tsurf(i,nsrf) = tsurf(i,1)
    490474!               ENDDO
     
    506490!               PRINT*, "phyetat0: Le champ <tsoil> est absent"
    507491!               PRINT*, "          Il prend donc la valeur de surface"
    508 !               DO i=1, ngridmx
     492!               DO i=1, ngrid
    509493!                  tsoil(i,isoil,nsrf)=tsurf(i,nsrf)
    510494!               ENDDO
     
    549533
    550534!
    551 ! surface roughness length (NB: z0 is a common in surfdat.h)
     535! surface roughness length (NB: z0 is a common in surfdat_h)
    552536!
    553537      ierr=nf90_inq_varid(nid,"z0",nvarid)
     
    577561      corner(2)=1
    578562      corner(3)=indextime
    579       edges(1)=ngridmx
    580       edges(2)=llm+1
     563      edges(1)=ngrid
     564      edges(2)=nlay+1
    581565      edges(3)=1
    582566      ierr=nf90_inq_varid(nid,"q2",nvarid)
     
    602586      corner(1)=1
    603587      corner(2)=indextime
    604       edges(1)=ngridmx
     588      edges(1)=ngrid
    605589      edges(2)=1
    606590      IF(nq.GE.1) THEN
     
    628612     &                  ' not found in file'
    629613             write(*,*) trim(txt), ' set to 0'
    630              do ig=1,ngridmx
     614             do ig=1,ngrid
    631615               qsurf(ig,iq)=0.
    632616             end do
     
    635619           !ierr=nf90_get_var(nid,nvarid,qsurf(1,iq))
    636620           ierr=nf90_get_var(nid,nvarid,surffield,corner,edges)
    637            qsurf(1:ngridmx,iq)=surffield(1:ngridmx)
     621           qsurf(1:ngrid,iq)=surffield(1:ngrid)
    638622             IF (ierr.NE.nf90_noerr) THEN
    639623               PRINT*, 'phyetat0: Lecture echouee pour <',trim(txt),'>'
     
    644628           xmin = 1.0E+20
    645629           xmax = -1.0E+20
    646            xmin = MINVAL(qsurf(1:ngridmx,iq))
    647            xmax = MAXVAL(qsurf(1:ngridmx,iq))
     630           xmin = MINVAL(qsurf(1:ngrid,iq))
     631           xmax = MAXVAL(qsurf(1:ngrid,iq))
    648632           PRINT*,'tracer on surface <',trim(txt),'>:',xmin,xmax
    649633         ENDDO
     
    660644!            if (yes.eq.'y') then
    661645!              write(*,*) 'OK, let s reindex qsurf'
    662 !                 do ig=1,ngridmx
     646!                 do ig=1,ngrid
    663647!                    do iq=nqtot,nqtot-nqold+1,-1
    664648!                       qsurf(ig,iq)=qsurf(ig,iq-nqtot+nqold)
     
    675659! as well as thermal inertia and volumetric heat capacity
    676660
    677       call soil_settings(nid,ngridmx,nsoil,tsurf,tsoil,indextime)
     661      call soil_settings(nid,ngrid,nsoil,tsurf,tsoil,indextime)
    678662c
    679663c Fermer le fichier:
Note: See TracChangeset for help on using the changeset viewer.