Ignore:
Timestamp:
Mar 19, 2012, 11:27:18 AM (13 years ago)
Author:
emillour
Message:

Generic GCM:
Some cleanup and bug fixing:

  • "cloudfrac" was not well written to restartfi (wrong size).
  • missing save attribute for "reffrad" in physiq.F90.
  • cleanup recomputation of surface pressure in newstart and change loop order in interp_horiz (which "fixes" an odd behaviour which fills some arrays with zeros, but only when using some versions of ifort!)

EM

Location:
trunk/LMDZ.GENERIC/libf/phystd
Files:
7 edited

Legend:

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

    r538 r588  
    8080         R     = 8.314511E+0 *1000.E+0/mugaz
    8181         rcp   = R/cpp
    82       elseif((cpp.ne.cpp_c) .or. (mugaz.ne.mugaz_c))then
     82!      elseif((cpp.ne.cpp_c) .or. (mugaz.ne.mugaz_c))then
     83         ! Ehouarn: tolerate a small mismatch between computed/stored values
     84       elseif((abs(1.-cpp/cpp_c).gt.1.e-6) .or.  &
     85              (abs(1.-mugaz/mugaz_c).gt.1.e-6)) then
    8386         if(check_cpp_match)then
    8487            print*,'Values do not match!'
  • trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90

    r586 r588  
    522522               temp=aerosol(ig,L_NLAYRAD-k,iaer)/QREFvis3d(ig,L_NLAYRAD-k,iaer)
    523523
    524                tauaero(2*k+2,iaer)=max(temp*pweight,0.0)
    525                tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.0)
     524               tauaero(2*k+2,iaer)=max(temp*pweight,0.d0)
     525               tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.d0)
    526526!
    527527            end do
  • trunk/LMDZ.GENERIC/libf/phystd/datareadnc.F

    r374 r588  
    4343
    4444      use datafile_mod, only: datadir
     45! to use  'getin'
     46      USE ioipsl_getincom
    4547      implicit none
    4648
     
    7779      REAL        zthe(imdp1*jmdp1)
    7880
    79       INTEGER   lnblnk, ierr
    80       EXTERNAL    lnblnk
     81      INTEGER   ierr
    8182
    8283      INTEGER   unit,nvarid
     
    110111c    Lecture NetCDF des donnees latitude et longitude
    111112c-----------------------------------------------------------------------
    112 
     113!    First get the corret value of dtadir, if not already done:
     114      ! default 'datadir' is set in "datadir_mod"
     115      call getin("datadir",datadir)
    113116      write(*,*) 'Choice of surface data is:'
    114117      filestring='ls '//trim(datadir)//'/*.nc'
  • trunk/LMDZ.GENERIC/libf/phystd/def_var.F90

    r135 r588  
    2727ierr=NF_REDEF(nid)
    2828
    29 print*,'in def_var.F90, dimids='
    30 print*,dimids
     29!print*,'in def_var.F90, dimids='
     30!print*,dimids
    3131
    3232! 2. Define the variable
  • trunk/LMDZ.GENERIC/libf/phystd/physdem1.F

    r253 r588  
    4848      REAL day_ini
    4949      INTEGER nsoil,nq
    50       integer ierr,idim1,idim2,idim3,idim4,idim5,nvarid
     50      integer ierr,idim1,idim2,idim3,idim4,idim5,idim6,nvarid
    5151
    5252      REAL phystep,time
     
    8282      real pzsig(ngridmx),pzgam(ngridmx),pzthe(ngridmx)
    8383      integer ig
    84       INTEGER lnblnk
    85       EXTERNAL lnblnk
    8684
    8785! flag which identifies if we are using old tracer names (qsurf01,...)
     
    144142      endif
    145143c
    146 !      ierr = NF_DEF_DIM (nid,"nlayer+1",llm+1,idim4)
    147144      ierr = NF_DEF_DIM (nid,"nlayer_plus_1",llm+1,idim4)
    148145      if (ierr.ne.NF_NOERR) then
     
    157154        WRITE(6,*)' nq = ',nq,' and ierr = ', ierr
    158155        write(6,*) NF_STRERROR(ierr)
     156      endif
     157
     158      ierr = NF_DEF_DIM (nid,"nlayer",llm,idim6)
     159      if (ierr.ne.NF_NOERR) then
     160        WRITE(6,*)'physdem1: Problem defining nlayer dimension'
     161        write(6,*) NF_STRERROR(ierr)
     162        call abort
    159163      endif
    160164
     
    448452c Write the physical fields
    449453
    450 ! CO2 Ice Cover
    451 !
    452 !      ierr = NF_REDEF (nid)
    453 !#ifdef NC_DOUBLE
    454 !      ierr = NF_DEF_VAR (nid, "co2ice", NF_DOUBLE, 1, idim2,nvarid)
    455 !#else
    456 !      ierr = NF_DEF_VAR (nid, "co2ice", NF_FLOAT, 1, idim2,nvarid)
    457 !#endif
    458 !      ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 13,
    459 !     .                        "CO2 ice cover")
    460 !      ierr = NF_ENDDEF(nid)
    461 !#ifdef NC_DOUBLE
    462 !      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,co2ice)
    463 !#else
    464 !      ierr = NF_PUT_VAR_REAL (nid,nvarid,co2ice)
    465 !#endif
    466 
    467 
    468 
    469 
    470454
    471455! Soil Thermal inertia
     
    574558      ierr = NF_REDEF (nid)
    575559#ifdef NC_DOUBLE
    576       ierr = NF_DEF_VAR (nid, "cloudfrac", NF_DOUBLE, 1, idim2,nvarid)
    577 #else
    578       ierr = NF_DEF_VAR (nid, "cloudfrac", NF_FLOAT, 1, idim2,nvarid)
     560      ierr = NF_DEF_VAR (nid, "cloudfrac", NF_DOUBLE, 2,
     561     & (/idim2,idim6/),nvarid)
     562#else
     563      ierr = NF_DEF_VAR (nid, "cloudfrac", NF_FLOAT, 2,
     564     & (/idim2,idim6/),nvarid)
    579565#endif
    580566      ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 14,
  • trunk/LMDZ.GENERIC/libf/phystd/physiq.F90

    r586 r588  
    209209      real latvl1,lonvl1             ! Viking Lander 1 point (for diagnostic)
    210210
    211       real reffrad(ngridmx,nlayermx,naerkind)                ! aerosol effective radius (m)
     211      real,save :: reffrad(ngridmx,nlayermx,naerkind) ! aerosol effective radius (m)
    212212
    213213!     Tendencies due to various processes:
    214214      real dqsurf(ngridmx,nqmx)
    215       real zdtlw(ngridmx,nlayermx)                            ! (K/s)
    216       real zdtsw(ngridmx,nlayermx)                            ! (K/s)
    217       save zdtlw, zdtsw
     215      real,save :: zdtlw(ngridmx,nlayermx)                    ! (K/s)
     216      real,save :: zdtsw(ngridmx,nlayermx)                    ! (K/s)
    218217
    219218      real cldtlw(ngridmx,nlayermx)                           ! (K/s) LW heating rate for clear areas
     
    17241723            if((ngridmx.gt.1) .or. (mod(icount-1,nint(ecritphy)).eq.0))then
    17251724               ! to record global radiative balance
    1726                open(92,file="rad_bal.out",form='formatted',access='append')
     1725               open(92,file="rad_bal.out",form='formatted',position='append')
    17271726               write(92,*) zday,ISR/totarea,ASR/totarea,OLR/totarea
    17281727               close(92)
    1729                open(93,file="tem_bal.out",form='formatted',access='append')
     1728               open(93,file="tem_bal.out",form='formatted',position='append')
    17301729               write(93,*) zday,Ts1,Ts2,Ts3,TsS
    17311730               close(93)
     
    17391738!     ------------------------------------------------------------------
    17401739      if(testradtimes)then
    1741          open(38,file="tau_phys.out",form='formatted',access='append')
     1740         open(38,file="tau_phys.out",form='formatted',position='append')
    17421741         ig=1
    17431742         do l=1,nlayer
     
    18131812            if((ngridmx.gt.1) .or. (mod(icount-1,nint(ecritphy)).eq.0))then
    18141813               ! to record global water balance
    1815                open(98,file="h2o_bal.out",form='formatted',access='append')
     1814               open(98,file="h2o_bal.out",form='formatted',position='append')
    18161815               write(98,*) zday,icesrf/totarea,liqsrf/totarea,icecol/totarea,vapcol/totarea
    18171816               close(98)
  • trunk/LMDZ.GENERIC/libf/phystd/tabfi.F

    r253 r588  
    176176      mugaz = tab_cntrl(tab0+8)
    177177      rcp = tab_cntrl(tab0+9)
     178      cpp=(8.314511/(mugaz/1000.0))/rcp
    178179      daysec = tab_cntrl(tab0+10)
    179180      dtphys = tab_cntrl(tab0+11)
     
    502503          write(*,*)
    503504          write(*,*) ' mugaz (new value):',mugaz
    504           r=8.314/(mugaz/1000.0)
     505          r=8.314511/(mugaz/1000.0)
    505506          write(*,*) ' R (new value):',r
    506507
     
    512513          write(*,*)
    513514          write(*,*) ' rcp (new value):',rcp
     515          r=8.314511/(mugaz/1000.0)
    514516          cpp=r/rcp
    515517          write(*,*) ' cpp (new value):',cpp
Note: See TracChangeset for help on using the changeset viewer.