Ignore:
Timestamp:
Feb 12, 2026, 9:09:12 AM (2 weeks ago)
Author:
jbclement
Message:

PEM:
Major refactor following the previous ones (r3989 and r3991) completing the large structural reorganization and cleanup of the PEM codebase. This revision introduces newly designed modules, standardizes interfaces with explicit ini/end APIs and adds native NetCDF I/O together with explicit PCM/PEM adapters. In detail:

  • Some PEM models were corrected or improved:
    • Frost/perennial ice semantics are clarified via renaming;
    • Soil temperature remapping clarified, notably by removing the rescaling of temperature deviation;
    • Geothermal flux for the PCM is computed based on the PEM state;
  • New explicit PEM/PCM adapters ("set_*"/"build4PCM_*") to decouple PEM internal representation from PCM file layouts and reconstruct consistent fields returned to the PCM;
  • New explicit build/teardown routines that centralize allocation and initialization ordering, reducing accidental use of uninitialized data and making the model lifecycle explicit;
  • Add native read/write helpers for NetCDF that centralize all low-level NetCDF interactions with major improvements (and more simplicity) compared to legacy PEM/PCM I/O (see the modules "io_netcdf" and "output"). They support reading, creation and writing of "diagevol.nc" (renamed from "diagpem.nc") and start/restart files;
  • Provide well-focused modules ("numerics"/"maths"/"utility"/"display") to host commonly-used primitives:
    • "numerics" defines numerical types and constants for reproducibility, portability across compilers and future transitions (e.g. quadruple precision experiments);
    • "display" provides a single controlled interface for runtime messages, status output and diagnostics, avoiding direct 'print'/'write' to enable silent mode, log redirection, and MPI-safe output in the future.
    • "utility" (new module) hosts generic helpers used throughout the code (e.g. "int2str" or "real2str");
  • Add modules "clim_state_init"/"clim_state_rec" which provide robust read/write logic for "start/startfi/startpem", including 1D fallbacks, mesh conversions and dimension checks. NetCDF file creation is centralized and explicit. Restart files are now self-consistent and future-proof, requiring changes only to affected variables;
  • Add module "atmosphere" which computes pressure fields, reconstructs potential temperature and air mass. It also holds the whole logic to define sigma or hybrid coordinates for altitudes;
  • Add module "geometry" to centrilize dimensions logic and grid conversions routines (including 2 new ones "dyngrd2vect"/"vect2dyngrd");
  • Add module "slopes" to isolate slopes handling;
  • Add module "surface" to isolate surface management. Notably, albedo and emissivity are now fully reconstructed following the PCM settings;
  • Add module "allocation" to check the dimension initialization and centrilize allocation/deallocation;
  • Finalize module decomposition and renaming to consolidate domain-based modules, purpose-based routines and physics/process-based variables;
  • The main program now drives a clearer sequence of conceptual steps (initialization / reading / evolution / update / build / writing) and fails explicitly instead of silently defaulting;
  • Ice table logic is made restart-safe;
  • 'Open'/'read' intrinsic logic is made safe and automatic;
  • Improve discoverability and standardize the data handling (private vs protected vs public);
  • Apply consistent documentation/header style already introduced;
  • Update deftank scripts to reflect new names and launch wrappers;

This revision is a structural milestone aiming to be behavior-preserving where possible. It has been tested via compilation and short integration runs. However, due to extensive renames, moves, and API changes, full validation is still ongoing.
Note: the revision includes one (possibly two) easter egg hidden in the code for future archaeologists of the PEM. No physics were harmed.
JBC

File:
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/evolution/physics.F90

    r4064 r4065  
    1 MODULE phys_constants
     1MODULE physics
    22!-----------------------------------------------------------------------
    33! NAME
    4 !     phys_constants
     4!     physics
    55!
    66! DESCRIPTION
    7 !     Physical constants used across PEM modules.
     7!     Fundamental constants read from file and used across PEM modules.
    88!
    99! AUTHORS & DATE
     
    1111!
    1212! NOTES
    13 !     Constants are initialized from NetCDF 'startfi.nc' via read_constants.
     13!
    1414!-----------------------------------------------------------------------
     15
     16! DEPENDENCIES
     17! ------------
     18use numerics, only: dp
    1519
    1620! DECLARATION
     
    1822implicit none
    1923
    20 ! MODULE PARAMETERS
    21 ! -----------------
    22 real, parameter :: pi = 4.*atan(1.) ! PI = 3.14159...
    23 
    24 real :: g     ! Gravity (m/s2)
    25 real :: r     ! Reduced gas constant,r = 8.314511/(mugaz/1000.0)
    26 real :: mugaz ! Molar mass of the atmosphere (g/mol)
    27 real :: rad   ! Radius of the planet (m)
    28 real :: cpp   ! Cp of the atmosphere
    29 real :: rcp   ! r/cpp
     24! PARAMETERS
     25! ----------
     26real(dp), protected :: g     ! Gravity [m/s2]
     27real(dp), protected :: r     ! Reduced gas constant,r = 8.314511/[mugaz/1000.0]
     28real(dp), protected :: mugaz ! Molar mass of the atmosphere [g/mol]
     29real(dp), protected :: rad   ! Radius of the planet [m]
     30real(dp), protected :: cpp   ! Cp of the atmosphere
     31real(dp), protected :: rcp   ! r/cpp (= kappa in the dynamics)
    3032
    3133contains
     
    3335
    3436!=======================================================================
    35 SUBROUTINE read_constants(filename)
     37SUBROUTINE init_physics(rad_in,g_in,mugaz_in,rcp_in)
    3638!-----------------------------------------------------------------------
    3739! NAME
    38 !     read_constants
     40!     init_physics
    3941!
    4042! DESCRIPTION
    41 !     Read physical constants from NetCDF file 'startfi.nc' in variable
    42 !     'controle' and initialize module-level constants.
     43!     Initialize physical constants.
    4344!
    4445! AUTHORS & DATE
    45 !     JB Clement, 12/2025
     46!     JB Clement, 01/2026
    4647!
    4748! NOTES
    48 !     Reads controle(5,7,8,9) for rad, g, mugaz, rcp and computes r.
     49!
    4950!-----------------------------------------------------------------------
    50 
    51 ! DEPENDENCIES
    52 ! ------------
    53 use netcdf
    5451
    5552! DECLARATION
     
    5956! ARGUMENTS
    6057! ---------
    61 character(*), intent(in) :: filename
    62 
    63 ! LOCAL VARIABLES
    64 ! ---------------
    65 real, dimension(:), allocatable :: controle
    66 integer :: ncid           ! File ID
    67 integer :: varid_controle ! Variable ID for 'controle'
    68 integer :: dimid_index    ! Dimension ID for 'index'
    69 integer :: nindex         ! Size of dimension 'index'
    70 integer :: ierr           ! Return codes
     58real(dp), intent(in) :: rad_in, g_in, mugaz_in, rcp_in
    7159
    7260! CODE
    7361! ----
    74 ! Open the NetCDF file
    75 ierr = nf90_open(trim(filename),NF90_NOWRITE,ncid)
    76 if (ierr /= nf90_noerr) then
    77     write(*,*) "Error opening file:", trim(nf90_strerror(ierr))
    78     error stop
    79 endif
     62rad = rad_in
     63g = g_in
     64mugaz = mugaz_in
     65rcp = rcp_in
     66r = 8.314511_dp*1000._dp/mugaz
     67cpp = r/rcp
    8068
    81 ! Get the dimension size of 'index'
    82 ierr = nf90_inq_dimid(ncid,"index",dimid_index)
    83 if (ierr /= nf90_noerr) then
    84     write(*,*) "Error getting dimid 'index':", trim(nf90_strerror(ierr))
    85     error stop
    86 endif
    87 ierr = nf90_inquire_dimension(ncid,dimid_index,len = nindex)
    88 if (ierr /= nf90_noerr) then
    89     write(*,*) "Error getting dimension length:", trim(nf90_strerror(ierr))
    90     error stop
    91 endif
    92 
    93 ! Get the variable ID for 'controle'
    94 allocate(controle(nindex))
    95 ierr = nf90_inq_varid(ncid,"controle",varid_controle)
    96 if (ierr /= nf90_noerr) then
    97     write(*,*) "Error getting variable ID 'controle':", trim(nf90_strerror(ierr))
    98     error stop
    99 endif
    100 ierr = nf90_get_var(ncid,varid_controle,controle)
    101 if (ierr /= nf90_noerr) then
    102     write(*,*) "Error reading 'controle':", trim(nf90_strerror(ierr))
    103     error stop
    104 endif
    105 
    106 ! Close the file
    107 ierr = nf90_close(ncid)
    108 if (ierr /= nf90_noerr) then
    109     write(*,*) "Error closing file:", trim(nf90_strerror(ierr))
    110     error stop
    111 endif
    112 
    113 ! Initialize the constants with 'controle' data
    114 rad   = controle(5)
    115 g     = controle(7)
    116 mugaz = controle(8)
    117 rcp   = controle(9)
    118 r     = 8.314511*1000./mugaz
    119 cpp   = r/rcp
    120 deallocate(controle)
    121 
    122 END SUBROUTINE read_constants
     69END SUBROUTINE init_physics
    12370!=======================================================================
    12471
    125 END MODULE phys_constants
     72END MODULE physics
Note: See TracChangeset for help on using the changeset viewer.