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

Last change on this file since 848 was 787, checked in by aslmd, 12 years ago

LMDZ.GENERIC. (Sorry for long text but this is a quite major commit)

Paving the path for parallel computations. And moving towards a more flexible code.

Automatic allocation is used within all routines in phystd. No further mention to ngridmx and nqmx.

  1. ngridmx and nqmx are still used in LMDZ.GENERIC in the dyn3d part
  2. if the LMDZ4/LMDZ5 dynamical core is used, there is no more fixed dimensions ngridmx and nqmx --> a fully flexible parallel implementation is now possible (e.g. no need to recompile when changing numbers of processors)

The important stuff :

  • Compilation checked with ifort. OK with and without debug mode. No errors. Checked for: gcm, newstart, rcm1d, kcm1d
  • RUN GCM: Running an Earth test case. Comparison with previous revision --> debug mode : perfect match. bit by bit (diff command). checked with plots --> O1 mode : close match (checked with plots) --> O2 mode : sometimes up to 0.5 K departure.... BUT in this new version O2 and O1 are quite close while in previous version O1 and O2 differed by about, well, typically 0.5 K (pictures available on request)
  • RUN NEWSTART : perfect match (bit-by-bit) in either debug or normal mode.
  • RUN RCM1D : perfect match in normal mode.
  • RUN KCM1D : not tested (I don't know what is the use of kcm1d)

List of main changes :

  • Additional arguments to some subroutines (ngrid and nq)
  • F77 include strategy is obsolete and replaced by F90 module strategy In this new strategy arrays are allocatable and allocated once at first use This has to be done for all common featuring arrays defined with ngridmx or nqmx

surfdat.h >> surfdat_h.F90
tracer.h >> tracer_h.F90
comsaison.h >> comsaison_h.F90
comgeomfi.h >> comgeomfi_h.F90
comsoil.h >> comsoil_h.F90
comdiurn.h >> comdiurn_h.F90
fisice.h >> DELETED. was not used. probably a fossil.
watercap.h >> DELETED. variable put in surfdat_h.F90

  • F77 'save' strategy is obsolete and replaced by F90 'allocatable save' strategy (see previous point and e.g. new version of physiq.F90)
  • Suppressing any mention to advtrac.h which is a common in the dynamics and needs nqmx This was easily solved by adding an argument with tracer names, coming from the dynamics This is probably not a definitive solution, ... but this allows for generic physics to work easily with either LMDZ.GENERIC or LMDZ dynamical cores
  • Removing consistency tests between nq and nqmx ; and ngrid and ngridmx. No use now!
  • Adaptation of rcm1d, kcm1d, newstart given above-mentioned changes

A note on phyetat0 and soil_setting:

  • Now written so that a slice of horizontal size 'ngrid' starting at grid point 'cursor' is read in startfi.nc 'cursor' is defined in dimphys.h and initialized by inifis (or in newstart) this is useful for parallel computations. default behavior is the usual one : sequential runs, cursor is 1, size ngrid is the whole global domain

A note on an additional change :

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