Ignore:
Timestamp:
May 5, 2023, 3:15:04 PM (19 months ago)
Author:
romain.vande
Message:

Mars PCM :
Correct start2archive to write watercaptag correctly.
Watercaptag will be set to false and correctly handle by the PCM in the case where we change resolution.
+ Correct inertiesoil writting in start2archive
RV

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/lect_start_archive.F

    r2943 r2959  
    2626     &             subslope_dist,end_comslope_h,ini_comslope_h
    2727      use netcdf
     28      use surfdat_h, ONLY: watercaptag
    2829      implicit none
    2930
     
    139140      real totcloudfracS(iip1,jjp1)
    140141      real watercapS(iip1,jjp1,nslope)
     142      real watercaptagS(iip1,jjp1)
    141143      real albedoS(iip1,jjp1,nslope)
    142144
     
    177179      real, dimension(:,:,:), allocatable :: watercapold
    178180      real, dimension(:,:), allocatable :: watercapold_noslope
     181      real, dimension(:,:), allocatable :: watercaptagold
     182      real, dimension(:), allocatable :: watercaptag_tmp
    179183      real, dimension(:,:,:), allocatable :: albedoold
    180184      real, dimension(:,:), allocatable :: albedoold_noslope
     
    426430      allocate(totcloudfracold(imold+1,jmold+1))
    427431      allocate(watercapold(imold+1,jmold+1,nslope))
     432      allocate(watercaptagold(imold+1,jmold+1))
     433      allocate(watercaptag_tmp(ngrid))
    428434      allocate(albedoold(imold+1,jmold+1,nslope))
    429435     
     
    908914        IF (ierr .NE. NF_NOERR) THEN
    909915           PRINT*, "lect_start_archive: Failed loading <watercap>"
     916           PRINT*, NF_STRERROR(ierr)
     917           CALL abort
     918        ENDIF
     919      ENDIF
     920
     921c
     922      ierr = NF_INQ_VARID (nid, "watercaptag", nvarid)
     923      IF (ierr .NE. NF_NOERR) THEN
     924         PRINT*, "lect_start_archive: <watercaptag> not in file"
     925         PRINT*, "watercaptag set to false, will be adapted in
     926     &            surfini of PCM"
     927         watercaptagold(:,:) = 0
     928      ELSE
     929        ierr = nf90_get_var(nid, nvarid,watercaptagold)
     930        IF (ierr .NE. NF_NOERR) THEN
     931           PRINT*, "lect_start_archive: Failed loading <watercaptag>"
    910932           PRINT*, NF_STRERROR(ierr)
    911933           CALL abort
     
    12301252     &       watercap(:,1))
    12311253
     1254c Watercaptag
     1255      if(imold.eq.iim .and. jmold.eq.jjm) then
     1256      else
     1257       print *, "We are doing an horizontal interpolation,
     1258     &    watercaptag will be set to false and redefined proprely
     1259     &    in the PCM(surfini)"
     1260       watercaptagold(:,:)=0
     1261      endif
     1262        call interp_horiz (watercaptagold(:,:),watercaptags(:,:),
     1263     &       imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv)
     1264        call gr_dyn_fi (1,iim+1,jjm+1,ngrid,watercaptags(:,:),
     1265     &       watercaptag_tmp(:))
     1266
     1267        do i=1,ngrid
     1268          if(watercaptag_tmp(i).gt. 0.5) then
     1269            watercaptag(i)=.true.
     1270          else
     1271            watercaptag(i)=.false.
     1272          endif
     1273        enddo
     1274
     1275
    12321276c Surface albedo
    12331277      call interp_horiz (albedoold(:,:,1),albedoS(:,:,1),
     
    16051649      deallocate(totcloudfracold)
    16061650      deallocate(watercapold)
     1651      deallocate(watercaptagold)
    16071652      deallocate(albedoold)
    16081653      deallocate(var,varp1)
Note: See TracChangeset for help on using the changeset viewer.