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

Last change on this file was 3519, checked in by slebonnois, 10 days ago

SL: some cleaning of rcm1d Venus

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