source: trunk/LMDZ.TITAN/libf/phytitan/rcm1d.F @ 898

Last change on this file since 898 was 894, checked in by slebonnois, 12 years ago

SL: small modifications on rcm1d (Venus, Titan)

File size: 14.1 KB
Line 
1
2      PROGRAM rcm1d
3     
4      USE infotrac
5      use control_mod
6      use comgeomphy
7      IMPLICIT NONE
8
9c=======================================================================
10c   subject:
11c   --------
12c   PROGRAM useful to run physical part of the venusian GCM in a 1D column
13c       
14c Can be compiled with a command like (e.g. for 55 layers)
15c  "makelmdz -p titan -d 55 rcm1d"
16
17c It requires the files "rcm1d.def" "physiq.def" "traceur.def"
18c      and a file describing the sigma layers (e.g. "z2sig.def")
19c
20c   author: Frederic Hourdin, R.Fournier,F.Forget (original Mars version)
21c   ------- Sebastien Lebonnois (Venus version)
22c   
23c=======================================================================
24
25c Version TITAN a tester et verifier 
26c  - verifier pour Ls...
27c  - faire un profile.F ...
28
29#include "dimensions.h"
30#include "dimsoil.h"
31#include "comcstfi.h"
32#include "comvert.h"
33#include "netcdf.inc"
34#include "logic.h"
35#include "clesphys.h"
36#include "iniprint.h"
37#include "tabcontrol.h"
38
39c --------------------------------------------------------------
40c  Declarations
41c --------------------------------------------------------------
42c
43      INTEGER unit           ! unite de lecture de "rcm1d.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,allocatable :: q(:,:) ! 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      REAL zmea(1), zstd(1)
66      REAL zsig(1), zgam(1), zthe(1)
67      REAL zpic(1), zval(1)
68
69c    Physical and dynamical tandencies (e.g.  m.s-2, K/s, Pa/s)
70      REAL du(llm),dv(llm),dtemp(llm)
71      REAL dudyn(llm),dvdyn(llm),dtempdyn(llm)
72      REAL dpsurf   
73      REAL,allocatable :: dq(:,:)
74
75c   Various intermediate variables
76      REAL zls
77      REAL phi(llm),s(llm)
78      REAL pk(llm),pks, w(llm)
79      INTEGER l, ierr, aslun
80      REAL tmp1(0:llm),tmp2(0:llm),tmp3(0:llm)                       
81
82      character*2 str2
83
84c normalement dans dyn3d/comconst.h
85      COMMON/cpdetvenus/cppdyn,nu_venus,t0_venus
86      REAL cppdyn,nu_venus,t0_venus
87      real pi
88
89c=======================================================================
90
91c=======================================================================
92c INITIALISATION
93c=======================================================================
94
95      lunout = 6
96
97c ------------------------------------------------------
98c  Constantes prescrites ICI
99c ------------------------------------------------------
100
101      pi=2.E+0*asin(1.E+0)
102
103c     Constante de Titan
104c     -------------------
105      planet_type = "titan"
106      rad=2575000.               ! rayon de Venus (m)
107      daysec=1.37889e6           ! duree du sol (s) 
108      omeg=2*pi/daysec           ! vitesse de rotation (rad.s-1)
109      g= 1.35                    ! gravite (m.s-2)
110      mugaz=28.                  ! Masse molaire de l'atm (g.mol-1)
111      cpp=1.039e3
112      cppdyn=cpp
113      r= 8.314511E+0 *1000.E+0/mugaz
114      rcp= r/cpp
115
116c-----------------------------------------------------------------------
117c   Initialisation des traceurs
118c   ---------------------------
119c  Choix du nombre de traceurs et du schema pour l'advection
120c  dans fichier traceur.def
121      call infotrac_init
122
123c Allocation de la tableau q : champs advectes   
124      allocate(q(llm,nqtot))
125      allocate(dq(llm,nqtot))
126
127c ------------------------------------------------------
128c  Lecture des parametres dans "rcm1d.def"
129c ------------------------------------------------------
130
131c   Opening parameters file "rcm1d.def"
132c   ---------------------------------------
133      unit =97
134      OPEN(unit,file='rcm1d.def',status='old',form='formatted'
135     .     ,iostat=ierr)
136
137      IF(ierr.ne.0) THEN
138        write(*,*) 'Problem to open "rcm1d.def'
139        write(*,*) 'Is it there ?'
140        stop
141      END IF
142
143c  Date et heure locale du debut du run
144c  ------------------------------------
145c    Date (en sols depuis le solstice de printemps) du debut du run
146      day0 = 0
147      PRINT *,'date de depart ?'
148      READ(unit,*) day0
149      day=REAL(day0)
150      PRINT *,day0
151c  Heure de demarrage
152      PRINT *,'heure de debut de simulation (entre 0 et 24) ?'
153      READ(unit,*) time
154      time=time/24.E+0
155
156c  Discretisation (Definition de la grille et des pas de temps)
157c  --------------
158c
159      nlayer=llm
160      nlevel=nlayer+1
161      nsoil=nsoilmx
162      PRINT *,'nombre de pas de temps par jour ?'
163      READ(unit,*) day_step
164      print*,day_step
165
166c     PRINT *,'nombre d appel au rayonnement par jour ?'
167c     READ(unit,*) nbapp_rad
168c     print*,nbapp_rad
169c LU DANS PHYSIQ.DEF...
170      nbapp_rad = 10.
171
172      PRINT *,'nombre de jours simules ?'
173      READ(unit,*) ndt
174      print*,ndt
175
176      ndt=ndt*day_step     
177      dtphys=daysec/day_step 
178
179c Pression de surface sur la planete
180c ------------------------------------
181c
182      PRINT *,'pression au sol'
183      READ(unit,*) psurf
184      PRINT *,psurf
185c Pression de reference 
186      pa     =  1.e4
187      preff  = 1.4e5
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 albedo
198c  ----------------------
199c ne sert pas ici...
200      albedo(1)=0.3
201c      PRINT *,'Albedo du sol nu ?'
202c      READ(unit,*) albedo(1)
203c      PRINT *,albedo(1)
204
205c   Initialisation speciales "physiq"
206c   ---------------------------------
207
208      CALL init_phys_lmdz(iim,jjm,llm,1,(/1/))
209      call initcomgeomphy
210      call infotrac_init
211      call ini_cpdet
212
213c   la surface de chaque maille est inutile en 1D --->
214      area(1)=1.E+0
215c de meme ?
216      cufi(1)=1.E+0
217      cvfi(1)=1.E+0
218
219c   le geopotentiel au sol est inutile en 1D car tout est controle
220c   par la pression de surface --->
221      phisfi(1)=0.E+0
222
223      CALL iniphysiq(1,llm,daysec,day0,dtphys,
224     .            lati,long,area,cufi,cvfi,rad,g,r,cpp)
225
226c   Initialisation pour prendre en compte les vents en 1-D
227c   ------------------------------------------------------
228 
229c    vent geostrophique
230      PRINT *,'composante vers l est du vent geostrophique (U) ?'
231      READ(unit,*) gru
232      PRINT *,'composante vers le nord du vent geostrophique (V) ?'
233      READ(unit,*) grv
234
235c     Initialisation des vents  au premier pas de temps
236      DO ilayer=1,nlayer
237         u(ilayer)=gru
238         v(ilayer)=grv
239      ENDDO
240
241c  calcul des pressions et altitudes en utilisant les niveaux sigma
242c  ----------------------------------------------------------------
243
244c    Vertical Coordinates  (hybrids)
245c    """"""""""""""""""""
246      CALL  disvert_noterre
247     
248c     Calcul au milieu des couches : Vient de la version Mars
249c     WARNING : le choix de placer le milieu des couches au niveau de
250c     pression intermédiaire est arbitraire et pourrait etre modifié.
251c     C'est fait de la meme facon dans disvert
252
253      DO l = 1, llm
254       aps(l) =  0.5 *( ap(l) +ap(l+1))
255       bps(l) =  0.5 *( bp(l) +bp(l+1))
256      ENDDO
257
258      DO ilevel=1,nlevel
259        plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
260      ENDDO
261
262      DO ilayer=1,nlayer
263        play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
264        pk(ilayer)  =cpp*(play(ilayer)/preff)**rcp
265c       write(120,*) ilayer,plev(ilayer),play(ilayer)
266      ENDDO
267c     write(120,*) nlevel,plev(nlevel)
268c     stop
269     
270      pks=cpp*(psurf/preff)**rcp
271
272c  profil de temperature et altitude au premier appel
273c  --------------------------------------------------
274
275c modif par rapport a Mars:
276c   on envoie dz/T=-log(play/psurf)*r/g dans profile
277      tmp1(0)=0.0
278      tmp1(1)= -log(play(1)/psurf)*r/g
279      DO ilayer=2,nlayer
280        tmp1(ilayer)=-log(play(ilayer)/play(ilayer-1))*r/g
281      ENDDO
282      DO ilayer=0,nlayer
283        tmp2(ilayer)=plev(ilayer+1)
284      ENDDO
285      call profile(unit,nlayer+1,tmp1,tmp2,tmp3)
286      CLOSE(unit)
287
288      print*,"               Pression        Altitude     Temperature"
289      ilayer=1
290      tsurf(1)=tmp3(0)
291       temp(1)=tmp3(1)
292       zlay(1)=tmp3(1)*tmp1(1)
293      print*,"           0",tsurf(1)
294      print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
295      DO ilayer=2,nlayer
296        temp(ilayer)=tmp3(ilayer)
297        zlay(ilayer)=zlay(ilayer-1)+tmp3(ilayer)*tmp1(ilayer)
298        print*,ilayer,play(ilayer),zlay(ilayer),temp(ilayer)
299      ENDDO
300     
301c     temperature du sous-sol
302c     ~~~~~~~~~~~~~~~~~~~~~~~
303      DO isoil=1,nsoil
304         tsoil(isoil)=93.
305      ENDDO
306
307c    Initialisation des traceurs
308c    ---------------------------
309
310      DO iq=1,nqtot
311        DO ilayer=1,nlayer
312           q(ilayer,iq) = 0.
313        ENDDO
314      ENDDO
315
316c    Initialisation des parametres d'oro
317c    -----------------------------------
318
319      zmea(1) = 0.
320      zstd(1) = 0.
321      zsig(1) = 0.
322      zgam(1) = 0.
323      zthe(1) = 0.
324      zpic(1) = 0.
325      zval(1) = 0.
326
327c  Initialisation Ls
328c ------------------
329         zls=0.
330         print*,'Ls=',zls*180./pi
331
332c  Ecriture de "startphy.nc"
333c  -------------------------
334c  (Ce fichier sera aussitot relu au premier
335c   appel de "physiq", mais il est necessaire pour passer
336c   les variables purement physiques a "physiq"...
337
338      solsw(1)    = 0.
339      sollwdown(1)= 0.
340      dlw(1)      = 0.
341      radsol(1)   = 0.
342     
343      radpas      = NINT(1.*day_step/nbapp_rad)
344      soil_model  = .true.
345
346      call phyredem("startphy.nc  ",
347     .              lati,long,
348     .              tsurf,tsoil,albedo,
349     .              solsw,sollwdown,dlw,radsol,
350     .    zmea, zstd, zsig, zgam, zthe, zpic, zval,
351     .              temp)
352
353c=======================================================================
354c  BOUCLE TEMPORELLE DU MODELE 1D
355c=======================================================================
356c
357      firstcall=.true.
358      lastcall=.false.
359
360      DO idt=1,ndt
361        IF (idt.eq.ndt) then
362         lastcall=.true.
363c        write(103,*) 'Ls=',zls*180./pi
364c        write(103,*) 'Lat=', lati(1)
365c        write(103,*) 'RunEnd - Atmos. Temp. File'
366c        write(103,*) 'RunEnd - Atmos. Temp. File'
367c        write(104,*) 'Ls=',zls*180./pi
368c        write(104,*) 'Lat=', lati(1)
369c        write(104,*) 'RunEnd - Atmos. Temp. File'
370        ENDIF
371
372c    calcul du geopotentiel
373c     ~~~~~~~~~~~~~~~~~~~~~
374! ADAPTATION GCM POUR CP(T)
375      DO ilayer=1,nlayer
376        s(ilayer)=(play(ilayer)/psurf)**rcp
377      ENDDO
378      phi(1)=cpp*temp(1)*(1.E+0-s(1))
379      DO ilayer=2,nlayer
380         phi(ilayer)=phi(ilayer-1)+
381     &     cpp*(temp(ilayer-1)/s(ilayer-1)+temp(ilayer)/s(ilayer))*0.5
382     &        *(s(ilayer-1)-s(ilayer))
383
384      ENDDO
385
386c       appel de la physique
387c       --------------------
388
389      CALL physiq (1,llm,nqtot,
390     ,     firstcall,lastcall,
391     ,     day,time,dtphys,
392     ,     plev,play,pk,phi,phisfi,
393     ,     presnivs,
394     ,     u,v,temp,q, 
395     ,     w,
396C - sorties
397     s     du,dv,dtemp,dq,dpsurf)
398
399c     print*,"DT APRES PHYSIQ=",day,time
400c     print*,dtemp
401c     print*,temp
402c     print*," "
403c     stop
404
405c       evolution du vent : modele 1D
406c       -----------------------------
407 
408c       la physique calcule les derivees temporelles de u et v.
409c       Pas de coriolis
410          DO ilayer=1,nlayer
411             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
412             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
413          ENDDO
414c     
415c       Calcul du temps au pas de temps suivant
416c       ---------------------------------------
417        firstcall=.false.
418        time=time+dtphys/daysec
419        IF (time.gt.1.E+0) then
420            time=time-1.E+0
421            day=day+1
422        ENDIF
423
424c       calcul des vitesses et temperature au pas de temps suivant
425c       ----------------------------------------------------------
426
427        DO ilayer=1,nlayer
428           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
429           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
430           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
431        ENDDO
432
433c       calcul des pressions au pas de temps suivant
434c       ----------------------------------------------------------
435
436           psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
437           DO ilevel=1,nlevel
438             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
439           ENDDO
440           DO ilayer=1,nlayer
441             play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
442           ENDDO
443
444      ENDDO   ! fin de la boucle temporelle
445
446c    ========================================================
447c    GESTION DES SORTIE
448c    ========================================================
449
450        print*,"Temperature finale:"
451        print*,temp
452       
453c stabilite
454      DO ilayer=1,nlayer
455        zlay(ilayer) = phi(ilayer)/g/1000.  !en km
456      ENDDO
457      DO ilayer=2,nlayer
458        tmp1(ilayer) =
459     .    (temp(ilayer)-temp(ilayer-1))/(zlay(ilayer)-zlay(ilayer-1))
460     .   + 1000.*g/cpp
461      ENDDO
462
463      OPEN(11,file='profile.new')
464      write (11,*) tsurf
465      DO ilayer=1,nlayer
466        write (11,*) zlay(ilayer),temp(ilayer),tmp1(ilayer)
467      ENDDO
468
469c    ========================================================
470      END
471 
472c***********************************************************************
473c***********************************************************************
474
475#include "../dyn3d/disvert_noterre.F"
476#include "../dyn3d/abort_gcm.F"
477#include "../dyn3d/dump2d.F"
478#include "../dyn3d/cpdet.F"
479
480c***********************************************************************
481      function ssum(n,sx,incx)
482c
483      IMPLICIT NONE
484c
485      integer n,incx,i,ix
486      real ssum,sx((n-1)*incx+1)
487c
488      ssum=0.
489      ix=1
490      do 10 i=1,n
491         ssum=ssum+sx(ix)
492         ix=ix+incx
49310    continue
494c
495      return
496      end
Note: See TracBrowser for help on using the repository browser.