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

Last change on this file since 2187 was 1726, checked in by slebonnois, 7 years ago

SL: adjustments after revision 1723, + some debug for cloud microphysics config and rcm1d

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