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

Last change on this file since 1576 was 1576, checked in by emillour, 8 years ago

All models:
Clean up the dynamics/physics interface to converge with LMDZ5;
get rid of tracerdyn parameter (which was supposed to send information
about wether tracers should be advected or not in the Mars and Generic
models).
EM

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