source: dynamico_lmdz/aquaplanet/LMDZ5/libf/dyn1d/infotrac.F90 @ 3815

Last change on this file since 3815 was 3815, checked in by ymipsl, 10 years ago

moved directories dyn1d and moved was empty.
YM

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