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

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

Several new files added as part of the climate evolution model
(main program kcm.F90). Some general cleanup in physiq.F90 and
callcorrk.F90. Bugs in dust radiative transfer and H2 Rayleigh
scattering corrected.

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