Ignore:
Timestamp:
Oct 24, 2024, 1:55:38 PM (3 days ago)
Author:
abarral
Message:

Replace F77 netcdf library by F90 netcdf library

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ6/trunk/libf/dyn3d/guide_mod.f90

    r5268 r5270  
    1414                    nf90_inq_dimid, nf90_inquire_dimension
    1515  use pres2lev_mod, only: pres2lev
     16  USE netcdf, ONLY: nf90_nowrite, nf90_open, nf90_inq_varid, nf90_close, &
     17          nf90_inq_dimid, nf90_inquire_dimension, nf90_float, nf90_def_var, &
     18          nf90_create, nf90_def_dim, nf90_open, nf90_unlimited, nf90_write, nf90_enddef, nf90_redef, &
     19          nf90_close, nf90_inq_varid, nf90_get_var, nf90_noerr, nf90_clobber, &
     20          nf90_64bit_offset, nf90_inq_dimid, nf90_inquire_dimension, nf90_put_var
    1621
    1722  IMPLICIT NONE
     
    7782    INCLUDE "dimensions.h"
    7883    INCLUDE "paramet.h"
    79     INCLUDE "netcdf.inc"
    8084
    8185    INTEGER                :: error,ncidpl,rid,rcod
     
    224228
    225229    endif
    226     error=NF_INQ_DIMID(ncidpl,'LEVEL',rid)
    227     IF (error.NE.NF90_NOERR) error=NF_INQ_DIMID(ncidpl,'PRESSURE',rid)
     230    error=nf90_inq_dimid(ncidpl,'LEVEL',rid)
     231    IF (error.NE.NF90_NOERR) error=nf90_inq_dimid(ncidpl,'PRESSURE',rid)
    228232    IF (error.NE.NF90_NOERR) THEN
    229233        CALL abort_gcm(modname,'Nudging: error reading pressure levels',1)
    230234    ENDIF
    231     error=NF_INQ_DIMLEN(ncidpl,rid,nlevnc)
     235    error=nf90_inquire_dimension(ncidpl,rid,len=nlevnc)
    232236    write(*,*)trim(modname)//' : number of vertical levels nlevnc', nlevnc
    233237    rcod = nf90_close(ncidpl)
     
    16601664    INCLUDE "dimensions.h"
    16611665    INCLUDE "paramet.h"
    1662     INCLUDE "netcdf.inc"
    16631666    INCLUDE "comgeom2.h"
    16641667   
     
    16861689! ----------------------------------------------
    16871690! Ouverture du fichier
    1688         ierr=NF_CREATE("guide_ins.nc",IOR(NF_CLOBBER,NF_64BIT_OFFSET),nid)
     1691        ierr=nf90_create("guide_ins.nc",IOR(nf90_clobber,nf90_64bit_offset),nid)
    16891692! Definition des dimensions
    1690         ierr=NF_DEF_DIM(nid,"LONU",iip1,id_lonu)
    1691         ierr=NF_DEF_DIM(nid,"LONV",iip1,id_lonv)
    1692         ierr=NF_DEF_DIM(nid,"LATU",jjp1,id_latu)
    1693         ierr=NF_DEF_DIM(nid,"LATV",jjm,id_latv)
    1694         ierr=NF_DEF_DIM(nid,"LEVEL",llm,id_lev)
    1695         ierr=NF_DEF_DIM(nid,"TIME",NF_UNLIMITED,id_tim)
     1693        ierr=nf90_def_dim(nid,"LONU",iip1,id_lonu)
     1694        ierr=nf90_def_dim(nid,"LONV",iip1,id_lonv)
     1695        ierr=nf90_def_dim(nid,"LATU",jjp1,id_latu)
     1696        ierr=nf90_def_dim(nid,"LATV",jjm,id_latv)
     1697        ierr=nf90_def_dim(nid,"LEVEL",llm,id_lev)
     1698        ierr=nf90_def_dim(nid,"TIME",nf90_unlimited,id_tim)
    16961699
    16971700! Creation des variables dimensions
     
    17101713             varid_alpha_q)
    17111714       
    1712         ierr=NF_ENDDEF(nid)
     1715        ierr=nf90_enddef(nid)
    17131716
    17141717! Enregistrement des variables dimensions
     
    17301733! Cr�ation des variables sauvegard�es
    17311734! --------------------------------------------------------------------
    1732         ierr = NF_REDEF(nid)
     1735        ierr = nf90_redef(nid)
    17331736! Pressure (GCM)
    17341737        dim4=(/id_lonv,id_latu,id_lev,id_tim/)
     
    17641767        ENDIF
    17651768       
    1766         ierr = NF_ENDDEF(nid)
    1767         ierr = NF_CLOSE(nid)
     1769        ierr = nf90_enddef(nid)
     1770        ierr = nf90_close(nid)
    17681771    ENDIF ! timestep=0
    17691772
     
    17711774! Enregistrement du champ
    17721775! --------------------------------------------------------------------
    1773     ierr=NF_OPEN("guide_ins.nc",NF_WRITE,nid)
     1776    ierr=nf90_open("guide_ins.nc",nf90_write,nid)
    17741777
    17751778    IF (varname=="SP") timestep=timestep+1
    17761779
    1777     ierr = NF_INQ_VARID(nid,varname,varid)
     1780    ierr = nf90_inq_varid(nid,varname,varid)
    17781781    SELECT CASE (varname)
    17791782    CASE ("SP","ps")
     
    18001803
    18011804    ierr = nf90_put_var(nid, varid, field2, start, count)
    1802     ierr = NF_CLOSE(nid)
     1805    ierr = nf90_close(nid)
    18031806
    18041807  END SUBROUTINE guide_out
Note: See TracChangeset for help on using the changeset viewer.