source: trunk/LMDZ.GENERIC/libf/phystd/radioactive_tracers.F90 @ 3746

Last change on this file since 3746 was 3724, checked in by mmaurice, 10 months ago

Generic PCM

Add radioactive tracers in order to compare tracers mixing between 3D and 1D.

MM

File size: 2.0 KB
Line 
1
2subroutine radioactive_tracers(ngrid,nlayer,nq,ptimestep,pq,zdqradio)
3
4    USE tracer_h, only: half_life, top_prod, bot_prod, noms
5    implicit none
6   
7    !==================================================================
8    !     
9    !     Purpose
10    !     -------
11    !     Calculates the decay of radioactive tracers
12    !     
13    !     Authors
14    !     -------
15    !     Maxime Maurice (14/04/2025)
16    !     
17    !==================================================================
18   
19    !     Arguments
20    integer,intent(in) :: ngrid                      ! number of atmospheric columns
21    integer,intent(in) :: nlayer                     ! number of atmospheric layers
22    integer,intent(in) :: nq                         ! number of tracers
23    real,intent(in)    :: ptimestep                  ! time interval
24    real,intent(in)    :: pq(ngrid,nlayer,nq)        ! tracers (kg/kg)
25    real,intent(out)   :: zdqradio(ngrid,nlayer,nq)  ! radioactive tracers tendencies
26
27    !     Local variables
28    integer iq
29
30        zdqradio(1:ngrid,1:nlayer,1:nq) = 0
31        do iq=1,nq
32        if (half_life(iq) .ne. 0) then
33            ! SOURCE
34            ! bottom
35            if (bot_prod(iq) .ne. 0) then
36                zdqradio(1:ngrid,1,iq) = zdqradio(1:ngrid,1,iq) + bot_prod(iq) ! <== bottom production
37            end if ! bot_prod != 0
38            ! top
39            if (top_prod(iq) .ne. 0) then
40                zdqradio(1:ngrid,nlayer,iq) = zdqradio(1:ngrid,nlayer,iq) + top_prod(iq) ! <== top production (prod = log(2)/half_life)
41            end if ! bot_prod != 0
42
43            ! DECAY
44            if (half_life(iq)/10 .le. ptimestep) then ! sanity test
45                write(*,*)"timestep is more than half-life / 10 of tracer ", noms(iq)
46                call abort_physic()
47            end if ! dt > 0.1 * half_life
48            zdqradio(1:ngrid,1:nlayer,iq) = zdqradio(1:ngrid,1:nlayer,iq) - LOG(2.)/half_life(iq)*pq(1:ngrid,1:nlayer,iq) ! <== decay
49        end if ! half_life(iq) != 0
50        end do ! iq=1,nq
51
52end subroutine radioactive_tracers
53
Note: See TracBrowser for help on using the repository browser.