source: LMDZ5/trunk/libf/dyn3d_common/infotrac.F90 @ 1945

Last change on this file since 1945 was 1945, checked in by lguez, 11 years ago

Duplicated files moved to dyn3d_common.

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