source: trunk/LMDZ.VENUS/libf/phyvenus/testphys1d.F @ 856

Last change on this file since 856 was 849, checked in by emillour, 13 years ago

Work on common dynamics and interfacing with different physics:

  • Put calls to PVtheta in dynamics between CPP_EARTH flags (because it calls tetalevel, which is supposed to be in the physics; only OK for Earth...).
  • Adapted makelmdz script so that one can compile main programs in dyn* or phy* (makelmdz_fcm already capable of doing that).
  • Moved start2archive-VT.F and start2archive-VT.F to phyvenus (as start2archive.F and newstart.F); leave adapting them to Titan for later.
  • Small correction to phyvenus/testphys1d.F (use module control_mod instead of control.h and call disvert_noterre).
  • removed "use histcom" in phyvenus/physiq.F and phytitan/physiq.F ; it is not needed since there is already a "use ioipsl" (and it moreover confused makelmdz_fcm...)
  • Had to add declaration of variable "zlsm1" in phytitan/physiq.F because it is used in "write_hist.h"; but note that it is used while not initialized (but what should it be initialized to?).

EM

File size: 14.5 KB
Line 
1
2      PROGRAM testphys1d
3     
4      use control_mod
5      IMPLICIT NONE
6
7c=======================================================================
8c   subject:
9c   --------
10c   PROGRAM useful to run physical part of the venusian GCM in a 1D column
11c       
12c Can be compiled with a command like (e.g. for 50 layers)
13c  "makegcm -p venus -d 50 testphys1d"
14
15c A REVOIR POUR VENUS....
16c It requires the files "testphys1d.def" "callphys.def"
17c      and a file describing the sigma layers (e.g. "z2sig.def")
18c
19c   author: Frederic Hourdin, R.Fournier,F.Forget (original Mars version)
20c   ------- Sebastien Lebonnois (Venus version)
21c   
22c=======================================================================
23
24c VENUS: revoir tous les include...... par rapport a calfis.F
25c        voir aussi les arguments de calfis: clesphy0
26c        voir aussi dans les include de physiq.F les variables qui
27c             doivent etre initialisees...
28#include "dimensions.h"
29#include "dimsoil.h"
30#include "comcstfi.h"
31!#include "control.h"
32#include "comvert.h"
33#include "netcdf.inc"
34#include "comg1d.h"
35#include "logic.h"
36#include "clesphys.h"
37#include "iniprint.h"
38
39c --------------------------------------------------------------
40c  Declarations
41c --------------------------------------------------------------
42c
43      INTEGER unit           ! unite de lecture de "testphys1d.def"
44      INTEGER unitstart      ! unite d'ecriture de "startphy.nc"
45      INTEGER nlayer,nlevel,nsoil,ndt
46      INTEGER ilayer,ilevel,isoil,idt,iq
47      LOGICAl firstcall,lastcall
48c
49      INTEGER day0          ! date initial (sol ; =0 a Ls=0)
50      REAL day              ! date durant le run
51      REAL time             ! time (0<time<1 ; time=0.5 a midi)
52      REAL play(llm)   ! Pressure at the middle of the layers (Pa)
53      REAL plev(llm+1) ! intermediate pressure levels (pa)
54      REAL psurf,tsurf(1)     
55      REAL u(llm),v(llm)  ! zonal, meridional wind
56      REAL gru,grv   ! prescribed "geostrophic" background wind
57      REAL temp(llm)   ! temperature at the middle of the layers
58      REAL q(llm,nqtot) ! tracer mixing ratio (e.g. kg/kg)
59      REAL tsoil(nsoilmx)   ! subsurface soik temperature (K)
60      REAL zlay(llm)   ! altitude estimee dans les couches (km)
61      REAL long(1),lati(1),area(1)
62      REAL cufi(1),cvfi(1)
63      REAL phisfi(1),albedo(1)
64      REAL solsw(1),sollwdown(1),dlw(1),radsol(1)
65
66c    Physical and dynamical tandencies (e.g.  m.s-2, K/s, Pa/s)
67      REAL du(llm),dv(llm),dtemp(llm)
68      REAL dudyn(llm),dvdyn(llm),dtempdyn(llm)
69      REAL dpsurf   
70      REAL dq(llm,nqtot)
71      REAL dqdyn(llm,nqtot)
72
73c   Various intermediate variables
74      REAL zls
75      REAL phi(llm),s(llm),aps(llm),bps(llm)
76      REAL pk(llm),pks, w(llm)
77      INTEGER l, ierr, aslun,radpas
78      REAL tmp1(0:llm),tmp2(0:llm)
79      INTEGER        longcles
80      PARAMETER    ( longcles = 20 )
81      REAL clesphy0( longcles      )
82                       
83
84      character*2 str2
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      rad=6051300.               ! rayon de Venus (m)  ~6051300 m
103      daysec= 1.0087e7           ! duree du sol (s)  ~1.e7 s
104      omeg=4.*asin(1.)/19.4141e6 ! vitesse de rotation (rad.s-1)
105      g= 8.87                    ! gravite (m.s-2) ~8.87
106      mugaz=43.44                ! Masse molaire de l'atm (g.mol-1) ~43.44
107! ADAPTATION GCM POUR CP(T)
108! VENUS: Cp(T) = cpp*(T/T0)^nu
109! avec T0=460. et nu=0.35
110!     cpp=1.0e3
111      cpp=9.0e2      ! version constante
112      r= 8.314511E+0 *1000.E+0/mugaz
113      rcp= r/cpp
114
115c ------------------------------------------------------
116c  Lecture des parametres dans "testphys1d.def"
117c ------------------------------------------------------
118
119c   Opening parameters file "testphys1d.def"
120c   ---------------------------------------
121      unit =97
122      OPEN(unit,file='testphys1d.def',status='old',form='formatted'
123     .     ,iostat=ierr)
124
125      IF(ierr.ne.0) THEN
126        write(*,*) 'Problem to open "testphys1d.def'
127        write(*,*) 'Is it there ?'
128        stop
129      END IF
130
131c  Date et heure locale du debut du run
132c  ------------------------------------
133c    Date (en sols depuis le solstice de printemps) du debut du run
134      day0 = 0
135      PRINT *,'date de depart ?'
136      READ(unit,*) day0
137      day=float(day0)
138      PRINT *,day0
139c  Heure de demarrage
140      PRINT *,'heure de debut de simulation (entre 0 et 24) ?'
141      READ(unit,*) time
142      time=time/24.E+0
143
144c  Discretisation (Definition de la grille et des pas de temps)
145c  --------------
146c
147      nlayer=llm
148      nlevel=nlayer+1
149      nsoil=nsoilmx
150      PRINT *,'nombre de pas de temps par jour ?'
151      READ(unit,*) day_step
152      print*,day_step
153
154      PRINT *,'nombre d appel au rayonnement par jour ?'
155      READ(unit,*) nbapp_rad
156      print*,nbapp_rad
157
158      PRINT *,'nombre de pas physique entre chaque sortie ?'
159      READ(unit,*) nbphypers
160      print*,nbphypers
161
162      PRINT *,'nombre de jours simules ?'
163      READ(unit,*) ndt
164      print*,ndt
165
166      ndt=ndt*day_step     
167      dtphys=daysec/day_step 
168
169c Pression de surface sur la planete
170c ------------------------------------
171c
172      PRINT *,'pression au sol'
173      READ(unit,*) psurf
174      PRINT *,psurf
175c Pression de reference  ! voir dyn3d/etat0_venus
176c     pa     =  5.e4
177      pa     =  1.e6
178      preff  = 9.2e6 ! 92 bars
179c     preff  = psurf
180 
181c  latitude/longitude
182c  -------------------
183      PRINT *,'latitude en degres ?'
184      READ(unit,*) lati(1)
185      PRINT *,lati(1)
186      long(1)=0.E+0
187
188c  Initialisation albedo
189c  ----------------------
190c
191      PRINT *,'Albedo du sol nu ?'
192      READ(unit,*) albedo(1)
193      PRINT *,albedo(1)
194
195c   Initialisation speciales "physiq"
196c   ---------------------------------
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   le geopotentiel au sol est inutile en 1D car tout est controle
204c   par la pression de surface --->
205      phisfi(1)=0.E+0
206
207      CALL iniphysiq(1,llm,daysec,day0,dtphys,
208     .            lati,long,area,cufi,cvfi,rad,g,r,cpp)
209
210c   Initialisation pour prendre en compte les vents en 1-D
211c   ------------------------------------------------------
212 
213c    vent geostrophique
214      PRINT *,'composante vers l est du vent geostrophique (U) ?'
215      READ(unit,*) gru
216      PRINT *,'composante vers le nord du vent geostrophique (V) ?'
217      READ(unit,*) grv
218
219c     Initialisation des vents  au premier pas de temps
220      DO ilayer=1,nlayer
221         u(ilayer)=gru
222         v(ilayer)=grv
223      ENDDO
224
225c  calcul des pressions et altitudes en utilisant les niveaux sigma
226c  ----------------------------------------------------------------
227
228c    Vertical Coordinates  (hybrids)
229c    """"""""""""""""""""
230      CALL  disvert_noterre
231     
232c     Calcul au milieu des couches : Vient de la version Mars
233c     WARNING : le choix de placer le milieu des couches au niveau de
234c     pression intermédiaire est arbitraire et pourrait etre modifié.
235c     C'est fait de la meme facon dans disvert
236
237      DO l = 1, llm
238       aps(l) =  0.5 *( ap(l) +ap(l+1))
239       bps(l) =  0.5 *( bp(l) +bp(l+1))
240      ENDDO
241
242      DO ilevel=1,nlevel
243        plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
244      ENDDO
245
246      DO ilayer=1,nlayer
247        play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
248        pk(ilayer)  =cpp*(play(ilayer)/preff)**rcp
249c       write(120,*) ilayer,plev(ilayer),play(ilayer)
250      ENDDO
251c     write(120,*) nlevel,plev(nlevel)
252c     stop
253     
254      pks=cpp*(psurf/preff)**rcp
255
256c  profil de temperature et altitude au premier appel
257c  --------------------------------------------------
258
259c modif par rapport a Mars:
260c   on envoie dz/T=-log(play/psurf)*r/g dans profile
261      tmp1(0)=0.0
262      tmp1(1)= -log(play(1)/psurf)*r/g
263      DO ilayer=2,nlayer
264        tmp1(ilayer)=-log(play(ilayer)/play(ilayer-1))*r/g
265      ENDDO
266      call profile(unit,nlayer+1,tmp1,tmp2)
267c file testphys1d.def closed in profile!
268
269      print*,"               Pression        Altitude     Temperature"
270      ilayer=1
271      tsurf(1)=tmp2(0)
272       temp(1)=tmp2(1)
273       zlay(1)=tmp2(1)*tmp1(1)
274      print*,"           0",tsurf(1)
275      print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
276      DO ilayer=2,nlayer
277        temp(ilayer)=tmp2(ilayer)
278        zlay(ilayer)=zlay(ilayer-1)+tmp2(ilayer)*tmp1(ilayer)
279        print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
280      ENDDO
281     
282c     temperature du sous-sol
283c     ~~~~~~~~~~~~~~~~~~~~~~~
284      DO isoil=1,nsoil
285         tsoil(isoil)=tsurf(1)
286      ENDDO
287
288c    Initialisation des traceurs
289c    ---------------------------
290
291      DO iq=1,nqtot
292        DO ilayer=1,nlayer
293           q(ilayer,iq) = 0.
294        ENDDO
295      ENDDO
296
297c    Initialisation pour les sorties GRADS dans "g1d.dat" et "g1d.ctl"
298c    ----------------------------------------------------------------
299c    (effectuee avec "writeg1d" a partir de "physiq.F"
300c    ou tout autre subroutine physique, y compris celle ci).
301
302        g1d_nlayer=nlayer
303        g1d_nomfich='g1d.dat'
304        g1d_unitfich=40
305        g1d_nomctl='g1d.ctl'
306        g1d_unitctl=41
307        g1d_premier=.true.
308        g2d_premier=.true.
309
310c  Ecriture de "startphy.nc"
311c  -------------------------
312c  (Ce fichier sera aussitot relu au premier
313c   appel de "physiq", mais il est necessaire pour passer
314c   les variables purement physiques a "physiq"...
315
316      clesphy0(1) = 0.
317      clesphy0(2) = nbapp_rad
318      clesphy0(3) = 0. ! cycle diurne
319      clesphy0(4) = 1. ! soil_model
320      clesphy0(5) = 1. ! new_oliq
321      clesphy0(6) = 0.
322      clesphy0(7) = 0.
323      clesphy0(8) = 0.
324
325      solsw(1)    = 0.
326      sollwdown(1)= 0.
327      dlw(1)      = 0.
328      radsol(1)   = 0.
329     
330      radpas      = NINT(day_step/clesphy0(2))
331      soil_model  = .true.
332      new_oliq    = .true.
333
334      call phyredem("startphy.nc  ",dtphys,radpas,
335     .              lati,long,
336     .              tsurf,tsoil,albedo,
337     .              solsw,sollwdown,dlw,radsol,
338     .              temp,q)
339
340c=======================================================================
341c  BOUCLE TEMPORELLE DU MODELE 1D
342c=======================================================================
343c
344      firstcall=.true.
345      lastcall=.false.
346
347      DO idt=1,ndt
348c        IF (idt.eq.ndt) lastcall=.true.
349        IF (idt.eq.ndt-day_step-1) then       !test
350         lastcall=.true.
351c toujours nulle dans le cas de Venus, pour l'instant...
352         zls = 0.0
353c        write(103,*) 'Ls=',zls*180./pi
354c        write(103,*) 'Lat=', lati(1)
355c        write(103,*) 'RunEnd - Atmos. Temp. File'
356c        write(103,*) 'RunEnd - Atmos. Temp. File'
357c        write(104,*) 'Ls=',zls*180./pi
358c        write(104,*) 'Lat=', lati(1)
359c        write(104,*) 'RunEnd - Atmos. Temp. File'
360        ENDIF
361
362c    calcul du geopotentiel
363c     ~~~~~~~~~~~~~~~~~~~~~
364! ADAPTATION GCM POUR CP(T)
365      DO ilayer=1,nlayer
366        s(ilayer)=(play(ilayer)/psurf)**rcp
367      ENDDO
368      phi(1)=cpp*temp(1)*(1.E+0-s(1))
369      DO ilayer=2,nlayer
370         phi(ilayer)=phi(ilayer-1)+
371     &     cpp*(temp(ilayer-1)/s(ilayer-1)+temp(ilayer)/s(ilayer))*0.5
372     &        *(s(ilayer-1)-s(ilayer))
373
374      ENDDO
375
376c       appel de la physique
377c       --------------------
378
379      CALL physiq (1,llm,nqtot,
380     ,     firstcall,lastcall,
381     ,     day,time,dtphys,
382     ,     plev,play,pk,phi,phisfi,
383     ,     presnivs,clesphy0,
384     ,     u,v,temp,q, 
385     ,     w,
386C - sorties
387     s     du,dv,dtemp,dq,dpsurf)
388
389c     print*,"DT APRES PHYSIQ=",day,time
390c     print*,dtemp
391c     print*,temp
392c     print*," "
393c     stop
394
395c       evolution du vent : modele 1D
396c       -----------------------------
397 
398c       la physique calcule les derivees temporelles de u et v.
399c       Pas de coriolis
400          DO ilayer=1,nlayer
401             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
402             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
403          ENDDO
404c     
405c         call writeg1d(1,nlayer,du,"du","du")
406c         call writeg1d(1,nlayer,dv,"dv","dv")
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
425c         call writeg1d(1,nlayer,u,"u","Vent zonal")
426c         call writeg1d(1,nlayer,v,"v","Vent meridien")
427
428c       calcul des pressions au pas de temps suivant
429c       ----------------------------------------------------------
430
431           psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
432           DO ilevel=1,nlevel
433             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
434           ENDDO
435           DO ilayer=1,nlayer
436             play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
437           ENDDO
438
439           call writeg1d(1,nlayer,play,"press","Pressure")
440           call writeg1d(1,nlayer,phi,"phi","Geopot")
441
442      ENDDO   ! fin de la boucle temporelle
443
444c    ========================================================
445c    GESTION DES SORTIE
446c    ========================================================
447
448c    fermeture pour conclure les sorties format grads dans "g1d.dat"
449c    et "g1d.ctl" (effectuee avec "writeg1d" a partir de "physiq.F"
450c    ou tout autre subroutine physique, y compris celle ci).
451
452        CALL endg1d(1,nlayer,zlay/1000.,ndt)
453        print*,"Temperature finale:"
454        print*,temp
455       
456c stabilite
457      DO ilayer=1,nlayer
458        zlay(ilayer) = phi(ilayer)/8870.  !en km
459      ENDDO
460      DO ilayer=2,nlayer
461        tmp1(ilayer) =
462     .    (temp(ilayer)-temp(ilayer-1))/(zlay(ilayer)-zlay(ilayer-1))
463     .   + 8870./cpp
464      ENDDO
465
466      OPEN(11,file='profile.new')
467      write (11,*) tsurf
468      DO ilayer=1,nlayer
469        write (11,*) zlay(ilayer),temp(ilayer),tmp1(ilayer)
470      ENDDO
471
472c    ========================================================
473      END
474 
475c***********************************************************************
476c***********************************************************************
477
478#include "../dyn3d/disvert_noterre.F"
479#include "../dyn3d/abort_gcm.F"
480#include "../dyn3d/dump2d.F"
481
Note: See TracBrowser for help on using the repository browser.