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

Last change on this file since 3585 was 3585, checked in by debatzbr, 12 days ago

Connecting microphysics to radiative transfer + miscellaneous cleans

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