Ignore:
Timestamp:
Apr 2, 2026, 9:15:00 AM (12 days ago)
Author:
emillour
Message:

Genric PCM:
Some changes following the "thermodynamics update":

  • bug fix on the size of pplev
  • renaming the module to something more explicit (thermo_mod.F90 => pcm_thermodynamics_mod.F90) and the Exner&potential temperature update routine "thermodynamics" => "thermodynamics_update"
  • some optimizations in "thermodynamics_update" : in the 'thermo_uni_ideal' case seting r(),cpp() and rcp() can be done at first call only.

EM

File:
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phygeneric/pcm_thermodynamics_mod.F90

    r4163 r4164  
    1       module thermo_mod
     1      module pcm_thermodynamics_mod
    22
    33      use comcstfi_mod, only: cppd_ref,cppv_ref,rd_ref, rcp_ref, mugaz_ref
     
    1212      contains
    1313     
    14       subroutine thermodynamics(thermo_phy, pplay, pplev, t, ngrid, nlayer, nq, q, iq, zh, zpopsk)
    15         CHARACTER(64), INTENT(IN) :: thermo_phy               ! flag
    16         REAL, INTENT(IN)    :: pplay(ngrid,nlayer)      ! pressure (Pa)
    17         REAL, INTENT(IN)    :: pplev(ngrid,nlayer)      ! pressure (Pa)
     14      subroutine thermodynamics_update(thermo_phy, pplay, pplev, t, ngrid, nlayer, nq, q, iq, zh, zpopsk)
     15        CHARACTER(LEN=*), INTENT(IN) :: thermo_phy               ! flag
     16        REAL, INTENT(IN)    :: pplay(ngrid,nlayer)    ! pressure (Pa) at mid-layer
     17        REAL, INTENT(IN)    :: pplev(ngrid,nlayer+1)  ! pressure (Pa) at layer interfaces
    1818        REAL, INTENT(IN)    :: t(ngrid,nlayer)      ! temperature (K)
    1919        INTEGER, INTENT(IN) :: ngrid, nlayer,nq     ! Number of cells, vertical layers and tracers
     
    2727        logical, save :: firstcall=.true.
    2828!$OMP THREADPRIVATE(firstcall)
     29        character(len=80),parameter :: myname = "thermodynamics_update"
    2930       
    3031        if (firstcall) then
     
    3233          ALLOCATE(cpp(ngrid,nlayer))
    3334          ALLOCATE(rcp(ngrid,nlayer))
     35         
     36          SELECT CASE (TRIM(thermo_phy))
     37            CASE('thermo_uni_ideal')
     38              ! Ideal gas, homogeneous
     39              r(:,:) = rd_ref
     40              cpp(:,:) = cppd_ref
     41              rcp(:,:) = rcp_ref
     42            CASE DEFAULT
     43            write(*,*) 'Bad selector for thermodynamics mod: <', TRIM(thermo_phy), '>'
     44            call abort_physic(trim(myname),'Bad selector for thermodynamics mod!',1)
     45          END SELECT
     46         
    3447          firstcall=.false.
    35         endif
     48        endif ! of if (firstcall)
    3649       
    3750        SELECT CASE (TRIM(thermo_phy))
     
    3952         CASE('thermo_uni_ideal')
    4053            ! Ideal gas
    41             r(:,:) = rd_ref
    42             cpp(:,:) = cppd_ref
    43             rcp(:,:) = rcp_ref
    4454            do l=1,nlayer
    4555              do ig=1,ngrid
     
    5161            write(*,*) 'Bad selector for thermodynamics mod: <', TRIM(thermo_phy), '>'
    5262            write(*,*) 'Option is <thermo_uni_ideal>'
    53             call abort
     63            call abort_physic(trim(myname),'Bad selector for thermodynamics mod!',1)
    5464        END SELECT
    5565       
    56         end subroutine thermodynamics
     66        end subroutine thermodynamics_update
    5767       
    58         end module thermo_mod
     68        end module pcm_thermodynamics_mod
Note: See TracChangeset for help on using the changeset viewer.