source: trunk/LMDZ.COMMON/libf/evolution/frost.F90 @ 4072

Last change on this file since 4072 was 4071, checked in by jbclement, 3 weeks ago

PEM:

  • Making the computation of ice tendencies more reliable by doing it after 'read_startpem' to know exactly where perennial ice is (no matter if there is a "startpem.nc" or not). Moreover, when there is no "startpem.nc", location of perennial ice depends now on 'watercaptag' and on huge amount of frost. This prevents negative ice tendency while there is no ice which can happen with weird PCM inputs (i.e. 'watercaptag' = F & 'watercap' /= 0 & no "stratpem.nc").
  • Few small safeguards throughout the code.

JBC

File size: 5.0 KB
RevLine 
[4065]1MODULE frost
[3991]2!-----------------------------------------------------------------------
3! NAME
[4065]4!     frost
[3991]5!
6! DESCRIPTION
7!     Module for managing frost variables.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
[3984]15
[4065]16! DEPENDENCIES
17! ------------
18use numerics, only: dp, di
19
[3991]20! DECLARATION
21! -----------
[3984]22implicit none
23
[4065]24! VARIABLES
25! ---------
[3989]26! Different types of frost retained by the PEM to give back to the PCM at the end
[4065]27real(dp), dimension(:,:), allocatable :: h2o_frost4PCM
28real(dp), dimension(:,:), allocatable :: co2_frost4PCM
[3984]29
[4065]30! PARAMETERS
31! ----------
32! Different types of frost in the PCM at the beginning [Pa]
33real(dp), dimension(:,:), allocatable, protected :: h2ofrost_PCM
34real(dp), dimension(:,:), allocatable, protected :: co2frost_PCM
[3984]35
36contains
[3991]37!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38
[3984]39!=======================================================================
[4065]40SUBROUTINE ini_frost()
[3991]41!-----------------------------------------------------------------------
42! NAME
[4065]43!     ini_frost
[3991]44!
45! DESCRIPTION
[4065]46!     Initialize the parameters of module 'frost'.
[3991]47!
48! AUTHORS & DATE
49!     JB Clement, 12/2025
50!
51! NOTES
52!
53!-----------------------------------------------------------------------
[3984]54
[4065]55! DEPENDENCIES
56! ------------
57use geometry, only: ngrid, nslope
58
[3991]59! DECLARATION
60! -----------
[3984]61implicit none
62
[3991]63! CODE
64! ----
[4065]65if (.not. allocated(h2ofrost_PCM)) allocate(h2ofrost_PCM(ngrid,nslope))
66if (.not. allocated(co2frost_PCM)) allocate(co2frost_PCM(ngrid,nslope))
67if (.not. allocated(h2o_frost4PCM)) allocate(h2o_frost4PCM(ngrid,nslope))
68if (.not. allocated(co2_frost4PCM)) allocate(co2_frost4PCM(ngrid,nslope))
[3984]69
[4065]70END SUBROUTINE ini_frost
[3984]71!=======================================================================
72
[3991]73!=======================================================================
[4065]74SUBROUTINE end_frost()
[3991]75!-----------------------------------------------------------------------
76! NAME
[4065]77!     end_frost
[3991]78!
79! DESCRIPTION
[4065]80!     Deallocate frost arrays.
[3991]81!
82! AUTHORS & DATE
83!     JB Clement, 12/2025
84!
85! NOTES
[4065]86!
[3991]87!-----------------------------------------------------------------------
[3984]88
[3991]89! DECLARATION
90! -----------
[3984]91implicit none
92
[3991]93! CODE
94! ----
[4065]95if (allocated(h2ofrost_PCM)) deallocate(h2ofrost_PCM)
96if (allocated(co2frost_PCM)) deallocate(co2frost_PCM)
97if (allocated(h2o_frost4PCM)) deallocate(h2o_frost4PCM)
98if (allocated(co2_frost4PCM)) deallocate(co2_frost4PCM)
[3984]99
[4065]100END SUBROUTINE end_frost
[3984]101!=======================================================================
102
[3991]103!=======================================================================
[4065]104SUBROUTINE set_h2ofrost_PCM(h2ofrost_PCM_in)
[3991]105!-----------------------------------------------------------------------
106! NAME
[4065]107!     set_h2ofrost_PCM
[3991]108!
109! DESCRIPTION
[4065]110!     Setter for 'h2ofrost_PCM'.
[3991]111!
112! AUTHORS & DATE
113!     JB Clement, 12/2025
114!
115! NOTES
116!
117!-----------------------------------------------------------------------
[3984]118
[3991]119! DECLARATION
120! -----------
[3984]121implicit none
122
[3991]123! ARGUMENTS
124! ---------
[4065]125real(dp), dimension(:,:), intent(in) :: h2ofrost_PCM_in
[3984]126
[3991]127! CODE
128! ----
[4065]129h2ofrost_PCM(:,:) = h2ofrost_PCM_in(:,:)
[3984]130
[4065]131END SUBROUTINE set_h2ofrost_PCM
[3984]132!=======================================================================
133
[3991]134!=======================================================================
[4065]135SUBROUTINE set_co2frost_PCM(co2frost_PCM_in)
[3991]136!-----------------------------------------------------------------------
137! NAME
[4065]138!     set_co2frost_PCM
[3991]139!
140! DESCRIPTION
[4065]141!     Setter for 'co2frost_PCM'.
[3991]142!
143! AUTHORS & DATE
144!     JB Clement, 12/2025
145!
146! NOTES
147!
148!-----------------------------------------------------------------------
[3984]149
[3991]150! DECLARATION
151! -----------
[3984]152implicit none
153
[3991]154! ARGUMENTS
155! ---------
[4065]156real(dp), dimension(:,:), intent(in) :: co2frost_PCM_in
[3984]157
[3991]158! CODE
159! ----
[4065]160co2frost_PCM(:,:) = co2frost_PCM_in(:,:)
[3984]161
[4065]162END SUBROUTINE set_co2frost_PCM
[3984]163!=======================================================================
164
[3991]165!=======================================================================
[4071]166SUBROUTINE compute_frost4PCM(minPCM_h2ofrost,minPCM_co2frost)
[3991]167!-----------------------------------------------------------------------
168! NAME
[4065]169!     compute_frost4PCM
[3991]170!
171! DESCRIPTION
[4065]172!     Compute the frost to give back to the PCM (metamorphism).
[3991]173!
174! AUTHORS & DATE
175!     JB Clement, 12/2025
176!
177! NOTES
[4065]178!     Frost for the PEM is the extra part of the PCM frost above the
179!     yearly minimum.
[3991]180!-----------------------------------------------------------------------
[3984]181
[4065]182! DEPENDENCIES
183! ------------
184use display, only: print_msg
185
[3991]186! DECLARATION
187! -----------
[3984]188implicit none
189
[4065]190! ARGUMENTS
191! ---------
[4071]192real(dp), dimension(:,:), intent(in) :: minPCM_h2ofrost, minPCM_co2frost
[4065]193
[3991]194! CODE
195! ----
[4065]196call print_msg('> Computing frost to give back to the PCM (metamorphism)')
[3984]197
[4065]198h2o_frost4PCM(:,:) = 0._dp
199co2_frost4PCM(:,:) = 0._dp
[4071]200where (h2ofrost_PCM(:,:) > 0._dp) h2o_frost4PCM(:,:) = h2ofrost_PCM(:,:) - minPCM_h2ofrost(:,:)
201where (co2frost_PCM(:,:) > 0._dp) co2_frost4PCM(:,:) = co2frost_PCM(:,:) - minPCM_co2frost(:,:)
[4065]202
203END SUBROUTINE compute_frost4PCM
[3991]204!=======================================================================
[3984]205
[4065]206END MODULE frost
Note: See TracBrowser for help on using the repository browser.