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

Last change on this file since 1542 was 1542, checked in by emillour, 9 years ago

Generic GCM:

  • fix for 1D in writediagfi to enable writing at "ecritphy" rate.
  • move iniprint.h to "misc"
  • Some code cleanup in anticipation of future updates:
    • changed variable names in comgeomphy.F90: give them more explicit names: rlond => longitude , rlatd => latitude, airephy => cell_area, cuphy => dx , cvphy => dy
    • removed long(), lati() and area() from comgeomfi_h.F90, use longitude(), latitude() and cell_are() from comgeomphy.F90 instead

EM

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