source: trunk/LMDZ.COMMON/libf/evolution/metamorphism.F90 @ 3989

Last change on this file since 3989 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: 3.2 KB
Line 
1MODULE metamorphism
2
3implicit none
4
5! Different types of frost retained by the PEM to give back to the PCM at the end
6real, dimension(:,:), allocatable :: h2o_frost4PCM
7real, dimension(:,:), allocatable :: co2_frost4PCM
8
9! Indices for frost taken from the PCM
10integer :: iPCM_h2ofrost, iPCM_co2frost
11
12!=======================================================================
13contains
14!=======================================================================
15
16SUBROUTINE ini_frost_id(nqtot,noms)
17
18implicit none
19
20! Arguments
21!----------
22integer,                        intent(in) :: nqtot
23character(*), dimension(nqtot), intent(in) :: noms
24
25! Local variables
26!----------------
27integer :: i
28
29! Code
30!-----
31! Initialization
32iPCM_h2ofrost = -1
33iPCM_co2frost = -1
34
35! Getting the index
36do i = 1,nqtot
37    if (trim(noms(i)) == "h2o_ice") iPCM_h2ofrost = i
38    if (trim(noms(i)) == "co2")     iPCM_co2frost = i
39enddo
40
41! Checking if everything has been found
42if (iPCM_h2ofrost < 0) error stop 'ini_frost_id: H2O frost index not found!'
43if (iPCM_co2frost < 0) error stop 'ini_frost_id: CO2 frost index not found!'
44
45END SUBROUTINE ini_frost_id
46!=======================================================================
47
48SUBROUTINE compute_frost(ngrid,nslope,h2ofrost_PCM,min_h2ofrost,co2frost_PCM,min_co2frost)
49
50implicit none
51
52! Arguments
53!----------
54integer,                       intent(in) :: ngrid, nslope
55real, dimension(ngrid,nslope), intent(in) :: h2ofrost_PCM, min_h2ofrost, co2frost_PCM, min_co2frost
56
57! Local variables
58!----------------
59
60! Code
61!-----
62write(*,*) '> Computing frost to give back to the PCM'
63
64! Allocation
65call ini_frost(ngrid,nslope)
66
67! Initialization
68h2o_frost4PCM(:,:) = 0.
69co2_frost4PCM(:,:) = 0.
70
71! Computation: frost for the PEM is the extra part of the PCM frost above the yearly minimum
72where (h2ofrost_PCM(:,:) > 0.) h2o_frost4PCM(:,:) = h2ofrost_PCM(:,:) - min_h2ofrost(:,:)
73where (co2frost_PCM(:,:) > 0.) co2_frost4PCM(:,:) = co2frost_PCM(:,:) - min_co2frost(:,:)
74
75END SUBROUTINE compute_frost
76!=======================================================================
77
78SUBROUTINE set_frost4PCM(PCMfrost)
79
80implicit none
81
82! Arguments
83!----------
84real, dimension(:,:,:), intent(inout) :: PCMfrost
85
86! Local variables
87!----------------
88
89! Code
90!-----
91write(*,*) '> Reconstructing frost for the PCM'
92PCMfrost(:,iPCM_h2ofrost,:) = h2o_frost4PCM(:,:)
93PCMfrost(:,iPCM_co2frost,:) = co2_frost4PCM(:,:)
94
95! Deallocation
96call end_frost()
97
98END SUBROUTINE set_frost4PCM
99!=======================================================================
100
101SUBROUTINE ini_frost(ngrid,nslope)
102
103implicit none
104
105! Arguments
106!----------
107integer, intent(in) :: ngrid, nslope
108
109! Local variables
110!----------------
111
112! Code
113!-----
114if (.not. allocated(h2o_frost4PCM)) allocate(h2o_frost4PCM(ngrid,nslope))
115if (.not. allocated(co2_frost4PCM)) allocate(co2_frost4PCM(ngrid,nslope))
116
117END SUBROUTINE ini_frost
118!=======================================================================
119
120SUBROUTINE end_frost()
121
122implicit none
123
124! Arguments
125!----------
126
127! Local variables
128!----------------
129
130! Code
131!-----
132if (allocated(h2o_frost4PCM)) deallocate(h2o_frost4PCM)
133if (allocated(co2_frost4PCM)) deallocate(co2_frost4PCM)
134
135END SUBROUTINE end_frost
136
137END MODULE metamorphism
Note: See TracBrowser for help on using the repository browser.