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/iniwritesoil.F90

    r38 r1047  
    1 subroutine iniwritesoil(nid)
     1subroutine iniwritesoil(nid,ngrid)
    22
    33! initialization routine for 'writediagoil'. Here we create/define
     
    55! (time-independent) parameters.
    66
     7use comsoil_h, only: mlayer, inertiedat, nsoilmx
     8
    79implicit none
    810
    911#include"dimensions.h"
    10 #include"dimphys.h"
     12!#include"dimphys.h"
    1113#include"paramet.h"
    1214#include"comcstfi.h"
    1315#include"comgeom.h"
    14 #include"comsoil.h"
     16!#include"comsoil.h"
    1517#include"netcdf.inc"
    1618
    1719! Arguments:
     20integer,intent(in) :: ngrid
    1821integer,intent(in) :: nid ! NetCDF output file ID
    1922
     
    154157ierr=NF_PUT_VAR_REAL(nid,varid,mlayer)
    155158#endif
    156 ! Note mlayer(0:nsoilmx-1) known from comsoil.h
     159! Note mlayer(0:nsoilmx-1) known from comsoil_h
    157160if (ierr.ne.NF_NOERR) then
    158161  write(*,*)"iniwritesoil: Error, could not write depth variable"
     
    240243
    241244! Recast data along 'dynamics' grid
    242 ! Note: inertiedat is known from comsoil.h
     245! Note: inertiedat is known from comsoil_h
    243246
    244247do l=1,nsoilmx
     
    246249  do i=1,iip1
    247250    data3(i,1,l)=inertiedat(1,l)
    248     data3(i,jjp1,l)=inertiedat(ngridmx,l)
     251    data3(i,jjp1,l)=inertiedat(ngrid,l)
    249252  enddo
    250253  ! rest of the grid
Note: See TracChangeset for help on using the changeset viewer.