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

Last change on this file since 486 was 486, checked in by rwordsworth, 13 years ago

Variable ice timestep added for ice evolution algorithm.
Sedimentation improved to allow consistent dust (L. Kerber).
Bug linked to allocatable matrices removed from rcm1d.F.
Treatment of initial aerosol radii in callcorrk.F90 improved.

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