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

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

18/07/2012 == JL

  • New water cycle scheme:
    • largescale now in F90. Robustness increased by i) including evap inside largescale ii) computing the

condensed water amount iteratively

  • same improvements in moistadj.
  • Water thermodynamical data and saturation curves centralized in module watercommn_h
    • The saturation curves used are now Tetens formula as they are analyticaly inversible (Ts(P)-> Ps(T)).

New saturation curve yields very good agreement with the former one.

  • Saturation curves are now generalized for arbitrary water amount (not just q<<1)
  • The old watersat should be removed soon.
  • The effect of water vapor on total (surface) pressure can be taken into account by setting

mass_redistrib=.true. in callphys.def (routine mass_redistribution inspired from co2_condense in martian
model but with a different scheme as many routines evaporate/condense water vapor).

  • New cloud and precipitation scheme (JL + BC):
    • The default recovery assumption for computing the total cloud fraction has been changed (total random gave too

large cloud fractions). See totalcloudfrac.F90 for details and to change this.

  • Totalcloudfraction now set the total cloud fraction to the fraction of the

optically thickest cloud and totalcloudfrac is thus called in aeropacity.

  • Only the total cloud fraction is used to compute optical depth in aeropacity (no more effective

optical depth with exponential formula).

  • 4 precipitation schemes are now available (see rain.F90 for details). The choice can be made using precip_scheme

in callphys.def. Usage of the more physically based model of Boucher et al 95 (precip_scheme=4) is recommended.
default behavior is set to the former "simple scheme" (precip_scheme=1).

  • See rain.f90 to determine the parameter to be defined in callphys.def as a function of the precipitation scheme used.
  • Physiq.F90 now written in a matricial (more F90) way.
  • Radii (H2O and CO2 cloud particles, aerosols, duts, ...) calculations now centralized in module radii_mod.F90

and work with the new aerosol scheme implemented by Laura K. Some inconsistency may remain in callsedim.


Implementation compiled with ifort and pgf90.
gcm.e runs in Earth and Early Mars case with CO2 and H2O cycle + dust.

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