source: trunk/LMDZ.COMMON/libf/evolution/tracers.F90 @ 4076

Last change on this file since 4076 was 4074, checked in by jbclement, 11 days ago

PEM:

  • Correct management of H2O ice tendency in 1D when there is not enough ice anymore.
  • Clean initialization of allocatable module arrays (especially needed when no slope)
  • One more renaming for consistency + few small updates thoughout the code.

JBC

File size: 9.3 KB
Line 
1MODULE tracers
2!-----------------------------------------------------------------------
3! NAME
4!     tracers
5!
6! DESCRIPTION
7!     Tracer species management for tracking.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
15
16! DEPENDENCIES
17! ------------
18use numerics, only: dp, di, k4
19
20! DECLARATION
21! -----------
22implicit none
23
24! PARAMETERS
25! ----------
26character(11),                                parameter :: traceurdef_name = 'traceur.def'
27integer(di),                                  protected :: nq        ! Number of tracers
28character(30), dimension(:),     allocatable, protected :: qnames    ! Names of tracers
29real(dp),      dimension(:,:,:), allocatable, protected :: q_PCM     ! Tracers in the "start.nc" at the beginning
30integer(di),                                  protected :: iPCM_qh2o ! Index for H2O vapor tracer from PCM
31real(dp),      dimension(:),     allocatable, protected :: mmol      ! Molar masses of tracers [g/mol]
32
33contains
34!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
35
36!=======================================================================
37SUBROUTINE ini_tracers()
38!-----------------------------------------------------------------------
39! NAME
40!     ini_tracers
41!
42! DESCRIPTION
43!     Allocate tracer molar mass array.
44!
45! AUTHORS & DATE
46!     JB Clement, 12/2025
47!
48! NOTES
49!
50!-----------------------------------------------------------------------
51
52! DEPENDENCIES
53! ------------
54use stoppage, only: stop_clean
55use geometry, only: ngrid, nlayer
56use display,  only: print_msg
57use utility,  only: int2str
58
59! DECLARATION
60! -----------
61implicit none
62
63! LOCAL VARIABLES
64! ---------------
65logical(k4) :: here
66integer(di) :: i, ierr, funit
67
68! CODE
69! ----
70inquire(file = traceurdef_name,exist = here)
71if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//traceurdef_name//'"!',1)
72call print_msg('> Reading "'//traceurdef_name//'"')
73open(newunit = funit,file = traceurdef_name,status = 'old',form = 'formatted',action = 'read',iostat = ierr)
74if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "'//traceurdef_name//'"!',ierr)
75read(funit,*) nq
76call print_msg('nq = '//int2str(nq))
77
78! Allocation
79if (.not. allocated(mmol)) allocate(mmol(nq))
80if (.not. allocated(qnames)) allocate(qnames(nq))
81if (.not. allocated(q_PCM)) allocate(q_PCM(ngrid,nlayer,nq))
82mmol(:) = 0._dp
83q_PCM(:,:,:) = 0._dp
84
85! Initialize the names of tracers
86do i = 1,nq
87    read(funit,*,iostat = ierr) qnames(i)
88    if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error reading the names of tracers in "'//traceurdef_name//'"!',1)
89end do
90close(funit)
91
92! Initialize the indices of tracers
93iPCM_qh2o = -1
94do i = 1,nq
95    if (qnames(i) == "h2o_vap") then
96        iPCM_qh2o = i
97        mmol(i) = 18._dp
98    end if
99end do
100
101! Checking if everything has been found
102if (iPCM_qh2o < 0) call stop_clean(__FILE__,__LINE__,'H2O vapour index not found!',1)
103
104END SUBROUTINE ini_tracers
105!=======================================================================
106
107!=======================================================================
108SUBROUTINE end_tracers()
109!-----------------------------------------------------------------------
110! NAME
111!     end_tracers
112!
113! DESCRIPTION
114!     Deallocate tracer molar mass array.
115!
116! AUTHORS & DATE
117!     JB Clement, 12/2025
118!
119! NOTES
120!
121!-----------------------------------------------------------------------
122
123! DECLARATION
124! -----------
125implicit none
126
127! CODE
128! ----
129if (allocated(mmol)) deallocate(mmol)
130if (allocated(qnames)) deallocate(qnames)
131if (allocated(q_PCM)) deallocate(q_PCM)
132
133END SUBROUTINE end_tracers
134!=======================================================================
135
136!=======================================================================
137SUBROUTINE set_q_PCM(q_PCM_in,iq)
138!-----------------------------------------------------------------------
139! NAME
140!     set_q_PCM
141!
142! DESCRIPTION
143!     Setter for 'q_PCM'.
144!
145! AUTHORS & DATE
146!     JB Clement, 12/2025
147!
148! NOTES
149!
150!-----------------------------------------------------------------------
151
152! DEPENDENCIES
153! ------------
154use geometry, only: dyngrd2vect, nlon, nlat, ngrid
155
156! DECLARATION
157! -----------
158implicit none
159
160! ARGUMENTS
161! ---------
162real(dp), dimension(:,:,:), intent(in) :: q_PCM_in
163integer(di),                intent(in) :: iq
164
165! LOCAL VARIABLES
166! ---------------
167integer(di) :: l
168
169! CODE
170! ----
171do l = 1,size(q_PCM_in,3)
172    call dyngrd2vect(nlon + 1,nlat,ngrid,q_PCM_in(:,:,l),q_PCM(:,l,iq))
173end do
174
175END SUBROUTINE set_q_PCM
176!=======================================================================
177
178!=======================================================================
179SUBROUTINE adapt_tracers2pressure(ps_avg_glob_old,ps_avg_glob,ps_ts,q_h2o_ts,q_co2_ts)
180!-----------------------------------------------------------------------
181! NAME
182!     adapt_tracers2pressure
183!
184! DESCRIPTION
185!     Adapt the timeseries of tracers VMR to the new pressure.
186!
187! AUTHORS & DATE
188!     JB Clement, 01/2026
189!
190! NOTES
191!
192!-----------------------------------------------------------------------
193
194! DEPENDENCIES
195! ------------
196use geometry,   only: ngrid, nday
197use atmosphere, only: ap, bp
198use display,    only: print_msg
199
200! DECLARATION
201! -----------
202implicit none
203
204! ARGUMENTS
205! ---------
206real(dp),                        intent(in)    :: ps_avg_glob_old, ps_avg_glob
207real(dp), dimension(ngrid,nday), intent(inout) :: ps_ts, q_h2o_ts, q_co2_ts
208
209! LOCAL VARIABLES
210! ---------------
211real(dp), dimension(ngrid,nday) :: plev1, plev2, plev1_old, plev2_old
212
213! CODE
214! ----
215call print_msg("> Adapting the timeseries of tracers VMR to the new pressure")
216
217! Build the first levels of pressure
218plev1_old(:,:) = ap(1) + bp(1)*ps_ts(:,:)
219plev2_old(:,:) = ap(2) + bp(2)*ps_ts(:,:)
220
221! Evolve timeseries of surface pressure
222ps_ts(:,:) = ps_ts(:,:)*ps_avg_glob/ps_avg_glob_old
223
224! Build the first levels of pressure
225plev1(:,:) = ap(1) + bp(1)*ps_ts(:,:)
226plev2(:,:) = ap(2) + bp(2)*ps_ts(:,:)
227
228! Adapt the timeseries of VMR tracers
229! H2O
230q_h2o_ts(:,:) = q_h2o_ts(:,:)*(plev1_old(:,:) - plev2_old(:,:))/(plev1(:,:) - plev2(:,:))
231where (q_h2o_ts < 0._dp) q_h2o_ts = 0._dp
232where (q_h2o_ts > 1._dp) q_h2o_ts = 1._dp
233! CO2
234q_co2_ts(:,:) = q_co2_ts(:,:)*(plev1_old(:,:) - plev2_old(:,:))/(plev1(:,:) - plev2(:,:))
235where (q_co2_ts < 0._dp) q_co2_ts = 0._dp
236where (q_co2_ts > 1._dp) q_co2_ts = 1._dp
237
238END SUBROUTINE adapt_tracers2pressure
239!=======================================================================
240
241!=======================================================================
242SUBROUTINE build4PCM_tracers(ps4PCM,q4PCM)
243!-----------------------------------------------------------------------
244! NAME
245!     build4PCM_tracers
246!
247! DESCRIPTION
248!     Reconstructs the tracer VMRs for the PCM.
249!
250! AUTHORS & DATE
251!     JB Clement, 12/2025
252!     C. Metz, 01/2026
253!
254! NOTES
255!
256!-----------------------------------------------------------------------
257
258! DEPENDENCIES
259! ------------
260use atmosphere, only: ap, bp, ps_PCM
261use geometry,   only: ngrid, nlayer
262use display,    only: print_msg
263
264! DECLARATION
265! -----------
266implicit none
267
268! ARGUMENTS
269! ---------
270real(dp), dimension(:),     intent(in)  :: ps4PCM
271real(dp), dimension(:,:,:), intent(out) :: q4PCM
272
273! LOCAL VARIABLES
274! ---------------
275integer(di)                           :: i, l, iq
276real(dp)                              :: extra_mass
277real(dp), dimension(:,:), allocatable :: dz_plev_PCM, dz_plev_new
278
279! CODE
280! ----
281call print_msg('> Building tracers for the PCM')
282allocate(dz_plev_PCM(ngrid,nlayer),dz_plev_new(ngrid,nlayer))
283
284! Compute thickness of pressure levels (= atmospheric layers)
285do l = 1,nlayer
286    dz_plev_PCM(:,l) = ap(l) + bp(l)*ps_PCM(:) - ap(l + 1) - bp(l + 1)*ps_PCM(:)
287    dz_plev_new(:,l) = ap(l) + bp(l)*ps_PCM(:) - ap(l + 1) - bp(l + 1)*ps4PCM(:)
288end do
289
290! Reconstruct tracers
291do iq = 1,nq
292    do l = 1,nlayer
293        do i = 1,ngrid
294            ! Standard scaling
295            q4PCM(i,l,iq) = q_PCM(i,l,iq)*dz_plev_PCM(i,l)/dz_plev_new(i,l)
296
297            ! CO2 logic: account for the dry air mass change
298            if (trim(qnames(iq)) == "co2") q4PCM(i,l,iq) = q4PCM(i,l,iq) + (dz_plev_new(i,l) - dz_plev_PCM(i,l)) / dz_plev_new(i,l)
299        end do
300    end do
301end do
302
303! Mass conservation and clipping
304do iq = 1,nq
305    do i = 1,ngrid
306        do l = 1,nlayer
307            ! Clipping negative values
308            if (q4PCM(i,l,iq) < 0._dp) q4PCM(i,l,iq) = 0._dp
309
310            ! VMR limit check (ignore number density tracers)
311            if (q4PCM(i,l,iq) > 1._dp .and.              &
312                (qnames(iq) /= "dust_number") .and.      &
313                (qnames(iq) /= "ccn_number") .and.       &
314                (qnames(iq) /= "stormdust_number") .and. &
315                (qnames(iq) /= "topdust_number")) then
316
317                if (l < nlayer) then
318                    extra_mass = (q4PCM(i,l,iq) - 1.)*dz_plev_new(i,l) ! extra_mass units: [VMR*pressure]
319                    q4PCM(i,l,iq) = 1._dp
320
321                    ! Transfer extra mass to the layer above
322                    ! Note: We divide by dz_plev_new(i,l + 1) to convert back to VMR
323                    q4PCM(i,l + 1,iq) = q4PCM(i,l + 1,iq) + extra_mass/dz_plev_new(i,l + 1)
324                else
325                    q4PCM(i,l,iq) = 1._dp
326                end if
327            end if
328        end do
329    end do
330end do
331
332deallocate(dz_plev_PCM,dz_plev_new)
333
334END SUBROUTINE build4PCM_tracers
335!=======================================================================
336
337END MODULE tracers
Note: See TracBrowser for help on using the repository browser.