Ignore:
Timestamp:
Apr 5, 2016, 10:51:51 AM (9 years ago)
Author:
emillour
Message:

Generic GCM: Towards a cleaner separation between physics and dynamics

  • Got rid of references to "dimensions.h" from physics packages: use nbp_lon (=iim), nbp_lat (=jjp1) and nbp_lev from module mod_grid_phy_lmdz (in phy_common) instead.
  • Removed module "comhdiff_mod.F90", as it is only used by module surf_heat_transp_mod.F90, moved module variables there.
  • Addedin "surf_heat_transp_mod" local versions of some arrays and routines (from dyn3d) required to compute gradient, divergence, etc. on the global dynamics grid. As before, the slab ocean only works in serial.

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/iniwritesoil.F90

    r1384 r1529  
    1 subroutine iniwritesoil(nid,ngrid)
     1subroutine iniwritesoil(nid,ngrid,inertia,area)
    22
    33! initialization routine for 'writediagoil'. Here we create/define
     
    55! (time-independent) parameters.
    66
    7 use comsoil_h, only: mlayer, inertiedat, nsoilmx
    8 use comcstfi_mod, only: pi
     7use comsoil_h, only: mlayer, nsoilmx
     8USE comcstfi_mod, only: pi
     9USE regular_lonlat_mod, ONLY: lon_reg, lat_reg
     10use mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
    911
    1012implicit none
    1113
    12 #include"dimensions.h"
    13 #include"paramet.h"
    14 #include"comgeom.h"
    15 #include"netcdf.inc"
     14include"netcdf.inc"
    1615
    1716! Arguments:
    1817integer,intent(in) :: ngrid
    1918integer,intent(in) :: nid ! NetCDF output file ID
     19real,intent(in) :: inertia(nbp_lon+1,nbp_lat,nsoilmx)
     20real,intent(in) :: area(nbp_lon+1,nbp_lat) ! mesh area (m2)
    2021
    2122! Local variables:
     
    3031integer,dimension(3) :: dimids ! to store IDs of dimensions of a variable
    3132character(len=60) :: text ! to store some text
    32 real,dimension(iip1,jjp1,nsoilmx) :: data3 ! to store 3D data
     33real,dimension(nbp_lon+1,nbp_lat,nsoilmx) :: data3 ! to store 3D data
    3334integer :: i,j,l,ig0
     35real :: lon_reg_ext(nbp_lon+1) ! extended longitudes
    3436
    3537! 1. Define the dimensions
     
    3840
    3941! Define the dimensions
    40 ierr=NF_DEF_DIM(nid,"longitude",iip1,idim_rlonv)
    41 ! iip1 known from paramet.h
     42ierr=NF_DEF_DIM(nid,"longitude",nbp_lon+1,idim_rlonv)
    4243if (ierr.ne.NF_NOERR) then
    4344  write(*,*)"iniwritesoil: Error, could not define longitude dimension"
    4445endif
    45 ierr=NF_DEF_DIM(nid,"latitude",jjp1,idim_rlatu)
    46 ! jjp1 known from paramet.h
     46ierr=NF_DEF_DIM(nid,"latitude",nbp_lat,idim_rlatu)
    4747if (ierr.ne.NF_NOERR) then
    4848  write(*,*)"iniwritesoil: Error, could not define latitude dimension"
    4949endif
    5050ierr=NF_DEF_DIM(nid,"depth",nsoilmx,idim_depth)
    51 ! nsoilmx known from dimphys.h
     51! nsoilmx known from comsoil_h
    5252if (ierr.ne.NF_NOERR) then
    5353  write(*,*)"iniwritesoil: Error, could not define depth dimension"
     
    8181ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
    8282
     83lon_reg_ext(1:nbp_lon)=lon_reg(1:nbp_lon)
     84!add extra redundant point (180 degrees, since lon_reg starts at -180
     85lon_reg_ext(nbp_lon+1)=-lon_reg_ext(1)
     86
    8387! Write longitude to file
    8488ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
    8589! Write
    8690#ifdef NC_DOUBLE
    87 ierr=NF_PUT_VAR_DOUBLE(nid,varid,rlonv*(180./pi))
    88 #else
    89 ierr=NF_PUT_VAR_REAL(nid,varid,rlonv*(180./pi))
    90 #endif
    91 ! Note: rlonv is known from comgeom.h
     91ierr=NF_PUT_VAR_DOUBLE(nid,varid,lon_reg_ext*(180./pi))
     92#else
     93ierr=NF_PUT_VAR_REAL(nid,varid,lon_reg_ext*(180./pi))
     94#endif
     95! Note: rlonv is known from comgeom.h and pi from comcstfi.h
    9296if (ierr.ne.NF_NOERR) then
    9397  write(*,*)"iniwritesoil: Error, could not write longitude variable"
     
    117121! Write
    118122#ifdef NC_DOUBLE
    119 ierr=NF_PUT_VAR_DOUBLE(nid,varid,rlatu*(180./pi))
    120 #else
    121 ierr=NF_PUT_VAR_REAL(nid,varid,rlatu*(180./pi))
    122 #endif
    123 ! Note: rlatu is known from comgeom.h
     123ierr=NF_PUT_VAR_DOUBLE(nid,varid,lat_reg*(180./pi))
     124#else
     125ierr=NF_PUT_VAR_REAL(nid,varid,lat_reg*(180./pi))
     126#endif
    124127if (ierr.ne.NF_NOERR) then
    125128  write(*,*)"iniwritesoil: Error, could not write longitude variable"
     
    209212! Write
    210213#ifdef NC_DOUBLE
    211 ierr=NF_PUT_VAR_DOUBLE(nid,varid,aire)
    212 #else
    213 ierr=NF_PUT_VAR_REAL(nid,varid,aire)
    214 #endif
    215 ! Note: aire is known from comgeom.h
     214ierr=NF_PUT_VAR_DOUBLE(nid,varid,area)
     215#else
     216ierr=NF_PUT_VAR_REAL(nid,varid,area)
     217#endif
    216218if (ierr.ne.NF_NOERR) then
    217219  write(*,*)"iniwritesoil: Error, could not write area variable"
     
    240242ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
    241243
    242 ! Recast data along 'dynamics' grid
    243 ! Note: inertiedat is known from comsoil_h
    244 
    245 do l=1,nsoilmx
    246   ! handle the poles
    247   do i=1,iip1
    248     data3(i,1,l)=inertiedat(1,l)
    249     data3(i,jjp1,l)=inertiedat(ngrid,l)
    250   enddo
    251   ! rest of the grid
    252   do j=2,jjm
    253     ig0=1+(j-2)*iim
    254     do i=1,iim
    255       data3(i,j,l)=inertiedat(ig0+i,l)
    256     enddo
    257     data3(iip1,j,l)=data3(1,j,l) ! extra (modulo) longitude
    258   enddo
    259 enddo ! of do l=1,nsoilmx
    260 
    261244! Write data2 to file
    262245ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
    263246! Write
    264247#ifdef NC_DOUBLE
    265 ierr=NF_PUT_VAR_DOUBLE(nid,varid,data3)
    266 #else
    267 ierr=NF_PUT_VAR_REAL(nid,varid,data3)
     248ierr=NF_PUT_VAR_DOUBLE(nid,varid,inertia)
     249#else
     250ierr=NF_PUT_VAR_REAL(nid,varid,inertia)
    268251#endif
    269252if (ierr.ne.NF_NOERR) then
Note: See TracChangeset for help on using the changeset viewer.