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

Last change on this file since 1458 was 1422, checked in by milmd, 10 years ago

In GENERIC, MARS and COMMON models replace some include files by modules (usefull for decoupling physics with dynamics).

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