source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/dyn3d_common/infotrac.F90 @ 222

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

Creating temporary dynamico/lmdz/saturn branche

YM

File size: 15.6 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, ierr2, ierr3
70   
71    character(len=80) :: line ! to store a line of text
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    IF (planet_type=='earth') THEN
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/='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     ! Test if config_inca is other then none for run without INCA
118     IF (type_trac/='inca' .AND. config_inca/='none') THEN
119       WRITE(lunout,*) 'config_inca will now be changed to none as you do not couple with INCA model'
120       config_inca='none'
121     END IF
122    ELSE
123     type_trac='plnt'  ! planets... May want to dissociate between each later.
124    ENDIF ! of IF (planet_type=='earth')
125
126!-----------------------------------------------------------------------
127!
128! 1) Get the true number of tracers + water vapor/liquid
129!    Here true tracers (nqtrue) means declared tracers (only first order)
130!
131!-----------------------------------------------------------------------
132    IF (planet_type=='earth') THEN
133     IF (type_trac == 'lmdz' .OR. type_trac == 'repr') THEN
134       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
135       IF(ierr.EQ.0) THEN
136          WRITE(lunout,*) trim(modname),': Open traceur.def : ok'
137          READ(90,*) nqtrue
138       ELSE
139          WRITE(lunout,*) trim(modname),': Problem in opening traceur.def'
140          WRITE(lunout,*) trim(modname),': WARNING using defaut values'
141          nqtrue=4 ! Defaut value
142       END IF
143       ! For Earth, water vapour & liquid tracers are not in the physics
144       nbtr=nqtrue-2
145     ELSE ! type_trac=inca
146       ! nbtr has been read from INCA by init_cont_lmdz() in gcm.F
147       nqtrue=nbtr+2
148     END IF
149
150     IF (nqtrue < 2) THEN
151       WRITE(lunout,*) trim(modname),': nqtrue=',nqtrue, ' is not allowded. 2 tracers is the minimum'
152       CALL abort_gcm('infotrac_init','Not enough tracers',1)
153     END IF
154
155! Transfert number of tracers to Reprobus
156     IF (type_trac == 'repr') THEN
157#ifdef REPROBUS
158       CALL Init_chem_rep_trac(nbtr)
159#endif
160     END IF
161
162    ELSE  ! not Earth
163       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
164       IF(ierr.EQ.0) THEN
165          WRITE(lunout,*) 'Open traceur.def : ok'
166          READ(90,*) nqtrue
167       ELSE
168          WRITE(lunout,*) 'Problem in opening traceur.def'
169          WRITE(lunout,*) 'ATTENTION using defaut values: nqtrue=1'
170          nqtrue=1 ! Defaut value
171       END IF
172       ! Other planets (for now); we have the same number of tracers
173       ! in the dynamics than in the physics
174       nbtr=nqtrue
175     
176    ENDIF  ! planet_type
177!
178! Allocate variables depending on nqtrue and nbtr
179!
180    ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue))
181    ALLOCATE(conv_flg(nbtr), pbl_flg(nbtr), tracnam(nbtr))
182    conv_flg(:) = 1 ! convection activated for all tracers
183    pbl_flg(:)  = 1 ! boundary layer activated for all tracers
184
185!-----------------------------------------------------------------------
186! 2)     Choix  des schemas d'advection pour l'eau et les traceurs
187!
188!     iadv = 1    schema  transport type "humidite specifique LMD"
189!     iadv = 2    schema   amont
190!     iadv = 14   schema  Van-leer + humidite specifique
191!                            Modif F.Codron
192!     iadv = 10   schema  Van-leer (retenu pour l'eau vapeur et liquide)
193!     iadv = 11   schema  Van-Leer pour hadv et version PPM (Monotone) pour vadv
194!     iadv = 12   schema  Frederic Hourdin I
195!     iadv = 13   schema  Frederic Hourdin II
196!     iadv = 16   schema  PPM Monotone(Collela & Woodward 1984)
197!     iadv = 17   schema  PPM Semi Monotone (overshoots autorisés)
198!     iadv = 18   schema  PPM Positif Defini (overshoots undershoots autorisés)
199!     iadv = 20   schema  Slopes
200!     iadv = 30   schema  Prather
201!
202!        Dans le tableau q(ij,l,iq) : iq = 1  pour l'eau vapeur
203!                                     iq = 2  pour l'eau liquide
204!       Et eventuellement             iq = 3,nqtot pour les autres traceurs
205!
206!        iadv(1): choix pour l'eau vap. et  iadv(2) : choix pour l'eau liq.
207!------------------------------------------------------------------------
208!
209!    Get choice of advection schema from file tracer.def or from INCA
210!---------------------------------------------------------------------
211    IF (planet_type=='earth') THEN
212     IF (type_trac == 'lmdz' .OR. type_trac == 'repr') THEN
213       IF(ierr.EQ.0) THEN
214          ! Continue to read tracer.def
215          DO iq=1,nqtrue
216             READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
217          END DO
218          CLOSE(90) 
219       ELSE ! Without tracer.def, set default values (for Earth!)
220         if ((nqtrue==4).and.(planet_type=="earth")) then
221          hadv(1) = 14
222          vadv(1) = 14
223          tnom_0(1) = 'H2Ov'
224          hadv(2) = 10
225          vadv(2) = 10
226          tnom_0(2) = 'H2Ol'
227          hadv(3) = 10
228          vadv(3) = 10
229          tnom_0(3) = 'RN'
230          hadv(4) = 10
231          vadv(4) = 10
232          tnom_0(4) = 'PB'
233         else
234           ! Error message, we need a traceur.def file
235           write(lunout,*) trim(modname),&
236           ': Cannot set default tracer names!'
237           write(lunout,*) trim(modname),' Make a traceur.def file!!!'
238           CALL abort_gcm('infotrac_init','Need a traceur.def file!',1)
239         endif ! of if (nqtrue==4)
240       END IF
241       
242       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
243       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
244       DO iq=1,nqtrue
245          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq)
246       END DO
247
248     ELSE  ! type_trac=inca : config_inca='aero' ou 'chem'
249! le module de chimie fournit les noms des traceurs
250! et les schemas d'advection associes.
251     
252#ifdef INCA
253       CALL init_transport( &
254            hadv, &
255            vadv, &
256            conv_flg, &
257            pbl_flg,  &
258            tracnam)
259#endif
260       tnom_0(1)='H2Ov'
261       tnom_0(2)='H2Ol'
262
263       DO iq =3,nqtrue
264          tnom_0(iq)=tracnam(iq-2)
265       END DO
266
267     END IF ! type_trac
268
269    ELSE  ! not Earth
270       IF(ierr.EQ.0) THEN
271          ! Continue to read tracer.def
272          DO iq=1,nqtrue
273             !READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
274            ! try to be smart when reading traceur.def
275            read(90,'(80a)') line ! store the line from traceur.def
276            ! assume format is hadv,vadv,tnom_0
277            read(line,*,iostat=ierr2) hadv(iq),vadv(iq),tnom_0(iq)
278            if (ierr2.ne.0) then
279              ! maybe format is tnom0,hadv,vadv
280              read(line,*,iostat=ierr3) tnom_0(iq),hadv(iq),vadv(iq)
281              if (ierr3.ne.0) then
282                ! assume only tnom0 is provided (havd and vad default to 10)
283                read(line,*) tnom_0(iq)
284                hadv(iq)=10
285                vadv(iq)=10
286              endif
287            endif ! of if(ierr2.ne.0)
288          END DO ! of DO iq=1,nqtrue
289          CLOSE(90) 
290       ELSE ! Without tracer.def
291          hadv(1) = 10
292          vadv(1) = 10
293          tnom_0(1) = 'dummy'
294       END IF
295       
296       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
297       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
298       DO iq=1,nqtrue
299          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq)
300       END DO
301
302    ENDIF  ! planet_type
303
304!-----------------------------------------------------------------------
305!
306! 3) Verify if advection schema 20 or 30 choosen
307!    Calculate total number of tracers needed: nqtot
308!    Allocate variables depending on total number of tracers
309!-----------------------------------------------------------------------
310    new_iq=0
311    DO iq=1,nqtrue
312       ! Add tracers for certain advection schema
313       IF (hadv(iq)<20 .AND. vadv(iq)<20 ) THEN
314          new_iq=new_iq+1  ! no tracers added
315       ELSE IF (hadv(iq)==20 .AND. vadv(iq)==20 ) THEN
316          new_iq=new_iq+4  ! 3 tracers added
317       ELSE IF (hadv(iq)==30 .AND. vadv(iq)==30 ) THEN
318          new_iq=new_iq+10 ! 9 tracers added
319       ELSE
320          WRITE(lunout,*) trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
321          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 1',1)
322       END IF
323    END DO
324   
325    IF (new_iq /= nqtrue) THEN
326       ! The choice of advection schema imposes more tracers
327       ! Assigne total number of tracers
328       nqtot = new_iq
329
330       WRITE(lunout,*) trim(modname),': The choice of advection schema for one or more tracers'
331       WRITE(lunout,*) 'makes it necessary to add tracers'
332       WRITE(lunout,*) trim(modname)//': ',nqtrue,' is the number of true tracers'
333       WRITE(lunout,*) trim(modname)//': ',nqtot, ' is the total number of tracers needed'
334
335    ELSE
336       ! The true number of tracers is also the total number
337       nqtot = nqtrue
338    END IF
339
340!
341! Allocate variables with total number of tracers, nqtot
342!
343    ALLOCATE(tname(nqtot), ttext(nqtot))
344    ALLOCATE(iadv(nqtot), niadv(nqtot))
345
346!-----------------------------------------------------------------------
347!
348! 4) Determine iadv, long and short name
349!
350!-----------------------------------------------------------------------
351    new_iq=0
352    DO iq=1,nqtrue
353       new_iq=new_iq+1
354
355       ! Verify choice of advection schema
356       IF (hadv(iq)==vadv(iq)) THEN
357          iadv(new_iq)=hadv(iq)
358       ELSE IF (hadv(iq)==10 .AND. vadv(iq)==16) THEN
359          iadv(new_iq)=11
360       ELSE
361          WRITE(lunout,*)trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
362
363          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 2',1)
364       END IF
365     
366       str1=tnom_0(iq)
367       tname(new_iq)= tnom_0(iq)
368       IF (iadv(new_iq)==0) THEN
369          ttext(new_iq)=trim(str1)
370       ELSE
371          ttext(new_iq)=trim(tnom_0(iq))//descrq(iadv(new_iq))
372       END IF
373
374       ! schemas tenant compte des moments d'ordre superieur
375       str2=ttext(new_iq)
376       IF (iadv(new_iq)==20) THEN
377          DO jq=1,3
378             new_iq=new_iq+1
379             iadv(new_iq)=-20
380             ttext(new_iq)=trim(str2)//txts(jq)
381             tname(new_iq)=trim(str1)//txts(jq)
382          END DO
383       ELSE IF (iadv(new_iq)==30) THEN
384          DO jq=1,9
385             new_iq=new_iq+1
386             iadv(new_iq)=-30
387             ttext(new_iq)=trim(str2)//txtp(jq)
388             tname(new_iq)=trim(str1)//txtp(jq)
389          END DO
390       END IF
391    END DO
392
393!
394! Find vector keeping the correspodence between true and total tracers
395!
396    niadv(:)=0
397    iiq=0
398    DO iq=1,nqtot
399       IF(iadv(iq).GE.0) THEN
400          ! True tracer
401          iiq=iiq+1
402          niadv(iiq)=iq
403       ENDIF
404    END DO
405
406
407    WRITE(lunout,*) trim(modname),': Information stored in infotrac :'
408    WRITE(lunout,*) trim(modname),': iadv  niadv tname  ttext :'
409    DO iq=1,nqtot
410       WRITE(lunout,*) iadv(iq),niadv(iq),&
411       ' ',trim(tname(iq)),' ',trim(ttext(iq))
412    END DO
413
414!
415! Test for advection schema.
416! This version of LMDZ only garantees iadv=10 and iadv=14 (14 only for water vapour) .
417!
418    DO iq=1,nqtot
419       IF (iadv(iq)/=10 .AND. iadv(iq)/=14 .AND. iadv(iq)/=0) THEN
420          WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
421          CALL abort_gcm('infotrac_init','In this version only iadv=10 and iadv=14 is tested!',1)
422       ELSE IF (iadv(iq)==14 .AND. iq/=1) THEN
423          WRITE(lunout,*)trim(modname),'STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
424          CALL abort_gcm('infotrac_init','In this version iadv=14 is only permitted for water vapour!',1)
425       END IF
426    END DO
427
428!-----------------------------------------------------------------------
429! Finalize :
430!
431    DEALLOCATE(tnom_0, hadv, vadv)
432    DEALLOCATE(tracnam)
433
434  END SUBROUTINE infotrac_init
435
436END MODULE infotrac
Note: See TracBrowser for help on using the repository browser.