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

Last change on this file since 3984 was 3984, checked in by jbclement, 7 days ago

PEM:
Addition of a module "tracers" to retain properties of atmospheric tracers.
JBC

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