source: trunk/LMDZ.GENERIC/libf/phystd/rcm1d.F @ 588

Last change on this file since 588 was 586, checked in by jleconte, 13 years ago

removed cpp3D and nonideal stuff.

  • Property svn:executable set to *
File size: 23.4 KB
Line 
1      program rcm1d
2
3      use radcommon_h, only: tauvis
4! to use  'getin'
5      use ioipsl_getincom
6
7      implicit none
8
9!==================================================================
10!     
11!     Purpose
12!     -------
13!     Run the physics package of the universal model in a 1D column.
14!     
15!     It can be compiled with a command like (e.g. for 25 layers):
16!                                  "makegcm -p std -d 25 rcm1d"
17!     It requires the files "callphys.def", "z2sig.def",
18!     "traceur.def", and "run.def" with a line "INCLUDEDEF=callphys.def"
19!
20!     Authors
21!     -------
22!     Frederic Hourdin
23!     R. Fournier
24!     F. Forget
25!     F. Montmessin (water ice added)
26!     R. Wordsworth
27!
28!==================================================================
29
30#include "dimensions.h"
31#include "paramet.h"
32#include "dimphys.h"
33#include "comgeomfi.h"
34#include "surfdat.h"
35#include "comsoil.h"
36#include "comdiurn.h"
37#include "callkeys.h"
38#include "comcstfi.h"
39#include "planete.h"
40#include "comsaison.h"
41#include "control.h"
42#include "comvert.h"
43#include "netcdf.inc"
44#include "watercap.h"
45#include "fisice.h"
46#include "logic.h"
47#include "advtrac.h"
48#include "comgeom.h"
49
50c --------------------------------------------------------------
51c  Declarations
52c --------------------------------------------------------------
53c
54      INTEGER unitstart      ! unite d'ecriture de "startfi"
55      INTEGER nlayer,nlevel,nsoil,ndt
56      INTEGER ilayer,ilevel,isoil,idt,iq
57      LOGICAl firstcall,lastcall
58c
59      INTEGER day0          ! date initial (sol ; =0 a Ls=0)
60      REAL day           ! date durant le run
61      REAL time             ! time (0<time<1 ; time=0.5 a midi)
62      REAL play(nlayermx)   ! Pressure at the middle of the layers (Pa)
63      REAL plev(nlayermx+1) ! intermediate pressure levels (pa)
64      REAL psurf,tsurf     
65      REAL u(nlayermx),v(nlayermx)  ! zonal, meridional wind
66      REAL gru,grv   ! prescribed "geostrophic" background wind
67      REAL temp(nlayermx)   ! temperature at the middle of the layers
68      REAL q(nlayermx,nqmx) ! tracer mixing ratio (e.g. kg/kg)
69      REAL qsurf(nqmx)      ! tracer surface budget (e.g. kg.m-2)
70      REAL tsoil(nsoilmx)   ! subsurface soik temperature (K)
71!      REAL co2ice           ! co2ice layer (kg.m-2) !not used anymore
72      integer :: i_co2_ice=0  ! tracer index of co2 ice
73      integer :: i_h2o_ice=0  ! tracer index of co2 ice
74      REAL emis             ! surface layer
75      REAL q2(nlayermx+1)   ! Turbulent Kinetic Energy
76      REAL zlay(nlayermx)   ! altitude estimee dans les couches (km)
77
78c    Physical and dynamical tandencies (e.g.  m.s-2, K/s, Pa/s)
79      REAL du(nlayermx),dv(nlayermx),dtemp(nlayermx)
80      REAL dudyn(nlayermx),dvdyn(nlayermx),dtempdyn(nlayermx)
81      REAL dpsurf   
82      REAL dq(nlayermx,nqmx)
83      REAL dqdyn(nlayermx,nqmx)
84
85c   Various intermediate variables
86!      INTEGER thermo
87      REAL zls
88      REAL phi(nlayermx),h(nlayermx),s(nlayermx)
89      REAL pks, ptif, w(nlayermx)
90      REAL qtotinit, mqtot(nqmx),qtot
91      INTEGER ierr, aslun
92      REAL tmp1(0:nlayermx),tmp2(0:nlayermx)
93      Logical  tracerdyn
94      integer :: nq=1 ! number of tracers
95 
96      character*2 str2
97      character (len=7) :: str7
98
99      logical oldcompare, earthhack,saveprofile
100
101!     added by RW for zlay computation
102      real Hscale, Hmax, rho, dz
103
104!     added by RW for autozlevs computation
105      real nu, xx, pMIN, zlev, Htop
106      real logplevs(nlayermx)
107
108!     added by BC
109      REAL cloudfrac(ngridmx,nlayermx)
110      REAL hice(ngridmx),totcloudfrac(ngridmx)
111
112c=======================================================================
113c INITIALISATION
114c=======================================================================
115
116      saveprofile=.false.
117
118c ------------------------------------------------------
119c  Default values for constants (corresponding to Mars)
120c ------------------------------------------------------
121
122      pi=2.E+0*asin(1.E+0)
123
124c     Constante de la Planete Mars
125c     ----------------------------
126      rad=3397200.               ! rayon de mars (m)  ~3397200 m
127      daysec=88775.              ! duree du sol (s)  ~88775 s
128      omeg=4.*asin(1.)/(daysec)  ! vitesse de rotation (rad.s-1)
129      g=3.72                     ! gravite (m.s-2) ~3.72 
130      mugaz=43.49                ! Masse molaire de l'atm (g.mol-1) ~43.49
131      rcp=.256793                ! = r/cp  ~0.256793
132      r= 8.314511E+0 *1000.E+0/mugaz
133      cpp= r/rcp
134      year_day = 669           ! duree de l'annee (sols) ~668.6
135      periastr = 206.66          ! min. dist. star-planet (Mkm) ~206.66
136      apoastr = 249.22           ! max. dist. star-planet (Mkm) ~249.22
137      peri_day =  485.           ! date du periastron (sols depuis printemps)
138      obliquit = 25.2           ! Obliquite de la planete (deg) ~25.2         
139 
140c     Parametres Couche limite et Turbulence
141c     --------------------------------------
142      z0 =  1.e-2                ! surface roughness (m) ~0.01
143      emin_turb = 1.e-6          ! energie minimale ~1.e-8
144      lmixmin = 30               ! longueur de melange ~100
145 
146c     propriete optiques des calottes et emissivite du sol
147c     ----------------------------------------------------
148      emissiv= 0.95              ! Emissivite du sol martien ~.95
149      emisice(1)=0.95            ! Emissivite calotte nord
150      emisice(2)=0.95            ! Emissivite calotte sud 
151
152      albedice(1)=0.5            ! Albedo calotte nord
153      albedice(2)=0.5            ! Albedo calotte sud
154      iceradius(1) = 100.e-6     ! mean scat radius of CO2 snow (north)
155      iceradius(2) = 100.e-6     ! mean scat radius of CO2 snow (south)
156      dtemisice(1) = 2.          ! time scale for snow metamorphism (north)
157      dtemisice(2) = 2.          ! time scale for snow metamorphism (south
158      hybrid=.false.
159
160c ------------------------------------------------------
161c  Load parameters from "run.def"
162c ------------------------------------------------------
163
164! also check if 'run1d.def' file is around (otherwise reading parameters
165! from callphys.def via getin() routine won't work.
166      open(90,file='run.def',status='old',form='formatted',
167     &     iostat=ierr)
168      if (ierr.ne.0) then
169        write(*,*) 'Cannot find required file "run.def"'
170        write(*,*) '  (which should contain some input parameters'
171        write(*,*) '   along with the following line:'
172        write(*,*) '   INCLUDEDEF=callphys.def'
173        write(*,*) '   )'
174        write(*,*) ' ... might as well stop here ...'
175        stop
176      else
177        close(90)
178      endif
179
180! check if we are going to run with or without tracers
181      write(*,*) "Run with or without tracer transport ?"
182      tracer=.false. ! default value
183      call getin("tracer",tracer)
184      write(*,*) " tracer = ",tracer
185
186! while we're at it, check if there is a 'traceur.def' file
187! and preocess it, if necessary. Otherwise initialize tracer names
188      if (tracer) then
189      ! load tracer names from file 'traceur.def'
190        open(90,file='traceur.def',status='old',form='formatted',
191     &       iostat=ierr)
192        if (ierr.eq.0) then
193          write(*,*) "rcm1d: Reading file traceur.def"
194          ! read number of tracers:
195          read(90,*,iostat=ierr) nq
196          if (ierr.ne.0) then
197            write(*,*) "rcm1d: error reading number of tracers"
198            write(*,*) "   (first line of traceur.def) "
199            stop
200          else
201            ! check that the number of tracers is indeed nqmx
202            if (nq.ne.nqmx) then
203              write(*,*) "rcm1d: error, wrong number of tracers:"
204              write(*,*) "nq=",nq," whereas nqmx=",nqmx
205              stop
206            endif
207          endif
208       
209          ! initialize advection schemes to Van-Leer for all tracers
210          do iq=1,nq
211            iadv(iq)=3 ! Van-Leer
212          enddo
213       
214          do iq=1,nq
215          ! minimal version, just read in the tracer names, 1 per line
216            read(90,*,iostat=ierr) tnom(iq)
217            if (ierr.ne.0) then
218              write(*,*) 'rcm1d: error reading tracer names...'
219              stop
220            endif
221          enddo !of do iq=1,nq
222! check for co2_ice / h2o_ice tracers:
223         i_co2_ice=0
224         i_h2o_ice=0
225         do iq=1,nq
226           if (tnom(iq)=="co2_ice") then
227             i_co2_ice=iq
228           elseif (tnom(iq)=="h2o_ice") then
229             i_h2o_ice=iq
230           endif
231         enddo
232         !if (i_co2_ice==0) then
233         !  write(*,*) "rcm1d: error, we need a 'co2_ice' tracer"
234         !  write(*,*) "   (add one to traceur.def)"
235         !  stop
236         !endif
237        else
238          write(*,*) 'Cannot find required file "traceur.def"'
239          write(*,*) ' If you want to run with tracers, I need it'
240          write(*,*) ' ... might as well stop here ...'
241          stop
242        endif
243        close(90)
244
245      else
246      ! we still need to set (dummy) tracer names for physdem1
247        nq=nqmx
248        do iq=1,nq
249          write(str7,'(a1,i2.2)')'q',iq
250          tnom(iq)=str7
251        enddo
252        ! actually, we'll need at least one "co2_ice" tracer
253        ! (for surface CO2 ice)
254        tnom(1)="co2_ice"
255        i_co2_ice=1
256      endif ! of if (tracer)
257
258
259c  Date et heure locale du debut du run
260c  ------------------------------------
261c    Date (en sols depuis le solstice de printemps) du debut du run
262      day0 = 0 ! default value for day0
263      write(*,*) 'Initial date (in martian sols ; =0 at Ls=0)?'
264      call getin("day0",day0)
265      day=float(day0)
266      write(*,*) " day0 = ",day0
267c  Heure de demarrage
268      time=0 ! default value for time
269      write(*,*)'Initial local time (in hours, between 0 and 24)?'
270      call getin("time",time)
271      write(*,*)" time = ",time
272      time=time/24.E+0 ! convert time (hours) to fraction of sol
273
274
275c  Discretisation (Definition de la grille et des pas de temps)
276c  --------------
277c
278      nlayer=nlayermx
279      nlevel=nlayer+1
280      nsoil=nsoilmx
281
282      day_step=48 ! default value for day_step
283      PRINT *,'Number of time steps per sol ?'
284      call getin("day_step",day_step)
285      write(*,*) " day_step = ",day_step
286
287       
288      ecritphy=day_step ! default value for ecritphy
289      PRINT *,'Nunber of steps between writediagfi ?'
290      call getin("ecritphy",ecritphy)
291      write(*,*) " ecritphy = ",ecritphy
292
293      ndt=10 ! default value for ndt
294      PRINT *,'Number of sols to run ?'
295      call getin("ndt",ndt)
296      write(*,*) " ndt = ",ndt
297
298      ndt=ndt*day_step     
299      dtphys=daysec/day_step 
300
301
302c output spectrum?
303      write(*,*)"Output spectral OLR?"
304      specOLR=.false.
305      call getin("specOLR",specOLR)
306      write(*,*)" specOLR = ",specOLR
307
308
309c Pression de surface sur la planete
310c ------------------------------------
311c
312      psurf=100000. ! default value for psurf
313      PRINT *,'Surface pressure (Pa) ?'
314      call getin("psurf",psurf)
315      write(*,*) " psurf = ",psurf
316c Pression de reference
317      pa=psurf/30.
318      preff=psurf ! these values are not needed in 1D  (are you sure JL12)
319
320c Gravity of planet
321c -----------------
322      g=10.0 ! default value for g
323      PRINT *,'Gravity ?'
324      call getin("g",g)
325      write(*,*) " g = ",g
326
327      mugaz=0.
328      cpp=0.
329
330      call su_gases
331      call calc_cpp_mugaz
332     
333
334c Proprietes des poussiere aerosol
335c --------------------------------
336      tauvis=0.2 ! default value for tauvis
337      print *,'Reference dust opacity at 700 Pa ?'
338      call getin("tauvis",tauvis)
339      write(*,*) " tauvis = ",tauvis
340 
341c  latitude/longitude
342c  ------------------
343      lati(1)=0 ! default value for lati(1)
344      PRINT *,'latitude (in degrees) ?'
345      call getin("latitude",lati(1))
346      write(*,*) " latitude = ",lati(1)
347      lati(1)=lati(1)*pi/180.E+0
348      long(1)=0.E+0
349      long(1)=long(1)*pi/180.E+0
350
351c  periastron/apoastron
352c  --------------------
353      periastr=227.94
354      apoastr=227.94 ! default = mean martian values
355      PRINT *,'periastron (in AU) ?'
356      call getin("periastr",periastr)
357      write(*,*) "periastron = ",periastr
358      periastr=periastr*149.598 ! AU to Mkm
359      PRINT *,'apoastron (in AU) ?'
360      call getin("apoastr",apoastr)
361      write(*,*) "apoastron = ",apoastr
362      apoastr=apoastr*149.598 ! AU to Mkm
363
364c  obliquity
365c  ---------
366      obliquit=0.0! default=zero in 1d
367      PRINT *,'obliquity (in AU) ?'
368      call getin("obliquit",obliquit)
369      write(*,*) "obliquity = ",obliquit
370
371c  Initialisation albedo / inertie du sol
372c  --------------------------------------
373      albedodat(1)=0.2 ! default value for albedodat
374      PRINT *,'Albedo of bare ground ?'
375      call getin("albedo",albedodat(1))
376      write(*,*) " albedo = ",albedodat(1)
377
378      inertiedat(1,1)=400 ! default value for inertiedat
379      PRINT *,'Soil thermal inertia (SI) ?'
380      call getin("inertia",inertiedat(1,1))
381      write(*,*) " inertia = ",inertiedat(1,1)
382c
383c  pour le schema d'ondes de gravite
384c  ---------------------------------
385c
386      zmea(1)=0.E+0
387      zstd(1)=0.E+0
388      zsig(1)=0.E+0
389      zgam(1)=0.E+0
390      zthe(1)=0.E+0
391
392c    Initialisation des traceurs
393c    ---------------------------
394
395      DO iq=1,nqmx
396        DO ilayer=1,nlayer
397           q(ilayer,iq) = 0.
398        ENDDO
399      ENDDO
400     
401      DO iq=1,nqmx
402        qsurf(iq) = 0.
403      ENDDO
404
405c   Initialisation speciales "physiq"
406c   ---------------------------------
407c   la surface de chaque maille est inutile en 1D --->
408      area(1)=1.E+0
409      aire(1)=area(1) !JL+EM to have access to the area in the diagfi.nc files. area in comgeomfi.h and aire in comgeom.h
410
411c   le geopotentiel au sol est inutile en 1D car tout est controle
412c   par la pression de surface --->
413      phisfi(1)=0.E+0
414
415c  "inifis" reproduit un certain nombre d'initialisations deja faites
416c  + lecture des clefs de callphys.def
417c  + calcul de la frequence d'appel au rayonnement
418c  + calcul des sinus et cosinus des longitude latitude
419
420!Mars possible issue with dtphys in input and include!!!
421      CALL inifis(1,llm,day0,daysec,dtphys,
422     .            lati,long,area,rad,g,r,cpp)
423c   Initialisation pour prendre en compte les vents en 1-D
424c   ------------------------------------------------------
425      ptif=2.E+0*omeg*sinlat(1)
426 
427
428c    vent geostrophique
429      gru=10. ! default value for gru
430      PRINT *,'zonal eastward component of the geostrophic wind (m/s) ?'
431      call getin("u",gru)
432      write(*,*) " u = ",gru
433      grv=0. !default value for grv
434      PRINT *,'meridional northward component of the geostrophic',
435     &' wind (m/s) ?'
436      call getin("v",grv)
437      write(*,*) " v = ",grv
438
439c     Initialisation des vents  au premier pas de temps
440      DO ilayer=1,nlayer
441         u(ilayer)=gru
442         v(ilayer)=grv
443      ENDDO
444
445c     energie cinetique turbulente
446      DO ilevel=1,nlevel
447         q2(ilevel)=0.E+0
448      ENDDO
449
450c  emissivity / surface co2 ice ( + h2o ice??)
451c  -------------------------------------------
452      emis=emissiv ! default value for emissivity
453      PRINT *,'Emissivity of bare ground ?'
454      call getin("emis",emis)
455      write(*,*) " emis = ",emis
456      emissiv=emis ! we do this so that condense_co2 sets things to the right
457                   ! value if there is no snow
458
459      if(i_co2_ice.gt.0)then
460         qsurf(i_co2_ice)=0 ! default value for co2ice
461         print*,'Initial CO2 ice on the surface (kg.m-2)'
462         call getin("co2ice",qsurf(i_co2_ice))
463         write(*,*) " co2ice = ",qsurf(i_co2_ice)
464         IF (qsurf(i_co2_ice).ge.1.E+0) THEN
465            ! if we have some CO2 ice on the surface, change emissivity
466            if (lati(1).ge.0) then ! northern hemisphere
467              emis=emisice(1)
468            else ! southern hemisphere
469              emis=emisice(2)
470            endif
471         ENDIF
472      endif
473
474c  calcul des pressions et altitudes en utilisant les niveaux sigma
475c  ----------------------------------------------------------------
476
477c    Vertical Coordinates
478c    """"""""""""""""""""
479      hybrid=.true.
480      PRINT *,'Hybrid coordinates ?'
481      call getin("hybrid",hybrid)
482      write(*,*) " hybrid = ", hybrid
483
484
485      autozlevs=.false.
486      PRINT *,'Auto-discretise vertical levels ?'
487      call getin("autozlevs",autozlevs)
488      write(*,*) " autozlevs = ", autozlevs
489
490      pceil=100.0 ! Pascals
491      PRINT *,'Ceiling pressure (Pa) ?'
492      call getin("pceil",pceil)
493      write(*,*) " pceil = ", pceil
494
495! Test of incompatibility:
496! if autozlevs used, cannot have hybrid too
497
498      if (autozlevs.and.hybrid) then
499         print*,'Cannot use autozlevs and hybrid together!'
500         call abort
501      endif
502
503      if(autozlevs)then
504           
505         open(91,file="z2sig.def",form='formatted')
506         read(91,*) Hscale
507         DO ilayer=1,nlayer-2
508            read(91,*) Hmax
509         enddo
510         close(91)
511 
512           
513         print*,'Hmax = ',Hmax,' km'
514         print*,'Auto-shifting Hscale to:'
515!     Hscale = Hmax / log(psurf/100.0)
516         Hscale = Hmax / log(psurf/pceil)
517         print*,'Hscale = ',Hscale,' km'
518         
519! none of this matters if we dont care about zlay
520         
521      endif
522
523      call disvert
524
525         if(.not.autozlevs)then
526            ! we want only the scale height from z2sig, in order to compute zlay
527            open(91,file="z2sig.def",form='formatted')
528            read(91,*) Hscale
529            close(91)
530         endif
531
532!         if(autozlevs)then
533!            open(94,file="Hscale.temp",form='formatted')
534!            read(94,*) Hscale
535!            close(94)
536!         endif
537
538         DO ilevel=1,nlevel
539            plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
540         ENDDO
541         
542         DO ilayer=1,nlayer
543            play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
544         ENDDO
545         
546
547
548         DO ilayer=1,nlayer
549!     zlay(ilayer)=-300.E+0 *r*log(play(ilayer)/plev(1))
550!     &   /g
551            zlay(ilayer)=-1000.0*Hscale*log(play(ilayer)/plev(1))
552         ENDDO
553
554!      endif
555
556c  profil de temperature au premier appel
557c  --------------------------------------
558      pks=psurf**rcp
559
560c altitude en km dans profile: on divise zlay par 1000
561      tmp1(0)=0.E+0
562      DO ilayer=1,nlayer
563        tmp1(ilayer)=zlay(ilayer)/1000.E+0
564      ENDDO
565      call profile(nlayer+1,tmp1,tmp2)
566
567      tsurf=tmp2(0)
568      DO ilayer=1,nlayer
569        temp(ilayer)=tmp2(ilayer)
570      ENDDO
571     
572
573
574! Initialize soil properties and temperature
575! ------------------------------------------
576      volcapa=1.e6 ! volumetric heat capacity
577      DO isoil=1,nsoil
578         inertiedat(1,isoil)=inertiedat(1,1) ! soil thermal inertia
579         tsoil(isoil)=tsurf  ! soil temperature
580      ENDDO
581
582! Initialize depths
583! -----------------
584      do isoil=0,nsoil-1
585        mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) ! mid-layer depth
586      enddo
587      do isoil=1,nsoil
588        layer(isoil)=2.e-4*(2.**(isoil-1)) ! layer depth
589      enddo
590
591
592c  Ecriture de "startfi"
593c  --------------------
594c  (Ce fichier sera aussitot relu au premier
595c   appel de "physiq", mais il est necessaire pour passer
596c   les variables purement physiques a "physiq"...
597
598      call physdem1("startfi.nc",long,lati,nsoilmx,nqmx,
599     &              dtphys,float(day0),
600     &              time,tsurf,tsoil,emis,q2,qsurf,
601     &              area,albedodat,inertiedat,zmea,zstd,zsig,zgam,zthe,
602     &              cloudfrac,totcloudfrac,hice)
603
604c=======================================================================
605c  BOUCLE TEMPORELLE DU MODELE 1D
606c=======================================================================
607
608      firstcall=.true.
609      lastcall=.false.
610
611      DO idt=1,ndt
612        IF (idt.eq.ndt) then       !test
613         lastcall=.true.
614         call stellarlong(day*1.0,zls)
615!         write(103,*) 'Ls=',zls*180./pi
616!         write(103,*) 'Lat=', lati(1)*180./pi
617!         write(103,*) 'Tau=', tauvis/700*psurf
618!         write(103,*) 'RunEnd - Atmos. Temp. File'
619!         write(103,*) 'RunEnd - Atmos. Temp. File'
620!         write(104,*) 'Ls=',zls*180./pi
621!         write(104,*) 'Lat=', lati(1)
622!         write(104,*) 'Tau=', tauvis/700*psurf
623!         write(104,*) 'RunEnd - Atmos. Temp. File'
624        ENDIF
625
626c    calcul du geopotentiel
627c     ~~~~~~~~~~~~~~~~~~~~~
628
629
630      DO ilayer=1,nlayer
631
632!              if(autozlevs)then
633!                 s(ilayer)=(play(ilayer)/psurf)**rcp
634!              else
635          s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
636!              endif
637              !s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
638          h(ilayer)=cpp*temp(ilayer)/(pks*s(ilayer))
639       ENDDO
640
641!      DO ilayer=1,nlayer
642!        s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
643!        h(ilayer)=cpp*temp(ilayer)/(pks*s(ilayer))
644!      ENDDO
645      phi(1)=pks*h(1)*(1.E+0-s(1))
646      DO ilayer=2,nlayer
647         phi(ilayer)=phi(ilayer-1)+
648     &               pks*(h(ilayer-1)+h(ilayer))*.5E+0
649     &                  *(s(ilayer-1)-s(ilayer))
650
651      ENDDO
652
653c       appel de la physique
654c       --------------------
655
656
657      CALL physiq (1,llm,nqmx,
658     ,     firstcall,lastcall,
659     ,     day,time,dtphys,
660     ,     plev,play,phi,
661     ,     u, v,temp, q, 
662     ,     w,
663C - sorties
664     s     du, dv, dtemp, dq,dpsurf,tracerdyn)
665
666
667c       evolution du vent : modele 1D
668c       -----------------------------
669 
670c       la physique calcule les derivees temporelles de u et v.
671c       on y rajoute betement un effet Coriolis.
672c
673c       DO ilayer=1,nlayer
674c          du(ilayer)=du(ilayer)+ptif*(v(ilayer)-grv)
675c          dv(ilayer)=dv(ilayer)+ptif*(-u(ilayer)+gru)
676c       ENDDO
677
678c       Pour certain test : pas de coriolis a l'equateur
679c       if(lati(1).eq.0.) then
680          DO ilayer=1,nlayer
681             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
682             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
683          ENDDO
684c       end if
685c     
686c
687c       Calcul du temps au pas de temps suivant
688c       ---------------------------------------
689        firstcall=.false.
690        time=time+dtphys/daysec
691        IF (time.gt.1.E+0) then
692            time=time-1.E+0
693            day=day+1
694        ENDIF
695
696c       calcul des vitesses et temperature au pas de temps suivant
697c       ----------------------------------------------------------
698
699        DO ilayer=1,nlayer
700           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
701           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
702           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
703        ENDDO
704
705c       calcul des pressions au pas de temps suivant
706c       ----------------------------------------------------------
707
708           psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
709           DO ilevel=1,nlevel
710              plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
711           ENDDO
712           DO ilayer=1,nlayer
713                 play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
714           ENDDO
715
716c       calcul traceur au pas de temps suivant
717c       --------------------------------------
718
719        DO iq = 1, nqmx
720          DO ilayer=1,nlayer
721             q(ilayer,iq)=q(ilayer,iq)+dtphys*dq(ilayer,iq)
722          ENDDO
723        END DO
724
725      ENDDO   ! fin de la boucle temporelle
726
727c    ========================================================
728c    GESTION DES SORTIE
729c    ========================================================
730
731! save temperature profile
732      if(saveprofile)then
733         OPEN(12,file='profile.out',form='formatted')
734         DO ilayer=1,nlayermx
735            write(12,*) temp(ilayer), play(ilayer)
736         ENDDO
737         CLOSE(12)
738      endif
739
740
741c    ========================================================
742      end                       !rcm1d
743 
744c***********************************************************************
745c***********************************************************************
746c     Subroutines Bidons utilise seulement en 3D, mais
747c     necessaire a la compilation de rcm1d en 1D
748
749      subroutine gr_fi_dyn
750      RETURN
751      END
752 
753c***********************************************************************
754c***********************************************************************
755
756#include "../dyn3d/disvert.F"
Note: See TracBrowser for help on using the repository browser.