source: trunk/LMDZ.VENUS/libf/phyvenus/dyn1d/rcm1d.F @ 1521

Last change on this file since 1521 was 1517, checked in by slebonnois, 10 years ago

SL: bugs en 1D pour Venus et Titan

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