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

Last change on this file since 1309 was 1308, checked in by emillour, 10 years ago

Generic GCM:
Some cleanup to simplify dynamics/physics interactions by getting rid
of dimphys.h (i.e. the nlayermx parameter) and minimizing use of
dimension.h in the physics.
EM

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