Ignore:
Timestamp:
Mar 26, 2025, 11:05:58 AM (5 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

File:
1 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
Note: See TracChangeset for help on using the changeset viewer.