source: LMDZ5/branches/LMDZ5V1.0-dev/libf/phylmd/radio_decay.F90 @ 5171

Last change on this file since 5171 was 1409, checked in by jghattas, 14 years ago
  • Added tracer "Age of stratospheric air", activated by putting Aga in tracer.def.
  • The logical rnpb is now controled during execution. rnpb is true only if both tracers RN and PB existe in tracer.def. RN and PB can now be removed from tracer.def
  • In tracer.def, the 2 water traceurs (H2Ov and H2Ol) must still be the first 2 tracers. The following tracers have no specific order(RN and PB can now change places). Still a minimum of 3 tracers are required.
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 1.8 KB
Line 
1!
2! $Id $
3!
4SUBROUTINE radio_decay(radio,rnpb,dtime,tautr,tr,d_tr)
5!
6! Caluclate radioactive decay for all tracers with radio(it)=true
7!
8  USE dimphy
9  USE infotrac, ONLY : nbtr
10  USE traclmdz_mod, ONLY : id_rn, id_pb
11  IMPLICIT NONE
12!-----------------------------------------------------------------------
13! Auteur(s): AA + CG (LGGE/CNRS) Date 24-06-94
14! Objet: Calcul de la tendance radioactive des traceurs type radioelements
15!        Cas particulier pour le couple radon-plomb : Le radon decroit en plomb
16!-----------------------------------------------------------------------
17!
18! Entrees
19!
20  LOGICAL,DIMENSION(nbtr),INTENT(IN)        :: radio ! .true. = traceur radioactif 
21  LOGICAL,INTENT(IN)                        :: rnpb  ! .true. = decroissance RN = source PB
22  REAL,INTENT(IN)                           :: dtime ! Pas de temps physique (secondes)
23  REAL,DIMENSION(nbtr),INTENT(IN)           :: tautr ! Constante de decroissance radioactive
24  REAL,DIMENSION(klon,klev,nbtr),INTENT(IN) :: tr    ! Concentrations traceurs U/kgA
25!
26! Sortie
27!
28  REAL,DIMENSION(klon,klev,nbtr),INTENT(OUT) :: d_tr  ! Tendance de decroissance radioactive
29!
30! Locales
31!
32  INTEGER  :: i,k,it
33
34
35  DO it = 1,nbtr
36     d_tr(:,:,it) = 0.
37     IF ( radio(it) ) THEN
38        IF (tautr(it) .GT. 0.) THEN
39           DO k = 1,klev
40              DO i = 1,klon
41                 d_tr(i,k,it) = - tr(i,k,it) * dtime / tautr(it)
42              END DO
43           END DO
44        END IF
45     END IF
46  END DO
47
48!-------------------------------------------------------
49! Cas particulier radon (id_rn) => plomb (id_pb)
50!-------------------------------------------------------
51  IF ( rnpb ) THEN
52     DO k = 1,klev
53        DO i = 1,klon
54           d_tr(i,k,id_pb) = d_tr(i,k,id_pb) - d_tr(i,k,id_rn)
55        ENDDO
56     ENDDO
57  ENDIF
58
59END SUBROUTINE radio_decay
Note: See TracBrowser for help on using the repository browser.