source: trunk/LMDZ.TITAN.old/libf/phytitan/dyn1d/rcm1d.F @ 3094

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

All GCMs:
Further adaptations to keep up with changes in LMDZ5 concerning
physics/dynamics separation (up to rev r2420 of LMDZ5)

  • all physics packages:
  • added module callphysiq_mod.F90 in dynphy_lonlat/phy* which contains the routine "call_physiq" which is called by calfis* and calls the physics. This way different "physiq" routine from different physics packages may be called: The calfis* routines now exposes all available fields that might be transmitted to physiq but which is actually send (ie: expected/needed by physiq) is decided in call_physiq.
  • turned "physiq.F[90]" into module "physiq_mod.F[90]" for better control of "physiq" arguments. for phyvenus/phytitan, extracted gr_fi_ecrit from physiq.F as gr_fi_ecrit.F90 (note that it can only work in serial).
  • misc:
  • updated wxios.F90 to keep up with LMDZ5 modifications.
  • dyn3d_common:
  • infotrac.F90 keep up with LMDZ5 modifications (cosmetics)
  • dyn3d:
  • gcm.F90: cosmetic cleanup.
  • leapfrog.F90: fix computation of date as function of itau.
  • dyn3dpar:
  • gcm.F: cosmetic cleanup.
  • leapfrog_p.F90: fix computation of date as function of itau.

NB: physics are given the date corresponding to the end of the
physics step.

  • dynphy_lonlat:
  • calfis.F : added computation of relative wind vorticity.
  • calfis_p.F: added computation of relative wind vorticity (input required by Earth physics)

EM

File size: 14.0 KB
Line 
1      PROGRAM rcm1d
2     
3      USE infotrac
4      use control_mod, only: planet_type, day_step
5!      use comgeomphy
6      USE phys_state_var_mod
7      USE comconst_mod, ONLY: cpp,t0_venus,nu_venus
8      use cpdet_mod, only: ini_cpdet
9      use moyzon_mod, only: plevmoy
10      USE comvert_mod, ONLY: ap,bp,presnivs,pa,preff,nivsigs,nivsig,
11     .                       aps,bps,scaleheight,pseudoalt,
12     .                       disvert_type,pressure_exner
13      USE iniphysiq_mod, ONLY: iniphysiq
14      USE mod_const_mpi, ONLY: comm_lmdz
15      USE physiq_mod, ONLY: physiq
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 55 layers)
24c  "makelmdz -p titan -d 55 rcm1d"
25
26c It requires the files "rcm1d.def" "physiq.def" "traceur.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(1)   
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),tmp3(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 Titan
99c     -------------------
100      planet_type = "titan"
101      rad=2575000.               ! rayon de Venus (m)
102      daysec=1.37889e6           ! duree du sol (s) 
103      omeg=2*pi/daysec           ! vitesse de rotation (rad.s-1)
104      g= 1.35                    ! gravite (m.s-2)
105      mugaz=28.                  ! Masse molaire de l'atm (g.mol-1)
106      cpp=1.039e3
107      r= 8.314511E+0 *1000.E+0/mugaz
108      rcp= r/cpp
109
110c-----------------------------------------------------------------------
111c   Initialisation des traceurs
112c   ---------------------------
113c  Choix du nombre de traceurs et du schema pour l'advection
114c  dans fichier traceur.def
115      call infotrac_init
116
117c Allocation de la tableau q : champs advectes   
118      allocate(q(llm,nqtot))
119      allocate(dq(llm,nqtot))
120
121c ------------------------------------------------------
122c  Lecture des parametres dans "rcm1d.def"
123c ------------------------------------------------------
124
125c   Opening parameters file "rcm1d.def"
126c   ---------------------------------------
127      unit =97
128      OPEN(unit,file='rcm1d.def',status='old',form='formatted'
129     .     ,iostat=ierr)
130
131      IF(ierr.ne.0) THEN
132        write(*,*) 'Problem to open "rcm1d.def'
133        write(*,*) 'Is it there ?'
134        stop
135      END IF
136
137c  Date et heure locale du debut du run
138c  ------------------------------------
139c    Date (en sols depuis le solstice de printemps) du debut du run
140      day0 = 0
141      PRINT *,'date de depart ?'
142      READ(unit,*) day0
143      day=REAL(day0)
144      PRINT *,day0
145c  Heure de demarrage
146      PRINT *,'heure de debut de simulation (entre 0 et 24) ?'
147      READ(unit,*) time
148      time=time/24.E+0
149
150c  Discretisation (Definition de la grille et des pas de temps)
151c  --------------
152c
153      nlayer=llm
154      nlevel=nlayer+1
155      nsoil=nsoilmx
156      PRINT *,'nombre de pas de temps par jour ?'
157      READ(unit,*) day_step
158      print*,day_step
159
160c     PRINT *,'nombre d appel au rayonnement par jour ?'
161c     READ(unit,*) nbapp_rad
162c     print*,nbapp_rad
163c LU DANS PHYSIQ.DEF...
164      nbapp_rad = 100.
165
166      PRINT *,'nombre de jours simules ?'
167      READ(unit,*) ndt
168      print*,ndt
169
170      ndt=ndt*day_step     
171      dtphys=daysec/day_step 
172      dtime=dtphys
173
174c Pression de surface sur la planete
175c ------------------------------------
176c
177      PRINT *,'pression au sol'
178      READ(unit,*) psurf
179      PRINT *,psurf
180c Pression de reference 
181      pa     =  1.e4
182      preff  = 1.4e5
183c     preff  = psurf
184 
185c  latitude/longitude
186c  -------------------
187      PRINT *,'latitude en degres ?'
188      READ(unit,*) lati(1)
189      PRINT *,lati(1)
190      long(1)=0.E+0
191
192c   Initialisation speciales "physiq"
193c   ---------------------------------
194
195!      CALL init_phys_lmdz(iim,jjm,llm,1,(/1/))
196
197c   la surface de chaque maille est inutile en 1D --->
198      area(1)=1.E+0
199c de meme ?
200      cufi(1)=1.E+0
201      cvfi(1)=1.E+0
202
203c Ehouarn: iniphysiq requires arrays related to (3D) dynamics grid,
204c e.g. for cell boundaries, which are meaningless in 1D; so pad these
205c with '0.' when necessary
206      CALL iniphysiq(1,1,llm,
207     &            1,comm_lmdz,
208     &            daysec,day0,dtphys,
209     &            (/lati(1),0./),(/0./),
210     &            (/0.,0./),(/long(1),0./),
211     &            (/ (/area,0./),(/0.,0./) /),
212     &            (/cufi,0.,0.,0./),
213     &            (/cvfi,0./),
214     &            rad,g,r,cpp,1)
215
216      call ini_cpdet
217
218c   le geopotentiel au sol est inutile en 1D car tout est controle
219c   par la pression de surface --->
220      phisfi(1)=0.E+0
221
222c   Initialisation pour prendre en compte les vents en 1-D
223c   ------------------------------------------------------
224 
225c    vent geostrophique
226      PRINT *,'composante vers l est du vent geostrophique (U) ?'
227      READ(unit,*) gru
228      PRINT *,'composante vers le nord du vent geostrophique (V) ?'
229      READ(unit,*) grv
230
231c     Initialisation des vents  au premier pas de temps
232      DO ilayer=1,nlayer
233         u(ilayer)=gru
234         v(ilayer)=grv
235      ENDDO
236
237c  calcul des pressions et altitudes en utilisant les niveaux sigma
238c  ----------------------------------------------------------------
239
240c    Vertical Coordinates  (hybrids)
241c    """"""""""""""""""""
242      CALL  disvert_noterre
243     
244c     Calcul au milieu des couches : Vient de la version Mars
245c     WARNING : le choix de placer le milieu des couches au niveau de
246c     pression intermédiaire est arbitraire et pourrait etre modifié.
247c     C'est fait de la meme facon dans disvert
248
249      DO l = 1, llm
250       aps(l) =  0.5 *( ap(l) +ap(l+1))
251       bps(l) =  0.5 *( bp(l) +bp(l+1))
252      ENDDO
253
254      DO ilevel=1,nlevel
255        plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
256      ENDDO
257      allocate(plevmoy(nlevel))
258      plevmoy(:)=plev(:)
259
260      DO ilayer=1,nlayer
261        play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
262        pk(ilayer)  =cpp*(play(ilayer)/preff)**rcp
263c       write(120,*) ilayer,plev(ilayer),play(ilayer)
264      ENDDO
265c     write(120,*) nlevel,plev(nlevel)
266c     stop
267     
268      pks=cpp*(psurf/preff)**rcp
269
270c  init des variables pour phyredem
271c  --------------------------------
272      call phys_state_var_init
273
274c  profil de temperature et altitude au premier appel
275c  --------------------------------------------------
276
277c modif par rapport a Mars:
278c   on envoie dz/T=-log(play/psurf)*r/g dans profile
279      tmp1(0)=0.0
280      tmp1(1)= -log(play(1)/psurf)*r/g
281      DO ilayer=2,nlayer
282        tmp1(ilayer)=-log(play(ilayer)/play(ilayer-1))*r/g
283      ENDDO
284      DO ilayer=0,nlayer
285        tmp2(ilayer)=plev(ilayer+1)
286      ENDDO
287      call profile(unit,nlayer+1,tmp1,tmp2,tmp3)
288      CLOSE(unit)
289
290      print*,"               Pression        Altitude     Temperature"
291      ilayer=1
292      ftsol(1)=tmp3(0)
293       temp(1)=tmp3(1)
294       zlay(1)=tmp3(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)=tmp3(ilayer)
299        zlay(ilayer)=zlay(ilayer-1)+tmp3(ilayer)*tmp1(ilayer)
300        print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
301      ENDDO
302     
303c     temperature du sous-sol
304c     ~~~~~~~~~~~~~~~~~~~~~~~
305      DO isoil=1,nsoil
306         ftsoil(1,isoil)=ftsol(1)
307      ENDDO
308
309c    Initialisation des traceurs
310c    ---------------------------
311
312      DO iq=1,nqtot
313        DO ilayer=1,nlayer
314           q(ilayer,iq) = 0.
315        ENDDO
316      ENDDO
317
318c    Initialisation des parametres d'oro
319c    -----------------------------------
320
321      zmea(1) = 0.
322      zstd(1) = 0.
323      zsig(1) = 0.
324      zgam(1) = 0.
325      zthe(1) = 0.
326      zpic(1) = 0.
327      zval(1) = 0.
328
329c  Initialisation Ls
330c ------------------
331         zls=0.
332         print*,'Ls=',zls*180./pi
333
334c  Initialisation albedo
335c  ----------------------
336
337      falbe(1)=0.3
338
339c  Ecriture de "startphy.nc"
340c  -------------------------
341c  (Ce fichier sera aussitot relu au premier
342c   appel de "physiq", mais il est necessaire pour passer
343c   les variables purement physiques a "physiq"...
344
345      solsw(1)    = 0.
346      sollw(1)    = 0.
347      fder(1)     = 0.
348      radsol(1)   = 0.
349     
350      radpas      = NINT(1.*day_step/nbapp_rad)
351      soil_model  = .true.
352
353      call phyredem("startphy.nc")
354
355c  deallocation des variables phyredem
356c  -----------------------------------
357      call phys_state_var_end
358
359c=======================================================================
360c  BOUCLE TEMPORELLE DU MODELE 1D
361c=======================================================================
362c
363      firstcall=.true.
364      lastcall=.false.
365
366      DO idt=1,ndt
367        IF (idt.eq.ndt) then
368         lastcall=.true.
369c        write(103,*) 'Ls=',zls*180./pi
370c        write(103,*) 'Lat=', lati(1)
371c        write(103,*) 'RunEnd - Atmos. Temp. File'
372c        write(103,*) 'RunEnd - Atmos. Temp. File'
373c        write(104,*) 'Ls=',zls*180./pi
374c        write(104,*) 'Lat=', lati(1)
375c        write(104,*) 'RunEnd - Atmos. Temp. File'
376        ENDIF
377
378c    calcul du geopotentiel
379c     ~~~~~~~~~~~~~~~~~~~~~
380! ADAPTATION GCM POUR CP(T)
381      DO ilayer=1,nlayer
382        s(ilayer)=(play(ilayer)/psurf)**rcp
383      ENDDO
384      phi(1)=cpp*temp(1)*(1.E+0-s(1))
385      DO ilayer=2,nlayer
386         phi(ilayer)=phi(ilayer-1)+
387     &     cpp*(temp(ilayer-1)/s(ilayer-1)+temp(ilayer)/s(ilayer))*0.5
388     &        *(s(ilayer-1)-s(ilayer))
389
390      ENDDO
391
392c       appel de la physique
393c       --------------------
394
395      CALL physiq (1,llm,nqtot,
396     ,     firstcall,lastcall,
397     ,     day,time,dtphys,
398     ,     plev,play,pk,phi,phisfi,
399     ,     presnivs,
400     ,     u,v,temp,q, 
401     ,     w,
402C - sorties
403     s     du,dv,dtemp,dq,dpsurf)
404
405c     print*,"DT APRES PHYSIQ=",day,time
406c     print*,dtemp
407c     print*,temp
408c     print*," "
409c     stop
410
411c       evolution du vent : modele 1D
412c       -----------------------------
413 
414c       la physique calcule les derivees temporelles de u et v.
415c       Pas de coriolis
416          DO ilayer=1,nlayer
417             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
418             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
419          ENDDO
420c     
421c       Calcul du temps au pas de temps suivant
422c       ---------------------------------------
423        firstcall=.false.
424        time=time+dtphys/daysec
425        IF (time.gt.1.E+0) then
426            time=time-1.E+0
427            day=day+1
428        ENDIF
429
430c       calcul des vitesses et temperature au pas de temps suivant
431c       ----------------------------------------------------------
432
433        DO ilayer=1,nlayer
434           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
435           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
436           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
437        ENDDO
438
439c       calcul des pressions au pas de temps suivant
440c       ----------------------------------------------------------
441
442           psurf=psurf+dtphys*dpsurf(1)   ! evolution de la pression de surface
443           DO ilevel=1,nlevel
444             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
445           ENDDO
446           DO ilayer=1,nlayer
447             play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
448           ENDDO
449
450      ENDDO   ! fin de la boucle temporelle
451
452c    ========================================================
453c    GESTION DES SORTIE
454c    ========================================================
455
456        print*,"Temperature finale:"
457        print*,temp
458       
459c stabilite
460      DO ilayer=1,nlayer
461        zlay(ilayer) = phi(ilayer)/g/1000.  !en km
462      ENDDO
463      DO ilayer=2,nlayer
464        tmp1(ilayer) =
465     .    (temp(ilayer)-temp(ilayer-1))/(zlay(ilayer)-zlay(ilayer-1))
466     .   + 1000.*g/cpp
467      ENDDO
468
469      OPEN(11,file='profile.new')
470      DO ilayer=1,nlayer
471        write (11,*) zlay(ilayer),temp(ilayer),tmp1(ilayer)
472      ENDDO
473
474c    ========================================================
475      END
476 
477c***********************************************************************
478c***********************************************************************
479
480!#include "../dyn3d_common/disvert_noterre.F"
481!#include "../dyn3d/abort_gcm.F"
Note: See TracBrowser for help on using the repository browser.