source: trunk/LMDZ.COMMON/libf/dyn3dpar/infotrac.F90 @ 832

Last change on this file since 832 was 492, checked in by emillour, 13 years ago

Common dynamics: updates to keep up with LMDZ5 Earth (rev 1605)
See file "DOC/chantiers/commit_importants.log" for details.
EM

File size: 14.9 KB
Line 
1! $Id$
2!
3MODULE infotrac
4
5! nqtot : total number of tracers and higher order of moment, water vapor and liquid included
6  INTEGER, SAVE :: nqtot
7
8! nbtr : number of tracers not including higher order of moment or water vapor or liquid
9!        number of tracers used in the physics
10  INTEGER, SAVE :: nbtr
11
12! Name variables
13  CHARACTER(len=20), ALLOCATABLE, DIMENSION(:), SAVE :: tname ! tracer short name for restart and diagnostics
14  CHARACTER(len=23), ALLOCATABLE, DIMENSION(:), SAVE :: ttext ! tracer long name for diagnostics
15
16! iadv  : index of trasport schema for each tracer
17  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: iadv
18
19! niadv : vector keeping the coorspondance between all tracers(nqtot) treated in the
20!         dynamic part of the code and the tracers (nbtr+2) used in the physics part of the code.
21  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: niadv ! equivalent dyn / physique
22
23! conv_flg(it)=0 : convection desactivated for tracer number it
24  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: conv_flg
25! pbl_flg(it)=0  : boundary layer diffusion desactivaded for tracer number it
26  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: pbl_flg
27
28  CHARACTER(len=4),SAVE :: type_trac
29 
30CONTAINS
31
32  SUBROUTINE infotrac_init
33    USE control_mod
34#ifdef REPROBUS
35    USE CHEM_REP, ONLY : Init_chem_rep_trac
36#endif
37    IMPLICIT NONE
38!=======================================================================
39!
40!   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
41!   -------
42!   Modif special traceur F.Forget 05/94
43!   Modif M-A Filiberti 02/02 lecture de traceur.def
44!
45!   Objet:
46!   ------
47!   GCM LMD nouvelle grille
48!
49!=======================================================================
50!   ... modification de l'integration de q ( 26/04/94 ) ....
51!-----------------------------------------------------------------------
52! Declarations
53
54    INCLUDE "dimensions.h"
55    INCLUDE "iniprint.h"
56
57! Local variables
58    INTEGER, ALLOCATABLE, DIMENSION(:) :: hadv  ! index of horizontal trasport schema
59    INTEGER, ALLOCATABLE, DIMENSION(:) :: vadv  ! index of vertical trasport schema
60
61    CHARACTER(len=15), ALLOCATABLE, DIMENSION(:) :: tnom_0  ! tracer short name
62    CHARACTER(len=8), ALLOCATABLE, DIMENSION(:) :: tracnam ! name from INCA
63    CHARACTER(len=3), DIMENSION(30) :: descrq
64    CHARACTER(len=1), DIMENSION(3)  :: txts
65    CHARACTER(len=2), DIMENSION(9)  :: txtp
66    CHARACTER(len=23)               :: str1,str2
67 
68    INTEGER :: nqtrue  ! number of tracers read from tracer.def, without higer order of moment
69    INTEGER :: iq, new_iq, iiq, jq, ierr
70 
71    character(len=*),parameter :: modname="infotrac_init"
72!-----------------------------------------------------------------------
73! Initialization :
74!
75    txts=(/'x','y','z'/)
76    txtp=(/'x ','y ','z ','xx','xy','xz','yy','yz','zz'/)
77
78    descrq(14)='VLH'
79    descrq(10)='VL1'
80    descrq(11)='VLP'
81    descrq(12)='FH1'
82    descrq(13)='FH2'
83    descrq(16)='PPM'
84    descrq(17)='PPS'
85    descrq(18)='PPP'
86    descrq(20)='SLP'
87    descrq(30)='PRA'
88   
89    IF (planet_type=='earth') THEN
90     ! Coherence test between parameter type_trac, config_inca and preprocessing keys
91     IF (type_trac=='inca') THEN
92       WRITE(lunout,*) 'You have choosen to couple with INCA chemestry model : type_trac=', &
93            type_trac,' config_inca=',config_inca
94       IF (config_inca/='aero' .AND. config_inca/='chem') THEN
95          WRITE(lunout,*) 'Incoherence between type_trac and config_inca. Model stops. Modify run.def'
96          CALL abort_gcm('infotrac_init','Incoherence between type_trac and config_inca',1)
97       END IF
98#ifndef INCA
99       WRITE(lunout,*) 'To run this option you must add cpp key INCA and compile with INCA code'
100       CALL abort_gcm('infotrac_init','You must compile with cpp key INCA',1)
101#endif
102     ELSE IF (type_trac=='repr') THEN
103       WRITE(lunout,*) 'You have choosen to couple with REPROBUS chemestry model : type_trac=', type_trac
104#ifndef REPROBUS
105       WRITE(lunout,*) 'To run this option you must add cpp key REPROBUS and compile with REPRPBUS code'
106       CALL abort_gcm('infotrac_init','You must compile with cpp key REPROBUS',1)
107#endif
108     ELSE IF (type_trac == 'lmdz') THEN
109       WRITE(lunout,*) 'Tracers are treated in LMDZ only : type_trac=', type_trac
110     ELSE
111       WRITE(lunout,*) 'type_trac=',type_trac,' not possible. Model stops'
112       CALL abort_gcm('infotrac_init','bad parameter',1)
113     END IF
114
115     ! Test if config_inca is other then none for run without INCA
116     IF (type_trac/='inca' .AND. config_inca/='none') THEN
117       WRITE(lunout,*) 'config_inca will now be changed to none as you do not couple with INCA model'
118       config_inca='none'
119     END IF
120    ELSE
121     type_trac='plnt'  ! planets... May want to dissociate between each later.
122    ENDIF ! of IF (planet_type=='earth')
123
124!-----------------------------------------------------------------------
125!
126! 1) Get the true number of tracers + water vapor/liquid
127!    Here true tracers (nqtrue) means declared tracers (only first order)
128!
129!-----------------------------------------------------------------------
130    IF (planet_type=='earth') THEN
131     IF (type_trac == 'lmdz' .OR. type_trac == 'repr') THEN
132       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
133       IF(ierr.EQ.0) THEN
134          WRITE(lunout,*) trim(modname),': Open traceur.def : ok'
135          READ(90,*) nqtrue
136       ELSE
137          WRITE(lunout,*) trim(modname),': Problem in opening traceur.def'
138          WRITE(lunout,*) trim(modname),': WARNING using defaut values'
139          nqtrue=4 ! Defaut value
140       END IF
141       ! For Earth, water vapour & liquid tracers are not in the physics
142       nbtr=nqtrue-2
143     ELSE ! type_trac=inca
144       ! nbtr has been read from INCA by init_cont_lmdz() in gcm.F
145       nqtrue=nbtr+2
146     END IF
147
148     IF (nqtrue < 2) THEN
149       WRITE(lunout,*) trim(modname),': nqtrue=',nqtrue, ' is not allowded. 2 tracers is the minimum'
150       CALL abort_gcm('infotrac_init','Not enough tracers',1)
151     END IF
152
153! Transfert number of tracers to Reprobus
154     IF (type_trac == 'repr') THEN
155#ifdef REPROBUS
156       CALL Init_chem_rep_trac(nbtr)
157#endif
158     END IF
159
160    ELSE  ! not Earth
161       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
162       IF(ierr.EQ.0) THEN
163          WRITE(lunout,*) 'Open traceur.def : ok'
164          READ(90,*) nqtrue
165       ELSE
166          WRITE(lunout,*) 'Problem in opening traceur.def'
167          WRITE(lunout,*) 'ATTENTION using defaut values: nqtrue=1'
168          nqtrue=1 ! Defaut value
169       END IF
170       ! Other planets (for now); we have the same number of tracers
171       ! in the dynamics than in the physics
172       nbtr=nqtrue
173     
174    ENDIF  ! planet_type
175!
176! Allocate variables depending on nqtrue and nbtr
177!
178    ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue))
179    ALLOCATE(conv_flg(nbtr), pbl_flg(nbtr), tracnam(nbtr))
180    conv_flg(:) = 1 ! convection activated for all tracers
181    pbl_flg(:)  = 1 ! boundary layer activated for all tracers
182
183!-----------------------------------------------------------------------
184! 2)     Choix  des schemas d'advection pour l'eau et les traceurs
185!
186!     iadv = 1    schema  transport type "humidite specifique LMD"
187!     iadv = 2    schema   amont
188!     iadv = 14   schema  Van-leer + humidite specifique
189!                            Modif F.Codron
190!     iadv = 10   schema  Van-leer (retenu pour l'eau vapeur et liquide)
191!     iadv = 11   schema  Van-Leer pour hadv et version PPM (Monotone) pour vadv
192!     iadv = 12   schema  Frederic Hourdin I
193!     iadv = 13   schema  Frederic Hourdin II
194!     iadv = 16   schema  PPM Monotone(Collela & Woodward 1984)
195!     iadv = 17   schema  PPM Semi Monotone (overshoots autorisés)
196!     iadv = 18   schema  PPM Positif Defini (overshoots undershoots autorisés)
197!     iadv = 20   schema  Slopes
198!     iadv = 30   schema  Prather
199!
200!        Dans le tableau q(ij,l,iq) : iq = 1  pour l'eau vapeur
201!                                     iq = 2  pour l'eau liquide
202!       Et eventuellement             iq = 3,nqtot pour les autres traceurs
203!
204!        iadv(1): choix pour l'eau vap. et  iadv(2) : choix pour l'eau liq.
205!------------------------------------------------------------------------
206!
207!    Get choice of advection schema from file tracer.def or from INCA
208!---------------------------------------------------------------------
209    IF (planet_type=='earth') THEN
210     IF (type_trac == 'lmdz' .OR. type_trac == 'repr') THEN
211       IF(ierr.EQ.0) THEN
212          ! Continue to read tracer.def
213          DO iq=1,nqtrue
214             READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
215          END DO
216          CLOSE(90) 
217       ELSE ! Without tracer.def, set default values (for Earth!)
218         if ((nqtrue==4).and.(planet_type=="earth")) then
219          hadv(1) = 14
220          vadv(1) = 14
221          tnom_0(1) = 'H2Ov'
222          hadv(2) = 10
223          vadv(2) = 10
224          tnom_0(2) = 'H2Ol'
225          hadv(3) = 10
226          vadv(3) = 10
227          tnom_0(3) = 'RN'
228          hadv(4) = 10
229          vadv(4) = 10
230          tnom_0(4) = 'PB'
231         else
232           ! Error message, we need a traceur.def file
233           write(lunout,*) trim(modname),&
234           ': Cannot set default tracer names!'
235           write(lunout,*) trim(modname),' Make a traceur.def file!!!'
236           CALL abort_gcm('infotrac_init','Need a traceur.def file!',1)
237         endif ! of if (nqtrue==4)
238       END IF
239       
240       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
241       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
242       DO iq=1,nqtrue
243          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq)
244       END DO
245
246     ELSE  ! type_trac=inca : config_inca='aero' ou 'chem'
247! le module de chimie fournit les noms des traceurs
248! et les schemas d'advection associes.
249     
250#ifdef INCA
251       CALL init_transport( &
252            hadv, &
253            vadv, &
254            conv_flg, &
255            pbl_flg,  &
256            tracnam)
257#endif
258       tnom_0(1)='H2Ov'
259       tnom_0(2)='H2Ol'
260
261       DO iq =3,nqtrue
262          tnom_0(iq)=tracnam(iq-2)
263       END DO
264
265     END IF ! type_trac
266
267    ELSE  ! not Earth
268       IF(ierr.EQ.0) THEN
269          ! Continue to read tracer.def
270          DO iq=1,nqtrue
271             READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
272          END DO
273          CLOSE(90) 
274       ELSE ! Without tracer.def
275          hadv(1) = 10
276          vadv(1) = 10
277          tnom_0(1) = 'dummy'
278       END IF
279       
280       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
281       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
282       DO iq=1,nqtrue
283          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq)
284       END DO
285
286    ENDIF  ! planet_type
287
288!-----------------------------------------------------------------------
289!
290! 3) Verify if advection schema 20 or 30 choosen
291!    Calculate total number of tracers needed: nqtot
292!    Allocate variables depending on total number of tracers
293!-----------------------------------------------------------------------
294    new_iq=0
295    DO iq=1,nqtrue
296       ! Add tracers for certain advection schema
297       IF (hadv(iq)<20 .AND. vadv(iq)<20 ) THEN
298          new_iq=new_iq+1  ! no tracers added
299       ELSE IF (hadv(iq)==20 .AND. vadv(iq)==20 ) THEN
300          new_iq=new_iq+4  ! 3 tracers added
301       ELSE IF (hadv(iq)==30 .AND. vadv(iq)==30 ) THEN
302          new_iq=new_iq+10 ! 9 tracers added
303       ELSE
304          WRITE(lunout,*) trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
305          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 1',1)
306       END IF
307    END DO
308   
309    IF (new_iq /= nqtrue) THEN
310       ! The choice of advection schema imposes more tracers
311       ! Assigne total number of tracers
312       nqtot = new_iq
313
314       WRITE(lunout,*) trim(modname),': The choice of advection schema for one or more tracers'
315       WRITE(lunout,*) 'makes it necessary to add tracers'
316       WRITE(lunout,*) trim(modname)//': ',nqtrue,' is the number of true tracers'
317       WRITE(lunout,*) trim(modname)//': ',nqtot, ' is the total number of tracers needed'
318
319    ELSE
320       ! The true number of tracers is also the total number
321       nqtot = nqtrue
322    END IF
323
324!
325! Allocate variables with total number of tracers, nqtot
326!
327    ALLOCATE(tname(nqtot), ttext(nqtot))
328    ALLOCATE(iadv(nqtot), niadv(nqtot))
329
330!-----------------------------------------------------------------------
331!
332! 4) Determine iadv, long and short name
333!
334!-----------------------------------------------------------------------
335    new_iq=0
336    DO iq=1,nqtrue
337       new_iq=new_iq+1
338
339       ! Verify choice of advection schema
340       IF (hadv(iq)==vadv(iq)) THEN
341          iadv(new_iq)=hadv(iq)
342       ELSE IF (hadv(iq)==10 .AND. vadv(iq)==16) THEN
343          iadv(new_iq)=11
344       ELSE
345          WRITE(lunout,*)trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
346
347          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 2',1)
348       END IF
349     
350       str1=tnom_0(iq)
351       tname(new_iq)= tnom_0(iq)
352       IF (iadv(new_iq)==0) THEN
353          ttext(new_iq)=trim(str1)
354       ELSE
355          ttext(new_iq)=trim(tnom_0(iq))//descrq(iadv(new_iq))
356       END IF
357
358       ! schemas tenant compte des moments d'ordre superieur
359       str2=ttext(new_iq)
360       IF (iadv(new_iq)==20) THEN
361          DO jq=1,3
362             new_iq=new_iq+1
363             iadv(new_iq)=-20
364             ttext(new_iq)=trim(str2)//txts(jq)
365             tname(new_iq)=trim(str1)//txts(jq)
366          END DO
367       ELSE IF (iadv(new_iq)==30) THEN
368          DO jq=1,9
369             new_iq=new_iq+1
370             iadv(new_iq)=-30
371             ttext(new_iq)=trim(str2)//txtp(jq)
372             tname(new_iq)=trim(str1)//txtp(jq)
373          END DO
374       END IF
375    END DO
376
377!
378! Find vector keeping the correspodence between true and total tracers
379!
380    niadv(:)=0
381    iiq=0
382    DO iq=1,nqtot
383       IF(iadv(iq).GE.0) THEN
384          ! True tracer
385          iiq=iiq+1
386          niadv(iiq)=iq
387       ENDIF
388    END DO
389
390
391    WRITE(lunout,*) trim(modname),': Information stored in infotrac :'
392    WRITE(lunout,*) trim(modname),': iadv  niadv tname  ttext :'
393    DO iq=1,nqtot
394       WRITE(lunout,*) iadv(iq),niadv(iq),&
395       ' ',trim(tname(iq)),' ',trim(ttext(iq))
396    END DO
397
398!
399! Test for advection schema.
400! This version of LMDZ only garantees iadv=10 and iadv=14 (14 only for water vapour) .
401!
402    DO iq=1,nqtot
403       IF (iadv(iq)/=10 .AND. iadv(iq)/=14 .AND. iadv(iq)/=0) THEN
404          WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
405          CALL abort_gcm('infotrac_init','In this version only iadv=10 and iadv=14 is tested!',1)
406       ELSE IF (iadv(iq)==14 .AND. iq/=1) THEN
407          WRITE(lunout,*)trim(modname),'STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
408          CALL abort_gcm('infotrac_init','In this version iadv=14 is only permitted for water vapour!',1)
409       END IF
410    END DO
411
412!-----------------------------------------------------------------------
413! Finalize :
414!
415    DEALLOCATE(tnom_0, hadv, vadv)
416    DEALLOCATE(tracnam)
417
418  END SUBROUTINE infotrac_init
419
420END MODULE infotrac
Note: See TracBrowser for help on using the repository browser.