source: trunk/LMDZ.TITAN/libf/phytitan/rcm1d.F @ 1243

Last change on this file since 1243 was 1121, checked in by slebonnois, 11 years ago

SL: bugs in rcm1d (venus, titan)

File size: 13.4 KB
Line 
1      PROGRAM rcm1d
2     
3      USE infotrac
4      use control_mod
5      use comgeomphy
6      USE phys_state_var_mod
7      use cpdet_mod, only: ini_cpdet
8      use moyzon_mod, only: plevmoy
9      IMPLICIT NONE
10
11c=======================================================================
12c   subject:
13c   --------
14c   PROGRAM useful to run physical part of the venusian GCM in a 1D column
15c       
16c Can be compiled with a command like (e.g. for 55 layers)
17c  "makelmdz -p titan -d 55 rcm1d"
18
19c It requires the files "rcm1d.def" "physiq.def" "traceur.def"
20c      and a file describing the sigma layers (e.g. "z2sig.def")
21c
22c   author: Frederic Hourdin, R.Fournier,F.Forget (original Mars version)
23c   ------- Sebastien Lebonnois (Venus version)
24c   
25c=======================================================================
26
27#include "dimensions.h"
28#include "dimsoil.h"
29#include "comcstfi.h"
30#include "comvert.h"
31#include "netcdf.inc"
32#include "logic.h"
33#include "clesphys.h"
34#include "iniprint.h"
35#include "tabcontrol.h"
36
37c --------------------------------------------------------------
38c  Declarations
39c --------------------------------------------------------------
40c
41      INTEGER unit           ! unite de lecture de "rcm1d.def"
42      INTEGER unitstart      ! unite d'ecriture de "startphy.nc"
43      INTEGER nlayer,nlevel,nsoil,ndt
44      INTEGER ilayer,ilevel,isoil,idt,iq
45      LOGICAl firstcall,lastcall
46c
47      INTEGER day0          ! date initial (sol ; =0 a Ls=0)
48      REAL day              ! date durant le run
49      REAL time             ! time (0<time<1 ; time=0.5 a midi)
50      REAL play(llm)   ! Pressure at the middle of the layers (Pa)
51      REAL plev(llm+1) ! intermediate pressure levels (pa)
52      REAL psurf     
53      REAL u(llm),v(llm)  ! zonal, meridional wind
54      REAL gru,grv   ! prescribed "geostrophic" background wind
55      REAL temp(llm)   ! temperature at the middle of the layers
56      REAL,allocatable :: q(:,:) ! tracer mixing ratio (e.g. kg/kg)
57      REAL zlay(llm)   ! altitude estimee dans les couches (km)
58      REAL long(1),lati(1),area(1)
59      REAL cufi(1),cvfi(1)
60      REAL phisfi(1)
61
62c    Physical and dynamical tandencies (e.g.  m.s-2, K/s, Pa/s)
63      REAL du(llm),dv(llm),dtemp(llm)
64      REAL dudyn(llm),dvdyn(llm),dtempdyn(llm)
65      REAL dpsurf   
66      REAL,allocatable :: dq(:,:)
67
68c   Various intermediate variables
69      REAL zls
70      REAL phi(llm),s(llm)
71      REAL pk(llm),pks, w(llm)
72      INTEGER l, ierr, aslun
73      REAL tmp1(0:llm),tmp2(0:llm),tmp3(0:llm)                       
74
75      character*2 str2
76
77c normalement dans dyn3d/comconst.h
78      COMMON/cpdetvenus/cppdyn,nu_venus,t0_venus
79      REAL cppdyn,nu_venus,t0_venus
80
81c=======================================================================
82
83c=======================================================================
84c INITIALISATION
85c=======================================================================
86
87      lunout = 6
88
89c ------------------------------------------------------
90c  Constantes prescrites ICI
91c ------------------------------------------------------
92
93      pi=2.E+0*asin(1.E+0)
94
95c     Constante de Titan
96c     -------------------
97      planet_type = "titan"
98      rad=2575000.               ! rayon de Venus (m)
99      daysec=1.37889e6           ! duree du sol (s) 
100      omeg=2*pi/daysec           ! vitesse de rotation (rad.s-1)
101      g= 1.35                    ! gravite (m.s-2)
102      mugaz=28.                  ! Masse molaire de l'atm (g.mol-1)
103      cpp=1.039e3
104      cppdyn=cpp
105      r= 8.314511E+0 *1000.E+0/mugaz
106      rcp= r/cpp
107
108c-----------------------------------------------------------------------
109c   Initialisation des traceurs
110c   ---------------------------
111c  Choix du nombre de traceurs et du schema pour l'advection
112c  dans fichier traceur.def
113      call infotrac_init
114
115c Allocation de la tableau q : champs advectes   
116      allocate(q(llm,nqtot))
117      allocate(dq(llm,nqtot))
118
119c ------------------------------------------------------
120c  Lecture des parametres dans "rcm1d.def"
121c ------------------------------------------------------
122
123c   Opening parameters file "rcm1d.def"
124c   ---------------------------------------
125      unit =97
126      OPEN(unit,file='rcm1d.def',status='old',form='formatted'
127     .     ,iostat=ierr)
128
129      IF(ierr.ne.0) THEN
130        write(*,*) 'Problem to open "rcm1d.def'
131        write(*,*) 'Is it there ?'
132        stop
133      END IF
134
135c  Date et heure locale du debut du run
136c  ------------------------------------
137c    Date (en sols depuis le solstice de printemps) du debut du run
138      day0 = 0
139      PRINT *,'date de depart ?'
140      READ(unit,*) day0
141      day=REAL(day0)
142      PRINT *,day0
143c  Heure de demarrage
144      PRINT *,'heure de debut de simulation (entre 0 et 24) ?'
145      READ(unit,*) time
146      time=time/24.E+0
147
148c  Discretisation (Definition de la grille et des pas de temps)
149c  --------------
150c
151      nlayer=llm
152      nlevel=nlayer+1
153      nsoil=nsoilmx
154      PRINT *,'nombre de pas de temps par jour ?'
155      READ(unit,*) day_step
156      print*,day_step
157
158c     PRINT *,'nombre d appel au rayonnement par jour ?'
159c     READ(unit,*) nbapp_rad
160c     print*,nbapp_rad
161c LU DANS PHYSIQ.DEF...
162      nbapp_rad = 100.
163
164      PRINT *,'nombre de jours simules ?'
165      READ(unit,*) ndt
166      print*,ndt
167
168      ndt=ndt*day_step     
169      dtphys=daysec/day_step 
170      dtime=dtphys
171
172c Pression de surface sur la planete
173c ------------------------------------
174c
175      PRINT *,'pression au sol'
176      READ(unit,*) psurf
177      PRINT *,psurf
178c Pression de reference 
179      pa     =  1.e4
180      preff  = 1.4e5
181c     preff  = psurf
182 
183c  latitude/longitude
184c  -------------------
185      PRINT *,'latitude en degres ?'
186      READ(unit,*) lati(1)
187      PRINT *,lati(1)
188      long(1)=0.E+0
189
190c   Initialisation speciales "physiq"
191c   ---------------------------------
192
193      CALL init_phys_lmdz(iim,jjm,llm,1,(/1/))
194      call initcomgeomphy
195      call ini_cpdet
196
197c   la surface de chaque maille est inutile en 1D --->
198      area(1)=1.E+0
199c de meme ?
200      cufi(1)=1.E+0
201      cvfi(1)=1.E+0
202
203c   le geopotentiel au sol est inutile en 1D car tout est controle
204c   par la pression de surface --->
205      phisfi(1)=0.E+0
206
207      CALL iniphysiq(1,llm,daysec,day0,dtphys,
208     .            lati,long,area,cufi,cvfi,rad,g,r,cpp)
209
210c   Initialisation pour prendre en compte les vents en 1-D
211c   ------------------------------------------------------
212 
213c    vent geostrophique
214      PRINT *,'composante vers l est du vent geostrophique (U) ?'
215      READ(unit,*) gru
216      PRINT *,'composante vers le nord du vent geostrophique (V) ?'
217      READ(unit,*) grv
218
219c     Initialisation des vents  au premier pas de temps
220      DO ilayer=1,nlayer
221         u(ilayer)=gru
222         v(ilayer)=grv
223      ENDDO
224
225c  calcul des pressions et altitudes en utilisant les niveaux sigma
226c  ----------------------------------------------------------------
227
228c    Vertical Coordinates  (hybrids)
229c    """"""""""""""""""""
230      CALL  disvert_noterre
231     
232c     Calcul au milieu des couches : Vient de la version Mars
233c     WARNING : le choix de placer le milieu des couches au niveau de
234c     pression intermédiaire est arbitraire et pourrait etre modifié.
235c     C'est fait de la meme facon dans disvert
236
237      DO l = 1, llm
238       aps(l) =  0.5 *( ap(l) +ap(l+1))
239       bps(l) =  0.5 *( bp(l) +bp(l+1))
240      ENDDO
241
242      DO ilevel=1,nlevel
243        plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
244      ENDDO
245      allocate(plevmoy(nlevel))
246      plevmoy(:)=plev(:)
247
248      DO ilayer=1,nlayer
249        play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
250        pk(ilayer)  =cpp*(play(ilayer)/preff)**rcp
251c       write(120,*) ilayer,plev(ilayer),play(ilayer)
252      ENDDO
253c     write(120,*) nlevel,plev(nlevel)
254c     stop
255     
256      pks=cpp*(psurf/preff)**rcp
257
258c  init des variables pour phyredem
259c  --------------------------------
260      call phys_state_var_init
261
262c  profil de temperature et altitude au premier appel
263c  --------------------------------------------------
264
265c modif par rapport a Mars:
266c   on envoie dz/T=-log(play/psurf)*r/g dans profile
267      tmp1(0)=0.0
268      tmp1(1)= -log(play(1)/psurf)*r/g
269      DO ilayer=2,nlayer
270        tmp1(ilayer)=-log(play(ilayer)/play(ilayer-1))*r/g
271      ENDDO
272      DO ilayer=0,nlayer
273        tmp2(ilayer)=plev(ilayer+1)
274      ENDDO
275      call profile(unit,nlayer+1,tmp1,tmp2,tmp3)
276      CLOSE(unit)
277
278      print*,"               Pression        Altitude     Temperature"
279      ilayer=1
280      ftsol(1)=tmp3(0)
281       temp(1)=tmp3(1)
282       zlay(1)=tmp3(1)*tmp1(1)
283      print*,"           0",ftsol(1)
284      print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
285      DO ilayer=2,nlayer
286        temp(ilayer)=tmp3(ilayer)
287        zlay(ilayer)=zlay(ilayer-1)+tmp3(ilayer)*tmp1(ilayer)
288        print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
289      ENDDO
290     
291c     temperature du sous-sol
292c     ~~~~~~~~~~~~~~~~~~~~~~~
293      DO isoil=1,nsoil
294         ftsoil(1,isoil)=ftsol(1)
295      ENDDO
296
297c    Initialisation des traceurs
298c    ---------------------------
299
300      DO iq=1,nqtot
301        DO ilayer=1,nlayer
302           q(ilayer,iq) = 0.
303        ENDDO
304      ENDDO
305
306c    Initialisation des parametres d'oro
307c    -----------------------------------
308
309      zmea(1) = 0.
310      zstd(1) = 0.
311      zsig(1) = 0.
312      zgam(1) = 0.
313      zthe(1) = 0.
314      zpic(1) = 0.
315      zval(1) = 0.
316
317c  Initialisation Ls
318c ------------------
319         zls=0.
320         print*,'Ls=',zls*180./pi
321
322c  Initialisation albedo
323c  ----------------------
324
325      falbe(1)=0.3
326
327c  Ecriture de "startphy.nc"
328c  -------------------------
329c  (Ce fichier sera aussitot relu au premier
330c   appel de "physiq", mais il est necessaire pour passer
331c   les variables purement physiques a "physiq"...
332
333      solsw(1)    = 0.
334      sollw(1)    = 0.
335      fder(1)     = 0.
336      radsol(1)   = 0.
337     
338      radpas      = NINT(1.*day_step/nbapp_rad)
339      soil_model  = .true.
340
341      call phyredem("startphy.nc")
342
343c  deallocation des variables phyredem
344c  -----------------------------------
345      call phys_state_var_end
346
347c=======================================================================
348c  BOUCLE TEMPORELLE DU MODELE 1D
349c=======================================================================
350c
351      firstcall=.true.
352      lastcall=.false.
353
354      DO idt=1,ndt
355        IF (idt.eq.ndt) then
356         lastcall=.true.
357c        write(103,*) 'Ls=',zls*180./pi
358c        write(103,*) 'Lat=', lati(1)
359c        write(103,*) 'RunEnd - Atmos. Temp. File'
360c        write(103,*) 'RunEnd - Atmos. Temp. File'
361c        write(104,*) 'Ls=',zls*180./pi
362c        write(104,*) 'Lat=', lati(1)
363c        write(104,*) 'RunEnd - Atmos. Temp. File'
364        ENDIF
365
366c    calcul du geopotentiel
367c     ~~~~~~~~~~~~~~~~~~~~~
368! ADAPTATION GCM POUR CP(T)
369      DO ilayer=1,nlayer
370        s(ilayer)=(play(ilayer)/psurf)**rcp
371      ENDDO
372      phi(1)=cpp*temp(1)*(1.E+0-s(1))
373      DO ilayer=2,nlayer
374         phi(ilayer)=phi(ilayer-1)+
375     &     cpp*(temp(ilayer-1)/s(ilayer-1)+temp(ilayer)/s(ilayer))*0.5
376     &        *(s(ilayer-1)-s(ilayer))
377
378      ENDDO
379
380c       appel de la physique
381c       --------------------
382
383      CALL physiq (1,llm,nqtot,
384     ,     firstcall,lastcall,
385     ,     day,time,dtphys,
386     ,     plev,play,pk,phi,phisfi,
387     ,     presnivs,
388     ,     u,v,temp,q, 
389     ,     w,
390C - sorties
391     s     du,dv,dtemp,dq,dpsurf)
392
393c     print*,"DT APRES PHYSIQ=",day,time
394c     print*,dtemp
395c     print*,temp
396c     print*," "
397c     stop
398
399c       evolution du vent : modele 1D
400c       -----------------------------
401 
402c       la physique calcule les derivees temporelles de u et v.
403c       Pas de coriolis
404          DO ilayer=1,nlayer
405             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
406             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
407          ENDDO
408c     
409c       Calcul du temps au pas de temps suivant
410c       ---------------------------------------
411        firstcall=.false.
412        time=time+dtphys/daysec
413        IF (time.gt.1.E+0) then
414            time=time-1.E+0
415            day=day+1
416        ENDIF
417
418c       calcul des vitesses et temperature au pas de temps suivant
419c       ----------------------------------------------------------
420
421        DO ilayer=1,nlayer
422           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
423           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
424           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
425        ENDDO
426
427c       calcul des pressions au pas de temps suivant
428c       ----------------------------------------------------------
429
430           psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
431           DO ilevel=1,nlevel
432             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
433           ENDDO
434           DO ilayer=1,nlayer
435             play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
436           ENDDO
437
438      ENDDO   ! fin de la boucle temporelle
439
440c    ========================================================
441c    GESTION DES SORTIE
442c    ========================================================
443
444        print*,"Temperature finale:"
445        print*,temp
446       
447c stabilite
448      DO ilayer=1,nlayer
449        zlay(ilayer) = phi(ilayer)/g/1000.  !en km
450      ENDDO
451      DO ilayer=2,nlayer
452        tmp1(ilayer) =
453     .    (temp(ilayer)-temp(ilayer-1))/(zlay(ilayer)-zlay(ilayer-1))
454     .   + 1000.*g/cpp
455      ENDDO
456
457      OPEN(11,file='profile.new')
458      DO ilayer=1,nlayer
459        write (11,*) zlay(ilayer),temp(ilayer),tmp1(ilayer)
460      ENDDO
461
462c    ========================================================
463      END
464 
465c***********************************************************************
466c***********************************************************************
467
468#include "../dyn3d/disvert_noterre.F"
469#include "../dyn3d/abort_gcm.F"
Note: See TracBrowser for help on using the repository browser.