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

Last change on this file since 1543 was 1543, checked in by emillour, 9 years ago

All models: Further adaptations to keep up with changes in LMDZ5 concerning
physics/dynamics separation:

  • dyn3d:
  • adapted gcm.F so that all physics initializations are now done in iniphysiq.
  • dyn3dpar:
  • adapted gcm.F so that all physics initializations are now done in iniphysiq.
  • updated calfis_p.F to follow up with changes.
  • copied over updated "bands.F90" from LMDZ5.
  • dynphy_lonlat:
  • calfis_p.F90, mod_interface_dyn_phys.F90, follow up of changes in phy_common/mod_* routines
  • phy_common:
  • added "geometry_mod.F90" to store information about the grid (replaces phy*/comgeomphy.F90) and give variables friendlier names: rlond => longitude , rlatd => latitude, airephy => cell_area, cuphy => dx , cvphy => dy
  • added "physics_distribution_mod.F90"
  • updated "mod_grid_phy_lmdz.F90", "mod_phys_lmdz_mpi_data.F90", "mod_phys_lmdz_para.F90", "mod_phys_lmdz_mpi_transfert.F90", "mod_grid_phy_lmdz.F90", "mod_phys_lmdz_omp_data.F90", "mod_phys_lmdz_omp_transfert.F90", "write_field_phy.F90" and "ioipsl_getin_p_mod.F90" to LMDZ5 versions.
  • phy[venus/titan/mars/std]:
  • removed "init_phys_lmdz.F90", "comgeomphy.F90"; adapted routines to use geometry_mod (longitude, latitude, cell_area, etc.)

EM

File size: 14.3 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      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
202c   la surface de chaque maille est inutile en 1D --->
203      area(1)=1.E+0
204c de meme ?
205      cufi(1)=1.E+0
206      cvfi(1)=1.E+0
207
208c Ehouarn: iniphysiq requires arrays related to (3D) dynamics grid,
209c e.g. for cell boundaries, which are meaningless in 1D; so pad these
210c with '0.' when necessary
211      CALL iniphysiq(1,1,llm,
212     &            1,comm_lmdz,
213     &            daysec,day0,dtphys,
214     &            (/lati(1),0./),(/0./),
215     &            (/0.,0./),(/long(1),0./),
216     &            (/ (/area,0./),(/0.,0./) /),
217     &            (/cufi,0.,0.,0./),
218     &            (/cvfi,0./),
219     &            rad,g,r,cpp,1)
220
221      call ini_cpdet
222
223c   le geopotentiel au sol est inutile en 1D car tout est controle
224c   par la pression de surface --->
225      phisfi(1)=0.E+0
226
227c   Initialisation pour prendre en compte les vents en 1-D
228c   ------------------------------------------------------
229 
230c    vent geostrophique
231      PRINT *,'composante vers l est du vent geostrophique (U) ?'
232      READ(unit,*) gru
233      PRINT *,'composante vers le nord du vent geostrophique (V) ?'
234      READ(unit,*) grv
235
236c     Initialisation des vents  au premier pas de temps
237      DO ilayer=1,nlayer
238         u(ilayer)=gru
239         v(ilayer)=grv
240      ENDDO
241
242c  calcul des pressions et altitudes en utilisant les niveaux sigma
243c  ----------------------------------------------------------------
244
245c    Vertical Coordinates  (hybrids)
246c    """"""""""""""""""""
247      CALL  disvert_noterre
248     
249c     Calcul au milieu des couches : Vient de la version Mars
250c     WARNING : le choix de placer le milieu des couches au niveau de
251c     pression intermédiaire est arbitraire et pourrait etre modifié.
252c     C'est fait de la meme facon dans disvert
253
254      DO l = 1, llm
255       aps(l) =  0.5 *( ap(l) +ap(l+1))
256       bps(l) =  0.5 *( bp(l) +bp(l+1))
257      ENDDO
258
259      DO ilevel=1,nlevel
260        plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
261      ENDDO
262
263      DO ilayer=1,nlayer
264        play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
265        pk(ilayer)  =cpp*(play(ilayer)/preff)**rcp
266c       write(120,*) ilayer,plev(ilayer),play(ilayer)
267      ENDDO
268c     write(120,*) nlevel,plev(nlevel)
269c     stop
270     
271      pks=cpp*(psurf/preff)**rcp
272
273c  init des variables pour phyredem
274c  --------------------------------
275      call phys_state_var_init
276
277c  profil de temperature et altitude au premier appel
278c  --------------------------------------------------
279
280c modif par rapport a Mars:
281c   on envoie dz/T=-log(play/psurf)*r/g dans profile
282      tmp1(0)=0.0
283      tmp1(1)= -log(play(1)/psurf)*r/g
284      DO ilayer=2,nlayer
285        tmp1(ilayer)=-log(play(ilayer)/play(ilayer-1))*r/g
286      ENDDO
287      call profile(unit,nlayer+1,tmp1,tmp2)
288      CLOSE(unit)
289
290      print*,"               Pression        Altitude     Temperature"
291      ilayer=1
292      ftsol(1)=tmp2(0)
293       temp(1)=tmp2(1)
294       zlay(1)=tmp2(1)*tmp1(1)
295      print*,"           0",ftsol(1)
296      print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
297      DO ilayer=2,nlayer
298        temp(ilayer)=tmp2(ilayer)
299        zlay(ilayer)=zlay(ilayer-1)+tmp2(ilayer)*tmp1(ilayer)
300        print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
301      ENDDO
302
303      allocate(tmoy(llm))
304      tmoy(:)=temp(:)
305     
306c     temperature du sous-sol
307c     ~~~~~~~~~~~~~~~~~~~~~~~
308      DO isoil=1,nsoil
309         ftsoil(1,isoil)=ftsol(1)
310      ENDDO
311
312c    Initialisation des traceurs
313c    ---------------------------
314
315      DO iq=1,nqtot
316        DO ilayer=1,nlayer
317           q(ilayer,iq) = 0.
318        ENDDO
319      ENDDO
320
321c FULL CHEMISTRY !! AJOUTER INIT AURELIEN...
322C Faudrait lire les cles avant pour mettre ca en option....
323c ou alors mettre ca dans physiq
324
325c    Initialisation des parametres d'oro
326c    -----------------------------------
327
328      zmea(1) = 0.
329      zstd(1) = 0.
330      zsig(1) = 0.
331      zgam(1) = 0.
332      zthe(1) = 0.
333      zpic(1) = 0.
334      zval(1) = 0.
335
336c  Initialisation albedo
337c  ----------------------
338
339      falbe(1)=0.1
340
341c  Ecriture de "startphy.nc"
342c  -------------------------
343c  (Ce fichier sera aussitot relu au premier
344c   appel de "physiq", mais il est necessaire pour passer
345c   les variables purement physiques a "physiq"...
346
347      solsw(1)    = 0.
348      sollw(1)    = 0.
349      fder(1)      = 0.
350      radsol(1)   = 0.
351     
352      radpas      = NINT(1.*day_step/nbapp_rad)
353      soil_model  = .true.
354
355      call phyredem("startphy.nc")
356
357c  deallocation des variables phyredem
358c  -----------------------------------
359      call phys_state_var_end
360
361c=======================================================================
362c  BOUCLE TEMPORELLE DU MODELE 1D
363c=======================================================================
364c
365      firstcall=.true.
366      lastcall=.false.
367
368      DO idt=1,ndt
369        IF (idt.eq.ndt) then
370         lastcall=.true.
371c toujours nulle dans le cas de Venus, pour l'instant...
372         zls = 0.0
373c        write(103,*) 'Ls=',zls*180./pi
374c        write(103,*) 'Lat=', lati(1)
375c        write(103,*) 'RunEnd - Atmos. Temp. File'
376c        write(103,*) 'RunEnd - Atmos. Temp. File'
377c        write(104,*) 'Ls=',zls*180./pi
378c        write(104,*) 'Lat=', lati(1)
379c        write(104,*) 'RunEnd - Atmos. Temp. File'
380        ENDIF
381
382c    calcul du geopotentiel
383c     ~~~~~~~~~~~~~~~~~~~~~
384! ADAPTATION GCM POUR CP(T)
385      DO ilayer=1,nlayer
386        s(ilayer)=(play(ilayer)/psurf)**rcp
387      ENDDO
388      phi(1)=cpp*temp(1)*(1.E+0-s(1))
389      DO ilayer=2,nlayer
390         phi(ilayer)=phi(ilayer-1)+
391     &     cpp*(temp(ilayer-1)/s(ilayer-1)+temp(ilayer)/s(ilayer))*0.5
392     &        *(s(ilayer-1)-s(ilayer))
393
394      ENDDO
395
396c       appel de la physique
397c       --------------------
398
399      CALL physiq (1,llm,nqtot,
400     ,     firstcall,lastcall,
401     ,     day,time,dtphys,
402     ,     plev,play,pk,phi,phisfi,
403     ,     presnivs,
404     ,     u,v,temp,q, 
405     ,     w,
406C - sorties
407     s     du,dv,dtemp,dq,dpsurf)
408
409c     calcul de rho
410       rho = 0.
411c     print*,rho
412
413
414c     print*,"DT APRES PHYSIQ=",day,time,dtime
415c     print*,dtemp
416c     print*,temp
417c     print*," "
418c     stop
419
420c       evolution du vent : modele 1D
421c       -----------------------------
422 
423c       la physique calcule les derivees temporelles de u et v.
424c       Pas de coriolis
425          DO ilayer=1,nlayer
426             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
427             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
428          ENDDO
429c     
430c       Calcul du temps au pas de temps suivant
431c       ---------------------------------------
432        firstcall=.false.
433        time=time+dtphys/daysec
434        IF (time.gt.1.E+0) then
435            time=time-1.E+0
436            day=day+1
437        ENDIF
438
439c       calcul des vitesses et temperature au pas de temps suivant
440c       ----------------------------------------------------------
441
442        DO ilayer=1,nlayer
443           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
444           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
445           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
446        ENDDO
447
448c       calcul des pressions au pas de temps suivant
449c       ----------------------------------------------------------
450
451           psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
452           DO ilevel=1,nlevel
453             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
454           ENDDO
455           DO ilayer=1,nlayer
456             play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
457           ENDDO
458
459      ENDDO   ! fin de la boucle temporelle
460
461c    ========================================================
462c    GESTION DES SORTIE
463c    ========================================================
464
465        print*,"Temperature finale:"
466        print*,temp
467       
468c stabilite
469      DO ilayer=1,nlayer
470        zlay(ilayer) = phi(ilayer)/g/1000.  !en km
471      ENDDO
472      DO ilayer=2,nlayer
473        tmp1(ilayer) =
474     .    (temp(ilayer)-temp(ilayer-1))/(zlay(ilayer)-zlay(ilayer-1))
475     .   + 1000.*g/cpp
476      ENDDO
477
478      OPEN(11,file='profile.new')
479      DO ilayer=1,nlayer
480        write (11,*) zlay(ilayer),temp(ilayer),tmp1(ilayer)
481      ENDDO
482
483c    ========================================================
484      END
485 
486c***********************************************************************
487c***********************************************************************
488
489!#include "../dyn3d_common/disvert_noterre.F"
490!#include "../dyn3d/abort_gcm.F"
491
Note: See TracBrowser for help on using the repository browser.