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

Last change on this file since 1017 was 1017, checked in by emillour, 11 years ago

Common dynamics: (and collateral adaptations in Venus physics)
Improved cpdet routines in and additional sponge mode:

  • Additionnal sponge mode (trigered with "callsponge" flag), in line with the one used in the Generic and Martian GCM. This sponge is called whenever there is a dissipation step.
  • Improvement of the cpdet routines : created routines tpot2t_glo_p and t2tpot_glo_p which handle fields on the whole dynamics (scaler) grid, which are more efficient than calling tpot2t_p or t2tpot_p with slabs of data (generated use of intermediate copies of these chunks of data at every call)
  • Turned cpdet.F into a module cpdet_mod.F90 (and correspondingly adapted all routines in the Venus physics).

EM

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