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

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

17/07/2012 == JL for LK

  • Generalization of aerosol scheme:
    • any number of aerosols can be used and id numbers are determined consistently by the code. Aerosol order not important anymore.
    • addition of a module with the id numbers for aerosols (aerosol_mod.F90).
    • initialization of aerosols id numbers in iniaerosol.F90
    • compile with -s x where x *must* be equal to the number of aerosols turned on in callphys.def (either by a flag or by dusttau>0 for dust). => may have to erase object files when compiling with s option for the first time.
  • For no aerosols, run with aeroco2=.true. and aerofixco2=.true (the default distribution for fixed co2

aerosols is 1.e-9; can be changed in aeropacity).

  • If starting from an old start file, recreate start file with the q=0 option in newstart.e.
  • update callphys.def with aeroXXX and aerofixXXX options (only XXX=co2,h2o supported for

now). Dust is activated by setting dusttau>0. See the early mars case in deftank.

  • To add other aerosols, see Laura Kerber.
  • Property svn:executable set to *
File size: 25.8 KB
Line 
1      program rcm1d
2
3! to use  'getin'
4      use ioipsl_getincom
5
6      implicit none
7
8!==================================================================
9!     
10!     Purpose
11!     -------
12!     Run the physics package of the universal model in a 1D column.
13!     
14!     It can be compiled with a command like (e.g. for 25 layers):
15!                                  "makegcm -p std -d 25 rcm1d"
16!     It requires the files "callphys.def", "z2sig.def",
17!     "traceur.def", and "run.def" with a line "INCLUDEDEF=callphys.def"
18!
19!     Authors
20!     -------
21!     Frederic Hourdin
22!     R. Fournier
23!     F. Forget
24!     F. Montmessin (water ice added)
25!     R. Wordsworth
26!     B. Charnay
27!     A. Spiga
28!
29!==================================================================
30
31#include "dimensions.h"
32#include "paramet.h"
33#include "dimphys.h"
34#include "comgeomfi.h"
35#include "surfdat.h"
36#include "comsoil.h"
37#include "comdiurn.h"
38#include "callkeys.h"
39#include "comcstfi.h"
40#include "planete.h"
41#include "comsaison.h"
42#include "control.h"
43#include "comvert.h"
44#include "netcdf.inc"
45#include "watercap.h"
46#include "fisice.h"
47#include "logic.h"
48#include "advtrac.h"
49#include "comgeom.h"
50
51c --------------------------------------------------------------
52c  Declarations
53c --------------------------------------------------------------
54c
55      INTEGER unitstart      ! unite d'ecriture de "startfi"
56      INTEGER nlayer,nlevel,nsoil,ndt
57      INTEGER ilayer,ilevel,isoil,idt,iq
58      LOGICAl firstcall,lastcall
59c
60      INTEGER day0          ! date initial (sol ; =0 a Ls=0)
61      REAL day           ! date durant le run
62      REAL time             ! time (0<time<1 ; time=0.5 a midi)
63      REAL play(nlayermx)   ! Pressure at the middle of the layers (Pa)
64      REAL plev(nlayermx+1) ! intermediate pressure levels (pa)
65      REAL psurf,tsurf     
66      REAL u(nlayermx),v(nlayermx)  ! zonal, meridional wind
67      REAL gru,grv   ! prescribed "geostrophic" background wind
68      REAL temp(nlayermx)   ! temperature at the middle of the layers
69      REAL q(nlayermx,nqmx) ! tracer mixing ratio (e.g. kg/kg)
70      REAL qsurf(nqmx)      ! tracer surface budget (e.g. kg.m-2)
71      REAL tsoil(nsoilmx)   ! subsurface soik temperature (K)
72!      REAL co2ice           ! co2ice layer (kg.m-2) !not used anymore
73      integer :: i_co2_ice=0  ! tracer index of co2 ice
74      integer :: i_h2o_ice=0  ! tracer index of h2o ice
75      integer :: i_h2o_vap=0  ! tracer index of h2o vapor
76      REAL emis             ! surface layer
77      REAL q2(nlayermx+1)   ! Turbulent Kinetic Energy
78      REAL zlay(nlayermx)   ! altitude estimee dans les couches (km)
79
80c    Physical and dynamical tandencies (e.g.  m.s-2, K/s, Pa/s)
81      REAL du(nlayermx),dv(nlayermx),dtemp(nlayermx)
82      REAL dudyn(nlayermx),dvdyn(nlayermx),dtempdyn(nlayermx)
83      REAL dpsurf   
84      REAL dq(nlayermx,nqmx)
85      REAL dqdyn(nlayermx,nqmx)
86
87c   Various intermediate variables
88!      INTEGER thermo
89      REAL zls
90      REAL phi(nlayermx),h(nlayermx),s(nlayermx)
91      REAL pks, ptif, w(nlayermx)
92      REAL qtotinit, mqtot(nqmx),qtot
93      INTEGER ierr, aslun
94      REAL tmp1(0:nlayermx),tmp2(0:nlayermx)
95      Logical  tracerdyn
96      integer :: nq=1 ! number of tracers
97 
98      character*2 str2
99      character (len=7) :: str7
100
101      logical oldcompare, earthhack,saveprofile
102
103!     added by RW for zlay computation
104      real Hscale, Hmax, rho, dz
105
106!     added by RW for autozlevs computation
107      real nu, xx, pMIN, zlev, Htop
108      real logplevs(nlayermx)
109
110!     added by BC
111      REAL cloudfrac(ngridmx,nlayermx)
112      REAL hice(ngridmx),totcloudfrac(ngridmx)
113      REAL qzero1D   !initial water amount on the ground
114
115c=======================================================================
116c INITIALISATION
117c=======================================================================
118
119      saveprofile=.false.
120      saveprofile=.true.
121
122c ----------------------------------------
123c  Default values  (corresponding to Mars)
124c ----------------------------------------
125
126      pi=2.E+0*asin(1.E+0)
127
128c     Parametres Couche limite et Turbulence
129c     --------------------------------------
130      z0 =  1.e-2                ! surface roughness (m) ~0.01
131      emin_turb = 1.e-6          ! energie minimale ~1.e-8
132      lmixmin = 30               ! longueur de melange ~100
133 
134c     propriete optiques des calottes et emissivite du sol
135c     ----------------------------------------------------
136      emissiv= 0.95              ! Emissivite du sol martien ~.95
137      emisice(1)=0.95            ! Emissivite calotte nord
138      emisice(2)=0.95            ! Emissivite calotte sud 
139
140      albedice(1)=0.5            ! Albedo calotte nord
141      albedice(2)=0.5            ! Albedo calotte sud
142      iceradius(1) = 100.e-6     ! mean scat radius of CO2 snow (north)
143      iceradius(2) = 100.e-6     ! mean scat radius of CO2 snow (south)
144      dtemisice(1) = 2.          ! time scale for snow metamorphism (north)
145      dtemisice(2) = 2.          ! time scale for snow metamorphism (south
146      hybrid=.false.
147
148c ------------------------------------------------------
149c  Load parameters from "run.def" and "gases.def"
150c ------------------------------------------------------
151
152! check if 'rcm1d.def' file is around
153      open(90,file='rcm1d.def',status='old',form='formatted',
154     &     iostat=ierr)
155      if (ierr.ne.0) then
156        write(*,*) 'Cannot find required file "rcm1d.def"'
157        write(*,*) 'which should contain some input parameters'
158        write(*,*) ' ... might as well stop here ...'
159        stop
160      else
161        close(90)
162      endif
163
164! now, run.def is needed anyway. so we create a dummy temporary one
165! for ioipsl to work. if a run.def is already here, stop the
166! program and ask the user to do a bit of cleaning
167      open(90,file='run.def',status='old',form='formatted',
168     &     iostat=ierr)
169      if (ierr.eq.0) then
170        close(90)
171        write(*,*) 'There is already a run.def file.'
172        write(*,*) 'This is not compatible with 1D runs.'
173        write(*,*) 'Please remove the file and restart the run.'
174        write(*,*) 'Runtime parameters are supposed to be in rcm1d.def'
175        stop
176      else
177        call system('touch run.def')
178        call system("echo 'INCLUDEDEF=callphys.def' >> run.def")
179        call system("echo 'INCLUDEDEF=rcm1d.def' >> run.def")
180      endif
181
182! check if we are going to run with or without tracers
183      write(*,*) "Run with or without tracer transport ?"
184      tracer=.false. ! default value
185      call getin("tracer",tracer)
186      write(*,*) " tracer = ",tracer
187
188! OK. now that run.def has been read once -- any variable is in memory.
189! so we can dump the dummy run.def
190      call system("rm -rf run.def")
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         i_h2o_vap=0
232         do iq=1,nq
233           if (tnom(iq)=="co2_ice") then
234             i_co2_ice=iq
235           elseif (tnom(iq)=="h2o_ice") then
236             i_h2o_ice=iq
237           elseif (tnom(iq)=="h2o_vap") then
238             i_h2o_vap=iq
239           endif
240         enddo
241         !if (i_co2_ice==0) then
242         !  write(*,*) "rcm1d: error, we need a 'co2_ice' tracer"
243         !  write(*,*) "   (add one to traceur.def)"
244         !  stop
245         !endif
246        else
247          write(*,*) 'Cannot find required file "traceur.def"'
248          write(*,*) ' If you want to run with tracers, I need it'
249          write(*,*) ' ... might as well stop here ...'
250          stop
251        endif
252        close(90)
253
254      else
255      ! we still need to set (dummy) tracer names for physdem1
256        nq=nqmx
257        do iq=1,nq
258          write(str7,'(a1,i2.2)')'q',iq
259          tnom(iq)=str7
260        enddo
261        ! actually, we'll need at least one "co2_ice" tracer
262        ! (for surface CO2 ice)
263        tnom(1)="co2_ice"
264        i_co2_ice=1
265      endif ! of if (tracer)
266
267!!! We have to check that check_cpp_match is F for 1D computations
268!!! We think this check is better than make a particular case for 1D in inifis or calc_cpp_mugaz
269      check_cpp_match = .false.
270      call getin("check_cpp_match",check_cpp_match)
271      if (check_cpp_match) then
272          print*,"In 1D modeling, check_cpp_match is supposed to be F"
273          print*,"Please correct callphys.def"
274          stop
275      endif
276
277!!! GEOGRAPHICAL INITIALIZATIONS
278     !!! AREA. useless in 1D
279      area(1)=1.E+0
280      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
281     !!! GEOPOTENTIAL. useless in 1D because control by surface pressure
282      phisfi(1)=0.E+0
283     !!! LATITUDE. can be set.
284      lati(1)=0 ! default value for lati(1)
285      PRINT *,'latitude (in degrees) ?'
286      call getin("latitude",lati(1))
287      write(*,*) " latitude = ",lati(1)
288      lati(1)=lati(1)*pi/180.E+0
289     !!! LONGITUDE. useless in 1D.
290      long(1)=0.E+0
291      long(1)=long(1)*pi/180.E+0
292
293!!! CALL INIFIS
294!!! - read callphys.def
295!!! - calculate sine and cosine of longitude and latitude
296!!! - calculate mugaz and cp
297!!! NB: some operations are being done dummily in inifis in 1D
298!!! - physical constants: nevermind, things are done allright below
299!!! - physical frequency: nevermind, in inifis this is a simple print
300      CALL inifis(1,llm,day0,daysec,dtphys,
301     .            lati,long,area,rad,g,r,cpp)
302!!! We check everything is OK.
303      PRINT *,"CHECK"
304      PRINT *,"--> mugaz = ",mugaz
305      PRINT *,"--> cpp = ",cpp
306      r = 8.314511E+0 * 1000.E+0 / mugaz
307      rcp = r / cpp
308
309
310!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
311!!!! PLANETARY CONSTANTS !!!!
312!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
313
314      g = -99999.
315      PRINT *,'GRAVITY in m s-2 ?'
316      call getin("g",g)
317      IF (g.eq.-99999.) THEN
318          PRINT *,"STOP. I NEED g IN RUN.DEF."
319          STOP
320      ELSE
321          PRINT *,"--> g = ",g
322      ENDIF
323
324      !rad = -99999.
325      !PRINT *,'PLANETARY RADIUS in m ?'
326      !call getin("rad",rad)
327      !IF (rad.eq.-99999.) THEN
328      !    PRINT *,"STOP. I NEED rad IN RUN.DEF."
329      !    STOP
330      !ELSE
331      !    PRINT *,"--> rad = ",rad
332      !ENDIF
333
334      daysec = -99999.
335      PRINT *,'LENGTH OF A DAY in s ?'
336      call getin("daysec",daysec)
337      IF (daysec.eq.-99999.) THEN
338          PRINT *,"STOP. I NEED daysec IN RUN.DEF."
339          STOP
340      ELSE
341          PRINT *,"--> daysec = ",daysec
342      ENDIF
343      omeg=4.*asin(1.)/(daysec)
344      PRINT *,"OK. FROM THIS I WORKED OUT:"
345      PRINT *,"--> omeg = ",omeg
346
347      year_day = -99999.
348      PRINT *,'LENGTH OF A YEAR in days ?'
349      call getin("year_day",year_day)
350      IF (year_day.eq.-99999.) THEN
351          PRINT *,"STOP. I NEED year_day IN RUN.DEF."
352          STOP
353      ELSE
354          PRINT *,"--> year_day = ",year_day
355      ENDIF
356
357      periastr = -99999.
358      PRINT *,'MIN DIST STAR-PLANET in AU ?'
359      call getin("periastr",periastr)
360      IF (periastr.eq.-99999.) THEN
361          PRINT *,"STOP. I NEED periastr IN RUN.DEF."
362          STOP
363      ELSE
364          PRINT *,"--> periastr = ",periastr
365      ENDIF
366
367      apoastr = -99999.
368      PRINT *,'MAX DIST STAR-PLANET in AU ?'
369      call getin("apoastr",apoastr)
370      IF (apoastr.eq.-99999.) THEN
371          PRINT *,"STOP. I NEED apoastr IN RUN.DEF."
372          STOP
373      ELSE
374          PRINT *,"--> apoastr = ",apoastr
375      ENDIF
376
377      peri_day = -99999.
378      PRINT *,'DATE OF PERIASTRON in days ?'
379      call getin("peri_day",peri_day)
380      IF (peri_day.eq.-99999.) THEN
381          PRINT *,"STOP. I NEED peri_day IN RUN.DEF."
382          STOP
383      ELSE IF (peri_day.gt.year_day) THEN
384          PRINT *,"STOP. peri_day.gt.year_day"
385          STOP
386      ELSE 
387          PRINT *,"--> peri_day = ", peri_day
388      ENDIF
389
390      obliquit = -99999.
391      PRINT *,'OBLIQUITY in deg ?'
392      call getin("obliquit",obliquit)
393      IF (obliquit.eq.-99999.) THEN
394          PRINT *,"STOP. I NEED obliquit IN RUN.DEF."
395          STOP
396      ELSE
397          PRINT *,"--> obliquit = ",obliquit
398      ENDIF
399
400      psurf = -99999.
401      PRINT *,'SURFACE PRESSURE in Pa ?'
402      call getin("psurf",psurf)
403      IF (psurf.eq.-99999.) THEN
404          PRINT *,"STOP. I NEED psurf IN RUN.DEF."
405          STOP
406      ELSE
407          PRINT *,"--> psurf = ",psurf
408      ENDIF
409      !! we need reference pressures
410      pa=psurf/30.
411      preff=psurf ! these values are not needed in 1D  (are you sure JL12)
412
413!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
414!!!! END PLANETARY CONSTANTS !!!!
415!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
416
417c  Date et heure locale du debut du run
418c  ------------------------------------
419c    Date (en sols depuis le solstice de printemps) du debut du run
420      day0 = 0 ! default value for day0
421      write(*,*) 'Initial date (in martian sols ; =0 at Ls=0)?'
422      call getin("day0",day0)
423      day=float(day0)
424      write(*,*) " day0 = ",day0
425c  Heure de demarrage
426      time=0 ! default value for time
427      write(*,*)'Initial local time (in hours, between 0 and 24)?'
428      call getin("time",time)
429      write(*,*)" time = ",time
430      time=time/24.E+0 ! convert time (hours) to fraction of sol
431
432
433c  Discretisation (Definition de la grille et des pas de temps)
434c  --------------
435c
436      nlayer=nlayermx
437      nlevel=nlayer+1
438      nsoil=nsoilmx
439
440      day_step=48 ! default value for day_step
441      PRINT *,'Number of time steps per sol ?'
442      call getin("day_step",day_step)
443      write(*,*) " day_step = ",day_step
444
445       
446      ecritphy=day_step ! default value for ecritphy
447      PRINT *,'Nunber of steps between writediagfi ?'
448      call getin("ecritphy",ecritphy)
449      write(*,*) " ecritphy = ",ecritphy
450
451      ndt=10 ! default value for ndt
452      PRINT *,'Number of sols to run ?'
453      call getin("ndt",ndt)
454      write(*,*) " ndt = ",ndt
455
456      ndt=ndt*day_step     
457      dtphys=daysec/day_step 
458      write(*,*)"-------------------------------------"
459      write(*,*)"-------------------------------------"
460      write(*,*)"--> Physical timestep is ",dtphys
461      write(*,*)"-------------------------------------"
462      write(*,*)"-------------------------------------"
463
464c output spectrum?
465      write(*,*)"Output spectral OLR?"
466      specOLR=.false.
467      call getin("specOLR",specOLR)
468      write(*,*)" specOLR = ",specOLR
469
470c
471c  pour le schema d'ondes de gravite
472c  ---------------------------------
473      zmea(1)=0.E+0
474      zstd(1)=0.E+0
475      zsig(1)=0.E+0
476      zgam(1)=0.E+0
477      zthe(1)=0.E+0
478
479c    Initialisation des traceurs
480c    ---------------------------
481
482      DO iq=1,nqmx
483        DO ilayer=1,nlayer
484           q(ilayer,iq) = 0.
485        ENDDO
486      ENDDO
487     
488      DO iq=1,nqmx
489        qsurf(iq) = 0.
490      ENDDO
491
492      if (water) then
493         qzero1D               = 0.0
494         qsurf(i_h2o_vap) = qzero1D
495      endif
496
497c   Initialisation pour prendre en compte les vents en 1-D
498c   ------------------------------------------------------
499      ptif=2.E+0*omeg*sinlat(1)
500 
501
502c    vent geostrophique
503      gru=10. ! default value for gru
504      PRINT *,'zonal eastward component of the geostrophic wind (m/s) ?'
505      call getin("u",gru)
506      write(*,*) " u = ",gru
507      grv=0. !default value for grv
508      PRINT *,'meridional northward component of the geostrophic',
509     &' wind (m/s) ?'
510      call getin("v",grv)
511      write(*,*) " v = ",grv
512
513c     Initialisation des vents  au premier pas de temps
514      DO ilayer=1,nlayer
515         u(ilayer)=gru
516         v(ilayer)=grv
517      ENDDO
518
519c     energie cinetique turbulente
520      DO ilevel=1,nlevel
521         q2(ilevel)=0.E+0
522      ENDDO
523
524c  emissivity / surface co2 ice ( + h2o ice??)
525c  -------------------------------------------
526      emis=emissiv ! default value for emissivity
527      PRINT *,'Emissivity of bare ground ?'
528      call getin("emis",emis)
529      write(*,*) " emis = ",emis
530      emissiv=emis ! we do this so that condense_co2 sets things to the right
531                   ! value if there is no snow
532
533      if(i_co2_ice.gt.0)then
534         qsurf(i_co2_ice)=0 ! default value for co2ice
535         print*,'Initial CO2 ice on the surface (kg.m-2)'
536         call getin("co2ice",qsurf(i_co2_ice))
537         write(*,*) " co2ice = ",qsurf(i_co2_ice)
538         IF (qsurf(i_co2_ice).ge.1.E+0) THEN
539            ! if we have some CO2 ice on the surface, change emissivity
540            if (lati(1).ge.0) then ! northern hemisphere
541              emis=emisice(1)
542            else ! southern hemisphere
543              emis=emisice(2)
544            endif
545         ENDIF
546      endif
547
548c  calcul des pressions et altitudes en utilisant les niveaux sigma
549c  ----------------------------------------------------------------
550
551c    Vertical Coordinates
552c    """"""""""""""""""""
553      hybrid=.true.
554      PRINT *,'Hybrid coordinates ?'
555      call getin("hybrid",hybrid)
556      write(*,*) " hybrid = ", hybrid
557
558
559      autozlevs=.false.
560      PRINT *,'Auto-discretise vertical levels ?'
561      call getin("autozlevs",autozlevs)
562      write(*,*) " autozlevs = ", autozlevs
563
564      pceil = psurf / 1000.0 ! Pascals
565      PRINT *,'Ceiling pressure (Pa) ?'
566      call getin("pceil",pceil)
567      write(*,*) " pceil = ", pceil
568
569! Test of incompatibility:
570! if autozlevs used, cannot have hybrid too
571
572      if (autozlevs.and.hybrid) then
573         print*,'Cannot use autozlevs and hybrid together!'
574         call abort
575      endif
576
577      if(autozlevs)then
578           
579         open(91,file="z2sig.def",form='formatted')
580         read(91,*) Hscale
581         DO ilayer=1,nlayer-2
582            read(91,*) Hmax
583         enddo
584         close(91)
585 
586           
587         print*,'Hmax = ',Hmax,' km'
588         print*,'Auto-shifting Hscale to:'
589!     Hscale = Hmax / log(psurf/100.0)
590         Hscale = Hmax / log(psurf/pceil)
591         print*,'Hscale = ',Hscale,' km'
592         
593! none of this matters if we dont care about zlay
594         
595      endif
596
597      call disvert
598
599         if(.not.autozlevs)then
600            ! we want only the scale height from z2sig, in order to compute zlay
601            open(91,file="z2sig.def",form='formatted')
602            read(91,*) Hscale
603            close(91)
604         endif
605
606!         if(autozlevs)then
607!            open(94,file="Hscale.temp",form='formatted')
608!            read(94,*) Hscale
609!            close(94)
610!         endif
611
612         DO ilevel=1,nlevel
613            plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
614         ENDDO
615         
616         DO ilayer=1,nlayer
617            play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
618         ENDDO
619         
620
621
622         DO ilayer=1,nlayer
623!     zlay(ilayer)=-300.E+0 *r*log(play(ilayer)/plev(1))
624!     &   /g
625            zlay(ilayer)=-1000.0*Hscale*log(play(ilayer)/plev(1))
626         ENDDO
627
628!      endif
629
630c  profil de temperature au premier appel
631c  --------------------------------------
632      pks=psurf**rcp
633
634c altitude en km dans profile: on divise zlay par 1000
635      tmp1(0)=0.E+0
636      DO ilayer=1,nlayer
637        tmp1(ilayer)=zlay(ilayer)/1000.E+0
638      ENDDO
639      call profile(nlayer+1,tmp1,tmp2)
640
641      tsurf=tmp2(0)
642      DO ilayer=1,nlayer
643        temp(ilayer)=tmp2(ilayer)
644      ENDDO
645     
646
647c  Initialisation albedo / inertie du sol
648c  --------------------------------------
649      albedodat(1)=0.2 ! default value for albedodat
650      PRINT *,'Albedo of bare ground ?'
651      call getin("albedo",albedodat(1))
652      write(*,*) " albedo = ",albedodat(1)
653
654      inertiedat(1,1)=400 ! default value for inertiedat
655      PRINT *,'Soil thermal inertia (SI) ?'
656      call getin("inertia",inertiedat(1,1))
657      write(*,*) " inertia = ",inertiedat(1,1)
658
659! Initialize soil properties and temperature
660! ------------------------------------------
661      volcapa=1.e6 ! volumetric heat capacity
662      DO isoil=1,nsoil
663         inertiedat(1,isoil)=inertiedat(1,1) ! soil thermal inertia
664         tsoil(isoil)=tsurf  ! soil temperature
665      ENDDO
666
667! Initialize depths
668! -----------------
669      do isoil=0,nsoil-1
670        mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) ! mid-layer depth
671      enddo
672      do isoil=1,nsoil
673        layer(isoil)=2.e-4*(2.**(isoil-1)) ! layer depth
674      enddo
675
676
677c  Ecriture de "startfi"
678c  --------------------
679c  (Ce fichier sera aussitot relu au premier
680c   appel de "physiq", mais il est necessaire pour passer
681c   les variables purement physiques a "physiq"...
682
683      call physdem1("startfi.nc",long,lati,nsoilmx,nqmx,
684     &              dtphys,float(day0),
685     &              time,tsurf,tsoil,emis,q2,qsurf,
686     &              area,albedodat,inertiedat,zmea,zstd,zsig,zgam,zthe,
687     &              cloudfrac,totcloudfrac,hice)
688
689c=======================================================================
690c  BOUCLE TEMPORELLE DU MODELE 1D
691c=======================================================================
692
693      firstcall=.true.
694      lastcall=.false.
695
696      DO idt=1,ndt
697        IF (idt.eq.ndt) then       !test
698         lastcall=.true.
699         call stellarlong(day*1.0,zls)
700!         write(103,*) 'Ls=',zls*180./pi
701!         write(103,*) 'Lat=', lati(1)*180./pi
702!         write(103,*) 'RunEnd - Atmos. Temp. File'
703!         write(103,*) 'RunEnd - Atmos. Temp. File'
704!         write(104,*) 'Ls=',zls*180./pi
705!         write(104,*) 'Lat=', lati(1)
706!         write(104,*) 'RunEnd - Atmos. Temp. File'
707        ENDIF
708
709c    calcul du geopotentiel
710c     ~~~~~~~~~~~~~~~~~~~~~
711
712
713      DO ilayer=1,nlayer
714
715!              if(autozlevs)then
716!                 s(ilayer)=(play(ilayer)/psurf)**rcp
717!              else
718          s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
719!              endif
720              !s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
721          h(ilayer)=cpp*temp(ilayer)/(pks*s(ilayer))
722       ENDDO
723
724!      DO ilayer=1,nlayer
725!        s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
726!        h(ilayer)=cpp*temp(ilayer)/(pks*s(ilayer))
727!      ENDDO
728      phi(1)=pks*h(1)*(1.E+0-s(1))
729      DO ilayer=2,nlayer
730         phi(ilayer)=phi(ilayer-1)+
731     &               pks*(h(ilayer-1)+h(ilayer))*.5E+0
732     &                  *(s(ilayer-1)-s(ilayer))
733
734      ENDDO
735
736c       appel de la physique
737c       --------------------
738
739
740      CALL physiq (1,llm,nqmx,
741     ,     firstcall,lastcall,
742     ,     day,time,dtphys,
743     ,     plev,play,phi,
744     ,     u, v,temp, q, 
745     ,     w,
746C - sorties
747     s     du, dv, dtemp, dq,dpsurf,tracerdyn)
748
749
750c       evolution du vent : modele 1D
751c       -----------------------------
752 
753c       la physique calcule les derivees temporelles de u et v.
754c       on y rajoute betement un effet Coriolis.
755c
756c       DO ilayer=1,nlayer
757c          du(ilayer)=du(ilayer)+ptif*(v(ilayer)-grv)
758c          dv(ilayer)=dv(ilayer)+ptif*(-u(ilayer)+gru)
759c       ENDDO
760
761c       Pour certain test : pas de coriolis a l'equateur
762c       if(lati(1).eq.0.) then
763          DO ilayer=1,nlayer
764             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
765             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
766          ENDDO
767c       end if
768c     
769c
770c       Calcul du temps au pas de temps suivant
771c       ---------------------------------------
772        firstcall=.false.
773        time=time+dtphys/daysec
774        IF (time.gt.1.E+0) then
775            time=time-1.E+0
776            day=day+1
777        ENDIF
778
779c       calcul des vitesses et temperature au pas de temps suivant
780c       ----------------------------------------------------------
781
782        DO ilayer=1,nlayer
783           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
784           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
785           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
786        ENDDO
787
788c       calcul des pressions au pas de temps suivant
789c       ----------------------------------------------------------
790
791           psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
792           DO ilevel=1,nlevel
793              plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
794           ENDDO
795           DO ilayer=1,nlayer
796                 play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
797           ENDDO
798
799c       calcul traceur au pas de temps suivant
800c       --------------------------------------
801
802        DO iq = 1, nqmx
803          DO ilayer=1,nlayer
804             q(ilayer,iq)=q(ilayer,iq)+dtphys*dq(ilayer,iq)
805          ENDDO
806        END DO
807
808c    ========================================================
809c    GESTION DES SORTIE
810c    ========================================================
811      if(saveprofile)then
812         OPEN(12,file='profile.out',form='formatted')
813         write(12,*) tsurf
814         DO ilayer=1,nlayermx
815            write(12,*) temp(ilayer) !, play(ilayer) !AS12 only temp so that iprofile=8 can be used
816         ENDDO
817         CLOSE(12)
818      endif
819
820
821      ENDDO   ! fin de la boucle temporelle
822
823c    ========================================================
824      end                       !rcm1d
825 
826c***********************************************************************
827c***********************************************************************
828c     Subroutines Bidons utilise seulement en 3D, mais
829c     necessaire a la compilation de rcm1d en 1D
830
831      subroutine gr_fi_dyn
832      RETURN
833      END
834 
835c***********************************************************************
836c***********************************************************************
837
838#include "../dyn3d/disvert.F"
Note: See TracBrowser for help on using the repository browser.