Changeset 3698 for trunk/LMDZ.PLUTO/libf


Ignore:
Timestamp:
Mar 26, 2025, 11:05:58 AM (10 months ago)
Author:
emillour
Message:

Pluto PCM:
Some fixes to enable runnnig with dynamico:

  • add "strictboundcorrk" flag in callcorrk_pluto to enable running even if outside of kmatrix temperatures (when strictboundcorrk=.true.)
  • add premature exiting of writediagsoil if not with lon-lat grid
  • while at it, turned surfprop.F90 into a module
  • in physiq, enforce the possibility to output subsurface-related field in most cases, not only when "fast=.true."
  • adapt reference xml files: subsurface quantities need to be defined on a dedicated grid, otherwise XIOS will generate misleading garbage values. Updated files are put in "deftank/dynamico" for now.

EM

Location:
trunk/LMDZ.PLUTO/libf/phypluto
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.PLUTO/libf/phypluto/callcorrk_pluto_mod.F90

    r3685 r3698  
    2828                              strobel,vmrch4_proffix,specOLR,vmrch4fix,&
    2929                              haze_radproffix,callmufi
     30      use callkeys_mod, only: strictboundcorrk
    3031      use optcv_pluto_mod, only: optcv_pluto
    3132      use optci_pluto_mod, only: optci_pluto
     
    210211      real,save :: levdat(Nfine),vmrdat(Nfine)
    211212      real :: vmrch4(ngrid,nlayer)              ! vmr ch4 from vmrch4_proffix
     213      character(len=100) :: message
     214      character(len=15),parameter :: subname="callcorrk_pluto"
    212215
    213216!=======================================================================
     
    634637          if(tlevrad(k).lt.tgasmin)then
    635638            print*,'Minimum temperature is outside the radiative'
    636             print*,'transfer kmatrix bounds, exiting.'
     639            print*,'transfer kmatrix bounds'
    637640            print*,'t(',k,')=',tlevrad(k),' < ',tgasmin
    638             ! call abort
     641            if (strictboundcorrk) then
     642              message="Minimum temperature outside of kmatrix bounds"
     643              call abort_physic(subname,message,1)
     644            else
     645              print*,'***********************************************'
     646              print*,'we allow model to continue with tlevrad<tgasmin'
     647              print*,'  ... we assume we know what you are doing ... '
     648              print*,'  ... but do not let this happen too often ... '
     649              print*,'***********************************************'
     650            endif
    639651          elseif(tlevrad(k).gt.tgasmax)then
    640652            print*,'Maximum temperature is outside the radiative'
    641             print*,'transfer kmatrix bounds, exiting.'
     653            print*,'transfer kmatrix bounds'
    642654            print*,'t(',k,')=',tlevrad(k),' > ',tgasmax
    643             ! call abort
     655            if (strictboundcorrk) then
     656              message="Maximum temperature outside of kmatrix bounds"
     657              call abort_physic(subname,message,1)
     658            else
     659              print*,'***********************************************'
     660              print*,'we allow model to continue with tlevrad>tgasmax'
     661              print*,'  ... we assume we know what you are doing ... '
     662              print*,'  ... but do not let this happen too often ... '
     663              print*,'***********************************************'
     664            endif
    644665          endif
    645666         enddo
  • trunk/LMDZ.PLUTO/libf/phypluto/physiq_mod.F90

    r3695 r3698  
    7070      use callcorrk_mod, only: callcorrk
    7171      use callcorrk_pluto_mod, only: callcorrk_pluto
     72      use surfprop_mod, only: surfprop
    7273      use vdifc_mod, only: vdifc
    7374      use vdifc_pluto_mod, only: vdifc_pluto
     
    19581959 !     Info about Ls, declin...
    19591960      IF (fast) THEN
    1960          if (is_master) write(*,*),'Ls=',zls*180./pi,' dec=',declin*180./pi
    1961          if (is_master) write(*,*),'zday=',zday,' ps=',globave
     1961         if (is_master) write(*,*)'Ls=',zls*180./pi,' dec=',declin*180./pi
     1962         if (is_master) write(*,*)'zday=',zday,' ps=',globave
    19621963         IF (lastcall) then
    1963             if (is_master) write(*,*),'lastcall'
     1964            if (is_master) write(*,*)'lastcall'
    19641965         ENDIF
    19651966      ELSE
    1966          if (is_master) write(*,*),'Ls=',zls*180./pi,'decli=',declin*180./pi,'zday=',zday
     1967         if (is_master) write(*,*)'Ls=',zls*180./pi,'decli=',declin*180./pi,'zday=',zday
    19671968      ENDIF
    19681969
     
    21392140      if (fast) then
    21402141         call write_output("fluxrad","fluxrad","W m-2",fluxrad)
    2141          call write_output("fluxgrd","fluxgrd","W m-2",fluxgrd)
    21422142         ! call write_output("dplanck","dplanck","W.s m-2 K-1",dplanck)
     2143      endif
     2144
     2145      if (callsoil) then
    21432146         ! "soil" variables
    2144          call write_output("capcal","capcal","W.s m-2 K-1",capcal)
     2147         call write_output("capcal","Surface Heat Capacity","W.s m-2 K-1",capcal)
    21452148         call write_output("tsoil","tsoil","K",tsoil)
    21462149         call write_output("therm_inertia","therm_inertia","S.I.",therm_inertia)
  • trunk/LMDZ.PLUTO/libf/phypluto/surfprop.F90

    r3539 r3698  
     1module surfprop_mod
     2
     3implicit none
     4
     5contains
     6
    17      subroutine surfprop(ngrid,nq,pfra,qsurf,tsurface,       &
    28                     capcal,adjust,dist,fluold,ptimestep,zls, &
     
    507513      end subroutine surfprop
    508514
     515end module surfprop_mod
  • trunk/LMDZ.PLUTO/libf/phypluto/writediagsoil.F90

    r3506 r3698  
    2323use mod_phys_lmdz_para, only : is_mpi_root, is_master, gather
    2424use mod_grid_phy_lmdz, only : klon_glo, Grid1Dto2D_glo, &
    25                               nbp_lon, nbp_lat
     25                              nbp_lon, nbp_lat, &
     26                              grid_type, unstructured
    2627
    2728implicit none
     
    7879real px2(ngrid)
    7980#endif
     81
     82! 0. This routine shoul only be used in lon-lat case
     83if (grid_type==unstructured) then
     84  return
     85endif
    8086
    8187! 1. Initialization step
Note: See TracChangeset for help on using the changeset viewer.