source: trunk/LMDZ.COMMON/libf/dyn3d/infotrac.F90 @ 253

Last change on this file since 253 was 66, checked in by emillour, 14 years ago

EM: Mise a niveau par rapport a la version terrestre (LMDZ5V2.0-dev, rev 1487)

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