source: trunk/LMDZ.COMMON/libf/evolution/phys_constants.F90

Last change on this file was 3989, checked in by jbclement, 5 days ago

PEM:
Massive structural refactor of the PEM codebase for improved readability, consistency and maintainability. The goal is to modernize, standardize and consolidate the code while removing legacy complexity. In detail, this change:

  • Performs large-scale cleanup removing unused code, obsolete routines, duplicated functionality and deprecated initialization logic;
  • Removes "*_mod" wrappers;
  • Replaces mixed naming conventions with clearer variable names, domain-based file/module names and purpose-based routine names;
  • Adds native reading/writing capabilities to the PEM removing the cumbersome dependency on Mars PCM subroutines;
  • Centralizes soil/ice/orbital/PEM configuration logic into dedicated modules;
  • Simplifies routines and updates calls/interfaces to match the new structure.

This refactor significantly clarifies the codebase and provides a cleaner foundation for forthcoming developments.
JBC

File size: 2.4 KB
Line 
1MODULE phys_constants
2
3implicit none
4
5real, parameter :: pi = 4.*atan(1.) ! PI = 3.14159...
6
7real :: g     ! Gravity (m/s2)
8real :: r     ! Reduced gas constant,r = 8.314511/(mugaz/1000.0)
9real :: mugaz ! Molar mass of the atmosphere (g/mol)
10real :: rad   ! Radius of the planet (m)
11real :: cpp   ! Cp of the atmosphere
12real :: rcp   ! r/cpp
13
14!=======================================================================
15contains
16!=======================================================================
17
18SUBROUTINE read_constants(filename)
19! To read the necessary constants in the variable 'controle' of the file "startfi.nc"
20
21use netcdf
22
23implicit none
24
25! Arguments
26!----------
27character(*), intent(in) :: filename
28
29! Local variables
30!----------------
31real, dimension(:), allocatable :: controle
32integer :: ncid           ! File ID
33integer :: varid_controle ! Variable ID for 'controle'
34integer :: dimid_index    ! Dimension ID for 'index'
35integer :: nindex         ! Size of dimension 'index'
36integer :: ierr           ! Return codes
37
38! Code
39!-----
40! Open the NetCDF file
41ierr = nf90_open(trim(filename),NF90_NOWRITE,ncid)
42if (ierr /= nf90_noerr) then
43    write(*,*) "Error opening file:", trim(nf90_strerror(ierr))
44    error stop
45endif
46
47! Get the dimension size of 'index'
48ierr = nf90_inq_dimid(ncid,"index",dimid_index)
49if (ierr /= nf90_noerr) then
50    write(*,*) "Error getting dimid 'index':", trim(nf90_strerror(ierr))
51    error stop
52endif
53ierr = nf90_inquire_dimension(ncid,dimid_index,len = nindex)
54if (ierr /= nf90_noerr) then
55    write(*,*) "Error getting dimension length:", trim(nf90_strerror(ierr))
56    error stop
57endif
58
59! Get the variable ID for 'controle'
60allocate(controle(nindex))
61ierr = nf90_inq_varid(ncid,"controle",varid_controle)
62if (ierr /= nf90_noerr) then
63    write(*,*) "Error getting variable ID 'controle':", trim(nf90_strerror(ierr))
64    error stop
65endif
66ierr = nf90_get_var(ncid,varid_controle,controle)
67if (ierr /= nf90_noerr) then
68    write(*,*) "Error reading 'controle':", trim(nf90_strerror(ierr))
69    error stop
70endif
71
72! Close the file
73ierr = nf90_close(ncid)
74if (ierr /= nf90_noerr) then
75    write(*,*) "Error closing file:", trim(nf90_strerror(ierr))
76    error stop
77endif
78
79! Initialize the constants with 'controle' data
80rad   = controle(5)
81g     = controle(7)
82mugaz = controle(8)
83rcp   = controle(9)
84r     = 8.314511*1000./mugaz
85cpp   = r/rcp
86deallocate(controle)
87
88END SUBROUTINE read_constants
89
90END MODULE phys_constants
Note: See TracBrowser for help on using the repository browser.