Changeset 4161


Ignore:
Timestamp:
Apr 1, 2026, 3:17:47 PM (8 days ago)
Author:
jmauxion
Message:

Mars PCM:

  • Fix improvedclouds_mod by adding a missing "else" condition on zdq.
  • Add the possibility to output watercaptag in physiq_mod.
  • Fix write_output to handle logical arrays with an unstructured grid.
  • Update some of the utils (localtime/zrecast) to handle phisfi.

JM

Location:
trunk/LMDZ.MARS
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/changelog.txt

    r4160 r4161  
    51415141== 31/03/2026 == JBC
    51425142Giving standard explicit names to variables related to "paleoclimate" module + adding yearly average of flux exchanged with subsurface ice in XIOS output file intended for the PEM.
     5143
     5144== 01/04/2026 == JM
     5145- Fix improvedclouds_mod by adding a missing "else" condition on zdq.
     5146- Add the possibility to output watercaptag in physiq_mod.
     5147- Fix write_output to handle logical arrays with an unstructured grid.
     5148- Update some of the utils (localtime/zrecast) to handle phisfi.
  • trunk/LMDZ.MARS/libf/phymars/improvedclouds_mod.F90

    r4044 r4161  
    443443        IF (spenttime.ne.0) then
    444444          zdq=(dMicetot/spenttime)!*(ptimestep-spenttime)
     445        ELSE
     446          ! Initialization for spenttime=0
     447          zdq=zpotcond(ig,l)*((ptimestep-spenttime)/ptimestep)
    445448        ENDIF
    446449        zdq=abs(zdq)
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r4151 r4161  
    33373337         call write_output("phisfi","Surface geopotential",
    33383338     &                    "m2s-2",phisfi(:))
     3339         call write_output("watercaptag","Watercap tag"
     3340     &         ,"Boolean",watercaptag(:))
    33393341         if (grid_type == regular_lonlat) then
    33403342           call write_output("area","Mesh area","m2",
  • trunk/LMDZ.MARS/libf/phymars/write_output_mod.F90

    r4058 r4161  
    285285logical, dimension(:), intent(in) :: field
    286286! Local argument used to convert logical to real
    287 real, dimension(ngrid) :: field_real
    288 logical :: is_active ! For XIOS, should this field be sent or not
    289 
     287real, allocatable, dimension(:) :: field_real
     288logical :: is_active ! For XIOS, should this field be sent or not
     289
     290allocate(field_real(size(field,1)))
    290291field_real = 0.
    291292where (field) field_real = 1.
  • trunk/LMDZ.MARS/util/localtime.F90

    r2577 r4161  
    10461046ierr=NF_INQ_VARID(infid,"phisinit",tmpvarid)
    10471047if (ierr.ne.NF_NOERR) then
    1048   write(*,*) "Failed to get phisinit ID. OK"
    1049   phisinit = 0.
    1050   phis = .false.
     1048  write(*,*) "Failed to get phisinit ID..."
     1049  write(*,*) "Try name phisfi (xios name)..."
     1050  ierr=NF_INQ_VARID(infid,"phisfi",tmpvarid)
     1051  if (ierr.ne.NF_NOERR) then
     1052    write(*,*) "Failed again. Ok, I skip phisinit."
     1053    phisinit = 0.
     1054    phis = .false.
     1055  else
     1056    ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
     1057    if (ierr.ne.NF_NOERR) then
     1058      write(*,*) "init2 ERRO: Failed reading phisinit"
     1059      stop
     1060    endif
     1061    phis = .true.
     1062  endif
    10511063else
    10521064  ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
  • trunk/LMDZ.MARS/util/zrecast.F90

    r3057 r4161  
    654654if (ierr.ne.NF_NOERR) then
    655655  write(*,*) "Warning: Failed to get phisinit ID from file ",trim(infile)
    656   infile2="diagfi.nc"
    657   write(*,*) "         Trying file ",trim(infile2)
    658   ierr=NF_OPEN(infile2,NF_NOWRITE,infid2)
     656  write(*,*) "Trying name phisfi (xios name)..."
     657  ierr=NF_INQ_VARID(infid,"phisfi",tmpvarid)
    659658  if (ierr.ne.NF_NOERR) then
    660     write(*,*) "Problem: Could not find/open that file"
    661     infile2="diagfi1.nc"
     659    write(*,*) "Failed to find phisfi..."
     660    infile2="diagfi.nc"
    662661    write(*,*) "         Trying file ",trim(infile2)
    663662    ierr=NF_OPEN(infile2,NF_NOWRITE,infid2)
    664663    if (ierr.ne.NF_NOERR) then
    665664      write(*,*) "Problem: Could not find/open that file"
    666       infile2="phisinit.nc"
     665      infile2="diagfi1.nc"
    667666      write(*,*) "         Trying file ",trim(infile2)
    668667      ierr=NF_OPEN(infile2,NF_NOWRITE,infid2)
    669668      if (ierr.ne.NF_NOERR) then
    670         write(*,*) "Error: Could not open that file either"
    671         write(*,*) "Might as well stop here"
    672         stop
     669        write(*,*) "Problem: Could not find/open that file"
     670        infile2="phisinit.nc"
     671        write(*,*) "         Trying file ",trim(infile2)
     672        ierr=NF_OPEN(infile2,NF_NOWRITE,infid2)
     673        if (ierr.ne.NF_NOERR) then
     674          write(*,*) "Error: Could not open that file either"
     675          write(*,*) "Might as well stop here"
     676          stop
     677        endif
    673678      endif
    674679    endif
    675   endif
    676 
    677   ! Get ID for phisinit
    678   ierr=NF_INQ_VARID(infid2,"phisinit",tmpvarid)
    679   if (ierr.ne.NF_NOERR) then
    680     write(*,*) "Error: Failed to get phisinit ID"
    681     stop
    682   endif
    683   ! Get physinit
    684   ierr=NF_GET_VAR_REAL(infid2,tmpvarid,phisinit)
    685   if (ierr.ne.NF_NOERR) then
    686     write(*,*) "Error: Failed reading phisinit"
    687     stop
    688   endif
    689   ! Close file
    690   write(*,*) 'OK, got phisinit'
    691   ierr=NF_CLOSE(infid2)
     680
     681    ! Get ID for phisinit
     682    ierr=NF_INQ_VARID(infid2,"phisinit",tmpvarid)
     683    if (ierr.ne.NF_NOERR) then
     684      write(*,*) "Error: Failed to get phisinit ID"
     685      stop
     686    endif
     687 
     688    ! Get physinit
     689    ierr=NF_GET_VAR_REAL(infid2,tmpvarid,phisinit)
     690    if (ierr.ne.NF_NOERR) then
     691      write(*,*) "Error: Failed reading phisinit"
     692      stop
     693    endif
     694    ! Close file
     695    write(*,*) 'OK, got phisinit'
     696    ierr=NF_CLOSE(infid2)
     697  else
     698    ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
     699    if (ierr.ne.NF_NOERR) then
     700      write(*,*) "Error: Failed reading phisinit"
     701      stop
     702    endif
     703  endif
    692704else
    693705  ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
     
    697709  endif
    698710endif
    699 
    700711!===============================================================================
    701712! 1.4 Choose and build the new vertical coordinate
Note: See TracChangeset for help on using the changeset viewer.