Ignore:
Timestamp:
Mar 26, 2026, 9:59:09 AM (6 weeks ago)
Author:
jbclement
Message:

Mars PCM:
r3903 made initialization choice of 'watercaptag' independent from its content read in 'startfi.nc' when 'icelocationmode = 5' (default).
This revision substitutes this obsolete content-based initialization when 'watercaptag' is not found in "start_archive.nc" or when doing horizontal interpolation for an explicit initialization with "locate_watercaptag".
JBC

File:
1 edited

Legend:

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

    r3902 r4155  
    1818c
    1919c   Objet:     Lecture des variables d'un fichier "start_archive"
    20 c              Plus besoin de régler ancienne valeurs grace
     20c              Plus besoin de regler ancienne valeurs grace
    2121c              a l'allocation dynamique de memoire (Yann Wanherdrick)
    2222c
     
    3232     &             subslope_dist,end_comslope_h,ini_comslope_h
    3333      use interp_line_mod, only: interp_line
     34      use surfini_mod, only: locate_watercaptag
     35      use geometry_mod, only: longitude, latitude
    3436      use netcdf
    3537      use surfdat_h, ONLY: watercaptag
     
    232234      integer :: nslope_read
    233235      logical :: no_slope=.false.
     236      logical :: no_watercaptag = .false.
    234237      integer :: ndims
    235238      integer, dimension(:), allocatable :: dimids
     
    955958      IF (ierr .NE. NF_NOERR) THEN
    956959         PRINT*, "lect_start_archive: <watercaptag> not in file"
    957          PRINT*, "watercaptag set to false, will be adapted in
    958      &            surfini of PCM"
    959          watercaptagold(:,:) = 0
     960         PRINT*, "'watercaptag' will be initialized with
     961     & ''locate_watercaptag''"
     962         no_watercaptag = .true.
    960963      ELSE
    961964        ierr = nf90_get_var(nid, nvarid,watercaptagold)
     
    12911294
    12921295c Watercaptag
    1293       if(imold.eq.iim .and. jmold.eq.jjm) then
     1296      if (no_watercaptag .or. imold /= iim .or. jmold /= jjm) then
     1297          if (imold /= iim .or. jmold /= jjm) then
     1298              print*, "We are doing an horizontal interpolation,
     1299     & 'watercaptag' will be initialized with ''locate_watercaptag''"
     1300          endif
     1301          call locate_watercaptag(ngrid,latitude,longitude,watercaptag)
    12941302      else
    1295        print *, "We are doing an horizontal interpolation,
    1296      &    watercaptag will be set to false and redefined proprely
    1297      &    in the PCM(surfini)"
    1298        watercaptagold(:,:)=0
    1299       endif
    1300         call interp_horiz (watercaptagold(:,:),watercaptags(:,:),
    1301      &       imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv)
    1302         call gr_dyn_fi (1,iim+1,jjm+1,ngrid,watercaptags(:,:),
    1303      &       watercaptag_tmp(:))
    1304 
    1305         do i=1,ngrid
    1306           if(watercaptag_tmp(i).gt. 0.5) then
    1307             watercaptag(i)=.true.
    1308           else
    1309             watercaptag(i)=.false.
    1310           endif
    1311         enddo
    1312 
     1303          call interp_horiz (watercaptagold(:,:),watercaptags(:,:),
     1304     &         imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv)
     1305          call gr_dyn_fi (1,iim+1,jjm+1,ngrid,watercaptags(:,:),
     1306     &         watercaptag_tmp(:))
     1307
     1308          do i = 1,ngrid
     1309              if (watercaptag_tmp(i) > 0.5) then
     1310                  watercaptag(i) = .true.
     1311              else
     1312                  watercaptag(i) = .false.
     1313              endif
     1314          enddo
     1315      endif
    13131316
    13141317c Surface albedo
Note: See TracChangeset for help on using the changeset viewer.