source: trunk/LMDZ.PLUTO/libf/phypluto/tracer_h.F90 @ 3945

Last change on this file since 3945 was 3936, checked in by tbertrand, 3 months ago

PLUTO PCM : correcting a bug in hazecloud (wrong lyman alpha fluxes due to mu0 being negative during nighttime) + cleaning routines
TB

File size: 7.9 KB
Line 
1
2       module tracer_h
3       !!------------------------------------------------------------------------------------------------------
4       !! Stores data related to physics tracers.
5       !!
6       !! The module provides additional methods:
7       !!   - indexoftracer : search for the index of a tracer in the global table (tracers_h:noms) by name.
8       !!   - nameoftracer  : get the name of tracer from a given index (of the global table).
9       !!   - dumptracers   : print the names of all tracers indexes given in argument.
10       !!------------------------------------------------------------------------------------------------------
11       implicit none
12
13       integer, save :: nqtot  ! total number of tracers
14!$OMP THREADPRIVATE(nqtot)
15
16       logical :: moderntracdef=.false. ! Standard or modern traceur.def
17!$OMP THREADPRIVATE(moderntracdef)
18
19       character*30, save, allocatable :: noms(:)! name of the tracer
20       real, save, allocatable :: mmol(:)        ! mole mass of tracer (g/mol)
21       real, save, allocatable :: aki(:)         ! to compute coefficient of thermal concduction if photochem
22       real, save, allocatable :: cpi(:)         ! to compute cpnew in concentration.F if photochem
23       real, save, allocatable :: radius(:)      ! dust and ice particle radius (m)
24       real, save, allocatable :: rho_q(:)       ! tracer densities (kg.m-3)
25       real, save, allocatable :: qext(:)        ! Single Scat. Extinction coeff at 0.67 um
26       real, save, allocatable :: qextrhor(:)    ! Intermediate for computing opt. depth from q
27
28       real,save :: varian      ! Characteristic variance of log-normal distribution
29       real,save :: r3n_q       ! used to compute r0 from number and mass mixing ratio
30       real,save :: rho_dust    ! Mars dust density (kg.m-3)
31       real,save :: rho_ice     ! Water ice density (kg.m-3)
32       real,save :: rho_ch4_ice ! ch4 ice density (kg.m-3)
33       real,save :: rho_co_ice  ! co ice density (kg.m-3)
34       real,save :: rho_n2      ! N2 ice density (kg.m-3)
35       real,save :: lw_ch4      ! Latent heat CH4 gas -> solid
36       real,save :: lw_co       ! Latent heat CO gas -> solid
37       real,save :: lw_n2       ! Latent heat N2 gas -> solid
38       integer,save :: nmono
39       real,save :: ref_r0        ! for computing reff=ref_r0*r0 (in log.n. distribution)
40!$OMP THREADPRIVATE(noms,mmol,aki,cpi,radius,rho_q,qext,qextrhor, &
41        !$OMP varian,r3n_q,rho_dust,rho_ice,rho_n2,lw_n2,ref_r0)
42
43       integer, save, allocatable :: is_rad(:)  ! 1 if   ""    ""  in radiative transfer, else 0
44!$OMP THREADPRIVATE(is_rad)
45
46       integer, save, allocatable :: is_recomb(:)      ! 1 if tracer used in recombining scheme, else 0 (if 1, must have is_rad=1)
47       integer, save, allocatable :: is_recomb_qset(:) ! 1 if tracer k-corr assume predefined vmr, else 0 (if 1, must have is_recomb=1)
48       integer, save, allocatable :: is_recomb_qotf(:) ! 1 if tracer recombination is done on-the-fly, else 0 (if 1, must have is_recomb_qset=0)
49!$OMP THREADPRIVATE(is_recomb,is_recomb_qset,is_recomb_qotf)
50       ! Lists of constants for condensable tracers
51
52! tracer indexes: these are initialized in initracer and should be 0 if the
53!                 corresponding tracer does not exist
54
55       ! Pluto chemistry
56       integer,save :: igcm_co_gas
57       integer,save :: igcm_n2
58       integer,save :: igcm_ar
59       integer,save :: igcm_ch4_gas ! methane gas
60!$OMP THREADPRIVATE(igcm_co_gas,igcm_n2,igcm_ar,igcm_ch4_gas)
61       ! Other tracers
62       integer,save :: igcm_ar_n2   ! for simulations using co2 + neutral gaz
63       integer,save :: igcm_ch4_ice ! methane ice
64       integer,save :: igcm_co_ice  ! CO ice
65!$OMP THREADPRIVATE(igcm_ar_n2,igcm_ch4_ice,igcm_co_ice)
66       integer,save :: igcm_prec_haze
67       integer,save :: igcm_haze
68       integer,save :: igcm_haze10
69       integer,save :: igcm_haze30
70       integer,save :: igcm_haze50
71       integer,save :: igcm_haze100
72!$OMP THREADPRIVATE(igcm_prec_haze,igcm_haze,igcm_haze10,igcm_haze30,igcm_haze50,igcm_haze100)
73       integer,save :: igcm_eddy1e6
74       integer,save :: igcm_eddy1e7
75       integer,save :: igcm_eddy5e7
76       integer,save :: igcm_eddy1e8
77       integer,save :: igcm_eddy5e8
78!$OMP THREADPRIVATE(igcm_eddy1e6,igcm_eddy1e7,igcm_eddy5e7,igcm_eddy1e8,igcm_eddy5e8)
79
80       ! Microphysical model
81       integer, save :: nmicro = 0                 !! Number of microphysics tracers.
82       integer, save, allocatable :: micro_indx(:) !! Indexes of all microphysical tracers
83!$OMP THREADPRIVATE(nmicro)
84
85       CONTAINS
86
87       FUNCTION indexoftracer(name, sensitivity) RESULT(idx)
88          !! Get the index of a tracer by name.
89          !!
90          !! The function searches in the global tracer table (tracer_h:noms)
91          !! for the given name and returns the first index matching "name".
92          !!
93          !! If no name in the table matches the given one, -1 is returned !
94          IMPLICIT NONE
95          CHARACTER(len=*), INTENT(in)  :: name         !! Name of the tracer to search.
96          LOGICAL, OPTIONAL, INTENT(in) :: sensitivity  !! Case sensitivity (true by default).
97          INTEGER                       :: idx          !! Index of the first tracer matching name or -1 if not found.
98          LOGICAL                       :: zsens
99          INTEGER                       :: j
100          CHARACTER(len=LEN(name))      :: zname
101          zsens = .true. ; IF(PRESENT(sensitivity)) zsens = sensitivity
102          idx = -1
103          IF (.NOT.ALLOCATED(noms)) RETURN
104          IF (zsens) THEN
105             DO j=1,SIZE(noms)
106                IF (TRIM(noms(j)) == TRIM(name)) THEN
107                   idx = j ; RETURN
108                ENDIF
109             ENDDO
110          ELSE
111             zname = to_lower(name)
112             DO j=1,SIZE(noms)
113                IF (TRIM(to_lower(noms(j))) == TRIM(zname)) THEN
114                   idx = j ; RETURN
115                ENDIF
116             ENDDO
117          ENDIF
118
119          CONTAINS
120
121          FUNCTION to_lower(istr) RESULT(ostr)
122             !! Lower case conversion function.
123             IMPLICIT NONE
124             CHARACTER(len=*), INTENT(in) :: istr
125             CHARACTER(len=LEN(istr))     :: ostr
126             INTEGER                      :: i,ic
127             ostr = istr
128             DO i = 1, LEN_TRIM(istr)
129                ic = ICHAR(istr(i:i))
130                IF (ic >= 65 .AND. ic < 90) ostr(i:i) = char(ic + 32)
131             ENDDO
132          END FUNCTION to_lower
133       END FUNCTION indexoftracer
134
135       FUNCTION nameoftracer(indx) RESULT(name)
136          !! Get the name of a tracer by index.
137          !!
138          !! The function searches in the global tracer table (tracer_h:noms)
139          !! and returns the name of the tracer at given index.
140          !!
141          !! If the index is out of range an empty string is returned.
142          IMPLICIT NONE
143          INTEGER, INTENT(in) :: indx   !! Index of the tracer name to retrieve.
144          CHARACTER(len=30)   :: name   !! Name of the tracer at given index.
145          name = ''
146          IF (.NOT.ALLOCATED(noms)) RETURN
147          IF (indx <= 0 .OR. indx > SIZE(noms)) RETURN
148          name = noms(indx)
149       END FUNCTION nameoftracer
150
151       SUBROUTINE dumptracers(indexes)
152          !! Print the names of the given list of tracers indexes.
153          INTEGER, DIMENSION(:), INTENT(in) :: indexes
154          INTEGER :: i,idx
155          CHARACTER(len=:), ALLOCATABLE :: suffix
156
157          IF (.NOT.ALLOCATED(noms)) THEN
158             WRITE(*,'(a)') "[tracers_h:dump_tracers] warning: 'noms' is not allocated, initracer has not be called yet"
159             RETURN
160          ENDIF
161
162          DO i=1,size(indexes)
163             idx = indexes(i)
164             IF (ANY(micro_indx == idx)) THEN
165                suffix = ' (micro)'
166             ELSE
167                suffix=" ()"
168             ENDIF
169             WRITE(*,'(I5,(a),I6,(a))') i," -> ",idx ," : "//TRIM(noms(idx))//suffix
170          ENDDO
171       END SUBROUTINE dumptracers
172
173       end module tracer_h
Note: See TracBrowser for help on using the repository browser.