source: trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F @ 2387

Last change on this file since 2387 was 2355, checked in by cmathe, 6 years ago

MARS GCM: improvement in the reading of input profile for 1D simulation

File size: 29.2 KB
Line 
1
2      PROGRAM testphys1d
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, nqperes,nqfils
8      use comsoil_h, only: volcapa, layer, mlayer, inertiedat, nsoilmx
9      use comgeomfi_h, only: sinlat, ini_fillgeom
10      use surfdat_h, only: albedodat, z0_default, emissiv, emisice,
11     &                     albedice, iceradius, dtemisice, z0,
12     &                     zmea, zstd, zsig, zgam, zthe, phisfi,
13     &                     watercaptag, watercap, hmons, summit, base
14      use slope_mod, only: theta_sl, psi_sl
15      use phyredem, only: physdem0,physdem1
16      use geometry_mod, only: init_geometry
17      use planete_h, only: year_day, periheli, aphelie, peri_day,
18     &                     obliquit, emin_turb, lmixmin
19      use comcstfi_h, only: pi, rad, omeg, g, mugaz, rcp, r, cpp
20      use time_phylmdz_mod, only: daysec, dtphys, day_step,
21     &                            ecritphy, iphysiq
22      use dimradmars_mod, only: tauscaling,tauvis,totcloudfrac
23      use co2cloud_mod, only: mem_Mccn_co2, mem_Mh2o_co2,
24     &                        mem_Nccn_co2
25      USE comvert_mod, ONLY: ap,bp,aps,bps,pa,preff,sig,
26     &                       presnivs,pseudoalt,scaleheight
27      USE vertical_layers_mod, ONLY: init_vertical_layers
28      USE logic_mod, ONLY: hybrid
29      use physics_distribution_mod, only: init_physics_distribution
30      use regular_lonlat_mod, only: init_regular_lonlat
31      use mod_interface_dyn_phys, only: init_interface_dyn_phys
32      USE phys_state_var_init_mod, ONLY: phys_state_var_init
33      USE physiq_mod, ONLY: physiq
34      USE read_profile_mod, ONLY: read_profile
35      IMPLICIT NONE
36
37c=======================================================================
38c   subject:
39c   --------
40c   PROGRAM useful to run physical part of the martian GCM in a 1D column
41c       
42c Can be compiled with a command like (e.g. for 25 layers)
43c  "makegcm -p mars -d 25 testphys1d"
44c It requires the files "testphys1d.def" "callphys.def"
45c   and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line)
46c      and a file describing the sigma layers (e.g. "z2sig.def")
47c
48c   author: Frederic Hourdin, R.Fournier,F.Forget
49c   -------
50c   
51c   update: 12/06/2003 including chemistry (S. Lebonnois)
52c                            and water ice (F. Montmessin)
53c
54c=======================================================================
55
56      include "dimensions.h"
57      integer, parameter :: ngrid = 1 !(2+(jjm-1)*iim - 1/jjm)
58      integer, parameter :: nlayer = llm
59!#include "dimradmars.h"
60!#include "comgeomfi.h"
61!#include "surfdat.h"
62!#include "slope.h"
63!#include "comsoil.h"
64!#include "comdiurn.h"
65      include "callkeys.h"
66!#include "comsaison.h"
67!#include "control.h"
68      include "netcdf.inc"
69      include "comg1d.h"
70!#include "advtrac.h"
71
72c --------------------------------------------------------------
73c  Declarations
74c --------------------------------------------------------------
75c
76      INTEGER unitstart      ! unite d'ecriture de "startfi"
77      INTEGER nlevel,nsoil,ndt
78      INTEGER ilayer,ilevel,isoil,idt,iq
79      LOGICAl firstcall,lastcall
80c
81      real,parameter :: odpref=610. ! DOD reference pressure (Pa)
82c
83      INTEGER day0          ! date initial (sol ; =0 a Ls=0)
84      REAL day           ! date durant le run
85      REAL time             ! time (0<time<1 ; time=0.5 a midi)
86      REAL play(nlayer)   ! Pressure at the middle of the layers (Pa)
87      REAL plev(nlayer+1) ! intermediate pressure levels (pa)
88      REAL psurf,tsurf(1)     
89      REAL u(nlayer),v(nlayer)  ! zonal, meridional wind
90      REAL gru,grv   ! prescribed "geostrophic" background wind
91      REAL temp(nlayer)   ! temperature at the middle of the layers
92      REAL,ALLOCATABLE :: q(:,:) ! tracer mixing ratio (e.g. kg/kg)
93      REAL,ALLOCATABLE :: qsurf(:) ! tracer surface budget (e.g. kg.m-2)
94      REAL tsoil(nsoilmx)   ! subsurface soik temperature (K)
95      REAL co2ice(1)        ! co2ice layer (kg.m-2)
96      REAL emis(1)          ! surface layer
97      REAL albedo(1,1)      ! surface albedo
98      REAL :: wstar(1)=0.    ! Thermals vertical velocity
99      REAL q2(nlayer+1)   ! Turbulent Kinetic Energy
100      REAL zlay(nlayer)   ! altitude estimee dans les couches (km)
101
102c    Physical and dynamical tandencies (e.g.  m.s-2, K/s, Pa/s)
103      REAL du(nlayer),dv(nlayer),dtemp(nlayer)
104      REAL dudyn(nlayer),dvdyn(nlayer),dtempdyn(nlayer)
105      REAL dpsurf(1)   
106      REAL,ALLOCATABLE :: dq(:,:)
107      REAL,ALLOCATABLE :: dqdyn(:,:)
108
109c   Various intermediate variables
110      INTEGER thermo
111      REAL zls
112      REAL phi(nlayer),h(nlayer),s(nlayer)
113      REAL pks, ptif, w(nlayer)
114      REAL qtotinit,qtot
115      real,allocatable :: mqtot(:)
116      INTEGER ierr, aslun
117      REAL tmp1(0:nlayer),tmp2(0:nlayer)
118      integer :: nq=1 ! number of tracers
119      real :: latitude(1), longitude(1), cell_area(1)
120
121      character*2 str2
122      character (len=7) :: str7
123      character(len=44) :: txt
124
125c   New flag to compute paleo orbital configurations + few variables JN
126      REAL halfaxe, excentric, Lsperi
127      Logical paleomars
128
129c   MVals: isotopes as in the dynamics (CRisi)
130      INTEGER :: ifils,ipere,generation
131      CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_transp ! transporting fluid short name
132      CHARACTER(len=80) :: line ! to store a line of text     
133      INTEGER ierr0
134      LOGICAL :: continu
135
136c=======================================================================
137
138c=======================================================================
139c INITIALISATION
140c=======================================================================
141! initialize "serial/parallel" related stuff
142!      CALL init_phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
143!      CALL init_phys_lmdz(1,1,llm,1,(/1/))
144!      call initcomgeomphy
145
146c ------------------------------------------------------
147c  Prescribed constants to be set here
148c ------------------------------------------------------
149
150      pi=2.E+0*asin(1.E+0)
151
152c     Mars planetary constants
153c     ----------------------------
154      rad=3397200.               ! mars radius (m)  ~3397200 m
155      daysec=88775.              ! length of a sol (s)  ~88775 s
156      omeg=4.*asin(1.)/(daysec)  ! rotation rate (rad.s-1)
157      g=3.72                     ! gravity (m.s-2) ~3.72 
158      mugaz=43.49                ! atmosphere mola mass (g.mol-1) ~43.49
159      rcp=.256793                ! = r/cp  ~0.256793
160      r= 8.314511E+0 *1000.E+0/mugaz
161      cpp= r/rcp
162      year_day = 669             ! lenght of year (sols) ~668.6
163      periheli = 206.66          ! minimum sun-mars distance (Mkm) ~206.66
164      aphelie = 249.22           ! maximum sun-mars distance (Mkm) ~249.22
165      halfaxe = 227.94           ! demi-grand axe de l'ellipse
166      peri_day =  485.           ! perihelion date (sols since N. Spring)
167      obliquit = 25.2            ! Obliquity (deg) ~25.2         
168      excentric = 0.0934         ! Eccentricity (0.0934)         
169 
170c     Planetary Boundary Layer and Turbulence parameters
171c     --------------------------------------------------
172      z0_default =  1.e-2        ! surface roughness (m) ~0.01
173      emin_turb = 1.e-6          ! minimal turbulent energy ~1.e-8
174      lmixmin = 30               ! mixing length ~100
175 
176c     cap properties and surface emissivities
177c     ----------------------------------------------------
178      emissiv= 0.95              ! Bare ground emissivity ~.95
179      emisice(1)=0.95            ! Northern cap emissivity
180      emisice(2)=0.95            ! Southern cap emisssivity
181      albedice(1)=0.5            ! Northern cap albedo
182      albedice(2)=0.5            ! Southern cap albedo
183      iceradius(1) = 100.e-6     ! mean scat radius of CO2 snow (north)
184      iceradius(2) = 100.e-6     ! mean scat radius of CO2 snow (south)
185      dtemisice(1) = 2.          ! time scale for snow metamorphism (north)
186      dtemisice(2) = 2.          ! time scale for snow metamorphism (south
187
188c     mesh surface (not a very usefull quantity in 1D)
189c     ----------------------------------------------------
190      cell_area(1)=1.E+0
191     
192c ------------------------------------------------------
193c  Loading run parameters from "run.def" file
194c ------------------------------------------------------
195
196
197! check if 'run.def' file is around (otherwise reading parameters
198! from callphys.def via getin() routine won't work.
199      open(99,file='run.def',status='old',form='formatted',
200     &     iostat=ierr)
201      if (ierr.ne.0) then
202        write(*,*) 'Cannot find required file "run.def"'
203        write(*,*) '  (which should contain some input parameters'
204        write(*,*) '   along with the following line:'
205        write(*,*) '   INCLUDEDEF=callphys.def'
206        write(*,*) '   )'
207        write(*,*) ' ... might as well stop here ...'
208        stop
209      else
210        close(99)
211      endif
212
213! check if we are going to run with or without tracers
214      write(*,*) "Run with or without tracer transport ?"
215      tracer=.false. ! default value
216      call getin("tracer",tracer)
217      write(*,*) " tracer = ",tracer
218
219! while we're at it, check if there is a 'traceur.def' file
220! and process it.
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.ne.0) then
226          write(*,*) 'Cannot find required file "traceur.def"'
227          write(*,*) ' If you want to run with tracers, I need it'
228          write(*,*) ' ... might as well stop here ...'
229          stop
230        else
231          write(*,*) "testphys1d: Reading file traceur.def"
232          ! read number of tracers:
233          read(90,*,iostat=ierr) nq
234          nqtot=nq ! set value of nqtot (in infotrac module) as nq
235          if (ierr.ne.0) then
236            write(*,*) "testphys1d: error reading number of tracers"
237            write(*,*) "   (first line of traceur.def) "
238            stop
239          endif
240        endif
241        ! allocate arrays:
242        allocate(tname(nq))
243        allocate(q(nlayer,nq))
244        allocate(qsurf(nq))
245        allocate(dq(nlayer,nq))
246        allocate(dqdyn(nlayer,nq))
247        allocate(mqtot(nq))
248        allocate(tnom_transp(nq))
249       
250        ! read tracer names from file traceur.def
251        do iq=1,nq
252          read(90,'(80a)',iostat=ierr) line ! store the line from traceur.def
253          if (ierr.ne.0) then
254            write(*,*) 'testphys1d: error reading tracer names...'
255            stop
256          endif
257          ! if format is tnom_0, tnom_transp (isotopes)
258          read(line,*,iostat=ierr0) tname(iq),tnom_transp(iq)
259          if (ierr0.ne.0) then
260            read(line,*) tname(iq)
261            tnom_transp(iq)='air'
262          endif
263
264        enddo
265        close(90)
266
267       ! Isotopes: as in the 3D case we have to determine father/son relations for isotopes and carrying fluid
268       ALLOCATE(nqfils(nqtot))
269       nqperes=0
270       nqfils(:)=0 
271       DO iq=1,nqtot
272       if (tnom_transp(iq) == 'air') then
273         ! ceci est un traceur père
274         WRITE(*,*) 'Le traceur',iq,', appele ',
275     &          trim(tname(iq)),', est un pere'
276         nqperes=nqperes+1
277       else !if (tnom_transp(iq) == 'air') then
278         ! ceci est un fils. Qui est son père?
279         WRITE(*,*) 'Le traceur',iq,', appele ',
280     &                trim(tname(iq)),', est un fils'
281         continu=.true.
282         ipere=1
283         do while (continu)           
284           if (tnom_transp(iq) .eq. tname(ipere)) then
285             ! Son père est ipere
286             WRITE(*,*) 'Le traceur',iq,'appele ',
287     &   trim(tname(iq)),' est le fils de ',
288     &   ipere,'appele ',trim(tname(ipere))
289             nqfils(ipere)=nqfils(ipere)+1         
290             continu=.false.
291           else !if (tnom_transp(iq) == tnom_0(ipere)) then
292             ipere=ipere+1
293             if (ipere.gt.nqtot) then
294                 WRITE(*,*) 'Le traceur',iq,'appele ',
295     &           trim(tname(iq)),', est orpelin.'
296                 CALL abort_gcm('infotrac_init',
297     &                  'Un traceur est orphelin',1)
298             endif !if (ipere.gt.nqtot) then
299           endif !if (tnom_transp(iq) == tnom_0(ipere)) then
300         enddo !do while (continu)
301       endif !if (tnom_transp(iq) == 'air') then
302       enddo !DO iq=1,nqtot
303       WRITE(*,*) 'nqperes=',nqperes   
304       WRITE(*,*) 'nqfils=',nqfils
305
306        ! initialize tracers here:
307        write(*,*) "testphys1d: initializing tracers"
308        call read_profile(nq, nlayer, qsurf, q)
309      else
310      ! we still need to set (dummy) tracer number and names for physdem1
311        nq=1
312        nqtot=nq ! set value of nqtot (in infotrac module) as nq
313        ! allocate arrays:
314        allocate(tname(nq))
315        allocate(q(nlayer,nq))
316        allocate(qsurf(nq))
317        allocate(dq(nlayer,nq))
318        allocate(dqdyn(nlayer,nq))
319        allocate(mqtot(nq))
320        do iq=1,nq
321          write(str7,'(a1,i2.2)')'t',iq
322          tname(iq)=str7
323        enddo
324      ! and just to be clean, also initialize tracers to zero for physdem1
325        q(:,:)=0
326        qsurf(:)=0     
327      endif ! of if (tracer)
328     
329      !write(*,*) "testphys1d q", q(1,:)
330      !write(*,*) "testphys1d qsurf", qsurf
331
332c  Date and local time at beginning of run
333c  ---------------------------------------
334c    Date (in sols since spring solstice) at beginning of run
335      day0 = 0 ! default value for day0
336      write(*,*) 'Initial date (in martian sols ; =0 at Ls=0)?'
337      call getin("day0",day0)
338      day=float(day0)
339      write(*,*) " day0 = ",day0
340c  Local time at beginning of run
341      time=0 ! default value for time
342      write(*,*)'Initial local time (in hours, between 0 and 24)?'
343      call getin("time",time)
344      write(*,*)" time = ",time
345      time=time/24.E+0 ! convert time (hours) to fraction of sol
346
347c  Discretization (Definition of grid and time steps)
348c  --------------
349c
350      nlevel=nlayer+1
351      nsoil=nsoilmx
352
353      day_step=48 ! default value for day_step
354      PRINT *,'Number of time steps per sol ?'
355      call getin("day_step",day_step)
356      write(*,*) " day_step = ",day_step
357
358      ecritphy=day_step ! default value for ecritphy, output every time step
359
360      ndt=10 ! default value for ndt
361      PRINT *,'Number of sols to run ?'
362      call getin("ndt",ndt)
363      write(*,*) " ndt = ",ndt
364
365      ndt=ndt*day_step     
366      dtphys=daysec/day_step 
367
368c Imposed surface pressure
369c ------------------------------------
370c
371      psurf=610. ! default value for psurf
372      PRINT *,'Surface pressure (Pa) ?'
373      call getin("psurf",psurf)
374      write(*,*) " psurf = ",psurf
375c Reference pressures
376      pa=20.   ! transition pressure (for hybrid coord.)
377      preff=610.      ! reference surface pressure
378 
379c Aerosol properties
380c --------------------------------
381      tauvis=0.2 ! default value for tauvis (dust opacity)
382      write(*,'("Reference dust opacity at ",f4.0," Pa ?")')odpref
383      call getin("tauvis",tauvis)
384      write(*,*) " tauvis = ",tauvis
385
386c Orbital parameters
387c ------------------
388      paleomars=.false. ! Default: no water ice reservoir
389      call getin("paleomars",paleomars)
390      if (paleomars.eqv..true.) then
391        write(*,*) "paleomars=", paleomars
392        write(*,*) "Orbital parameters from callphys.def"
393        write(*,*) "Enter eccentricity & Lsperi"
394        print *, 'Martian eccentricity (0<e<1) ?'
395        call getin('excentric ',excentric)
396        write(*,*)"excentric =",excentric
397        print *, 'Solar longitude of perihelion (0<Ls<360) ?'
398        call getin('Lsperi',Lsperi )
399        write(*,*)"Lsperi=",Lsperi
400        Lsperi = Lsperi*pi/180.0 ! Put it in rad for peri_day
401        periheli = halfaxe*(1-excentric)
402        aphelie = halfaxe*(1+excentric)
403        call call_dayperi(Lsperi,excentric,peri_day,year_day)
404        write(*,*) "Corresponding orbital params for GCM"
405        write(*,*) " periheli = ",periheli
406        write(*,*) " aphelie = ",aphelie
407        write(*,*) "date of perihelion (sol)",peri_day
408      else
409        write(*,*) "paleomars=", paleomars
410        write(*,*) "Default present-day orbital parameters"
411        write(*,*) "Unless specified otherwise"
412        print *,'Min. distance Sun-Mars (Mkm)?'
413        call getin("periheli",periheli)
414        write(*,*) " periheli = ",periheli
415
416        print *,'Max. distance Sun-Mars (Mkm)?'
417        call getin("aphelie",aphelie)
418        write(*,*) " aphelie = ",aphelie
419
420        print *,'Day of perihelion?'
421        call getin("periday",peri_day)
422        write(*,*) " periday = ",peri_day
423
424        print *,'Obliquity?'
425        call getin("obliquit",obliquit)
426        write(*,*) " obliquit = ",obliquit
427      end if
428 
429c  latitude/longitude
430c  ------------------
431      latitude(1)=0 ! default value for latitude
432      PRINT *,'latitude (in degrees) ?'
433      call getin("latitude",latitude(1))
434      write(*,*) " latitude = ",latitude
435      latitude=latitude*pi/180.E+0
436      longitude=0.E+0
437      longitude=longitude*pi/180.E+0
438
439!  some initializations (some of which have already been
440!  done above!) and loads parameters set in callphys.def
441!  and allocates some arrays
442!Mars possible matter with dtphys in input and include!!!
443! Initializations below should mimick what is done in iniphysiq for 3D GCM
444      call init_physics_distribution(regular_lonlat,4,
445     &                               1,1,1,nlayer,1)
446      call init_interface_dyn_phys
447      call init_regular_lonlat(1,1,longitude,latitude,
448     &                   (/0.,0./),(/0.,0./))
449      call init_geometry(1,longitude,latitude,
450     &                   (/0.,0.,0.,0./),(/0.,0.,0.,0./),
451     &                   cell_area)
452! Ehouarn: init_vertial_layers called later (because disvert not called yet)
453!      call init_vertical_layers(nlayer,preff,scaleheight,
454!     &                      ap,bp,aps,bps,presnivs,pseudoalt)
455      call init_dimphy(1,nlayer) ! Initialize dimphy module
456      call phys_state_var_init(1,llm,nq,tname,
457     .          day0,time,daysec,dtphys,rad,g,r,cpp,
458     .          nqperes,nqfils)! MVals: variables isotopes
459      call ini_fillgeom(1,latitude,longitude,(/1.0/))
460      call conf_phys(1,llm,nq)
461
462      ! in 1D model physics are called every time step
463      ! ovverride iphysiq value that has been set by conf_phys
464      if (iphysiq/=1) then
465        write(*,*) "testphys1d: setting iphysiq=1"
466        iphysiq=1
467      endif
468
469c  Initialize albedo / soil thermal inertia
470c  ----------------------------------------
471c
472      albedodat(1)=0.2 ! default value for albedodat
473      PRINT *,'Albedo of bare ground ?'
474      call getin("albedo",albedodat(1))
475      write(*,*) " albedo = ",albedodat(1)
476      albedo(1,1)=albedodat(1)
477
478      inertiedat(1,1)=400 ! default value for inertiedat
479      PRINT *,'Soil thermal inertia (SI) ?'
480      call getin("inertia",inertiedat(1,1))
481      write(*,*) " inertia = ",inertiedat(1,1)
482
483      z0(1)=z0_default ! default value for roughness
484      write(*,*) 'Surface roughness length z0 (m)?'
485      call getin("z0",z0(1))
486      write(*,*) " z0 = ",z0(1)
487
488! Initialize local slope parameters (only matters if "callslope"
489! is .true. in callphys.def)
490      ! slope inclination angle (deg) 0: horizontal, 90: vertical
491      theta_sl(1)=0.0 ! default: no inclination
492      call getin("slope_inclination",theta_sl(1))
493      ! slope orientation (deg)
494      ! 0 == Northward, 90 == Eastward, 180 == Southward, 270 == Westward
495      psi_sl(1)=0.0 ! default value
496      call getin("slope_orientation",psi_sl(1))
497     
498c
499c  for the gravity wave scheme
500c  ---------------------------------
501c
502      zmea(1)=0.E+0
503      zstd(1)=0.E+0
504      zsig(1)=0.E+0
505      zgam(1)=0.E+0
506      zthe(1)=0.E+0
507c
508c  for the slope wind scheme
509c  ---------------------------------
510
511      hmons(1)=0.E+0
512      PRINT *,'hmons is initialized to ',hmons(1)
513      summit(1)=0.E+0
514      PRINT *,'summit is initialized to ',summit(1)
515      base(1)=0.E+0
516c
517c  Default values initializing the coefficients calculated later
518c  ---------------------------------
519
520      tauscaling(1)=1. ! calculated in aeropacity_mod.F
521      totcloudfrac(1)=1. ! calculated in watercloud_mod.F     
522
523c   Specific initializations for "physiq"
524c   -------------------------------------
525c   surface geopotential is not used (or useful) since in 1D
526c   everything is controled by surface pressure
527      phisfi(1)=0.E+0
528
529c   Initialization to take into account prescribed winds
530c   ------------------------------------------------------
531      ptif=2.E+0*omeg*sinlat(1)
532 
533c    geostrophic wind
534      gru=10. ! default value for gru
535      PRINT *,'zonal eastward component of the geostrophic wind (m/s) ?'
536      call getin("u",gru)
537      write(*,*) " u = ",gru
538      grv=0. !default value for grv
539      PRINT *,'meridional northward component of the geostrophic',
540     &' wind (m/s) ?'
541      call getin("v",grv)
542      write(*,*) " v = ",grv
543
544c     Initialize winds  for first time step
545      DO ilayer=1,nlayer
546         u(ilayer)=gru
547         v(ilayer)=grv
548         w(ilayer)=0 ! default: no vertical wind
549      ENDDO
550
551c     Initialize turbulente kinetic energy
552      DO ilevel=1,nlevel
553         q2(ilevel)=0.E+0
554      ENDDO
555
556c  CO2 ice on the surface
557c  -------------------
558      co2ice(1)=0.E+0 ! default value for co2ice
559      PRINT *,'Initial CO2 ice on the surface (kg.m-2)'
560      call getin("co2ice",co2ice)
561      write(*,*) " co2ice = ",co2ice
562! Initialization for CO2 clouds (could be improved to read initial profiles)
563      mem_Mccn_co2(:,:)=0
564      mem_Mh2o_co2(:,:)=0
565      mem_Nccn_co2(:,:)=0
566c
567c  emissivity
568c  ----------
569      emis=emissiv
570      IF (co2ice(1).eq.1.E+0) THEN
571         emis=emisice(1) ! northern hemisphere
572         IF(latitude(1).LT.0) emis=emisice(2) ! southern hemisphere
573      ENDIF
574
575 
576
577c  Compute pressures and altitudes of atmospheric levels
578c  ----------------------------------------------------------------
579
580c    Vertical Coordinates
581c    """"""""""""""""""""
582      hybrid=.true.
583      PRINT *,'Hybrid coordinates ?'
584      call getin("hybrid",hybrid)
585      write(*,*) " hybrid = ", hybrid
586
587      CALL  disvert_noterre
588      ! now that disvert has been called, initialize module vertical_layers_mod
589      call init_vertical_layers(nlayer,preff,scaleheight,
590     &                      ap,bp,aps,bps,presnivs,pseudoalt)
591
592      DO ilevel=1,nlevel
593        plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
594      ENDDO
595
596      DO ilayer=1,nlayer
597        play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
598      ENDDO
599
600      DO ilayer=1,nlayer
601        zlay(ilayer)=-200.E+0 *r*log(play(ilayer)/plev(1))
602     &   /g
603      ENDDO
604
605
606c  Initialize temperature profile
607c  --------------------------------------
608      pks=psurf**rcp
609
610c altitude in km in profile: divide zlay by 1000
611      tmp1(0)=0.E+0
612      DO ilayer=1,nlayer
613        tmp1(ilayer)=zlay(ilayer)/1000.E+0
614      ENDDO
615
616      call profile(nlayer+1,tmp1,tmp2)
617
618      tsurf=tmp2(0)
619      DO ilayer=1,nlayer
620        temp(ilayer)=tmp2(ilayer)
621      ENDDO
622     
623
624
625! Initialize soil properties and temperature
626! ------------------------------------------
627      volcapa=1.e6 ! volumetric heat capacity
628      DO isoil=1,nsoil
629         inertiedat(1,isoil)=inertiedat(1,1) ! soil thermal inertia
630         tsoil(isoil)=tsurf(1)  ! soil temperature
631      ENDDO
632
633! Initialize depths
634! -----------------
635      do isoil=0,nsoil-1
636        mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) ! mid-layer depth
637      enddo
638      do isoil=1,nsoil
639        layer(isoil)=2.e-4*(2.**(isoil-1)) ! layer depth
640      enddo
641
642c    Initialize traceurs
643c    ---------------------------
644
645      if (photochem.or.callthermos) then
646         write(*,*) 'Initializing chemical species'
647         ! thermo=0: initialize over all atmospheric layers
648         thermo=0
649         call inichim_newstart(q,psurf,sig,nq,latitude,longitude,
650     $   cell_area,thermo,qsurf)
651      endif
652
653c Check if the surface is a water ice reservoir
654c --------------------------------------------------
655      watercap(1)=0 ! Initialize watercap
656      watercaptag(1)=.false. ! Default: no water ice reservoir
657      print *,'Water ice cap on ground ?'
658      call getin("watercaptag",watercaptag)
659      write(*,*) " watercaptag = ",watercaptag
660     
661
662c    Initialization for GRADS outputs in "g1d.dat" and "g1d.ctl"
663c    ----------------------------------------------------------------
664c    (output done in "writeg1d", typically called by "physiq.F")
665
666        g1d_nlayer=nlayer
667        g1d_nomfich='g1d.dat'
668        g1d_unitfich=40
669        g1d_nomctl='g1d.ctl'
670        g1d_unitctl=41
671        g1d_premier=.true.
672        g2d_premier=.true.
673
674c  Write a "startfi" file
675c  --------------------
676c  This file will be read during the first call to "physiq".
677c  It is needed to transfert physics variables to "physiq"...
678
679      call physdem0("startfi.nc",longitude,latitude,nsoilmx,ngrid,llm,
680     &              nq,dtphys,float(day0),0.,cell_area,
681     &              albedodat,inertiedat,zmea,zstd,zsig,zgam,zthe,
682     &              hmons,summit,base)
683      call physdem1("startfi.nc",nsoilmx,ngrid,llm,nq,
684     &              dtphys,time,
685     &              tsurf,tsoil,co2ice,albedo,emis,q2,qsurf,tauscaling,
686     &              totcloudfrac,wstar,
687     &              mem_Mccn_co2,mem_Nccn_co2,
688     &              mem_Mh2o_co2,watercap)
689
690c=======================================================================
691c  1D MODEL TIME STEPPING LOOP
692c=======================================================================
693c
694      firstcall=.true.
695      lastcall=.false.
696
697      DO idt=1,ndt
698c        IF (idt.eq.ndt) lastcall=.true.
699        IF (idt.eq.ndt-day_step-1) then       !test
700         lastcall=.true.
701         call solarlong(day*1.0,zls)
702         write(103,*) 'Ls=',zls*180./pi
703         write(103,*) 'Lat=', latitude(1)*180./pi
704         write(103,*) 'Tau=', tauvis/odpref*psurf
705         write(103,*) 'RunEnd - Atmos. Temp. File'
706         write(103,*) 'RunEnd - Atmos. Temp. File'
707         write(104,*) 'Ls=',zls*180./pi
708         write(104,*) 'Lat=', latitude(1)
709         write(104,*) 'Tau=', tauvis/odpref*psurf
710         write(104,*) 'RunEnd - Atmos. Temp. File'
711        ENDIF
712
713c     compute geopotential
714c     ~~~~~~~~~~~~~~~~~~~~~
715      DO ilayer=1,nlayer
716        s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
717        h(ilayer)=cpp*temp(ilayer)/(pks*s(ilayer))
718      ENDDO
719      phi(1)=pks*h(1)*(1.E+0-s(1))
720      DO ilayer=2,nlayer
721         phi(ilayer)=phi(ilayer-1)+
722     &               pks*(h(ilayer-1)+h(ilayer))*.5E+0
723     &                  *(s(ilayer-1)-s(ilayer))
724
725      ENDDO
726
727c       call physics
728c       --------------------
729!      write(*,*) "testphys1d avant q", q(1,:)
730      CALL physiq (1,llm,nq,
731     ,     firstcall,lastcall,
732     ,     day,time,dtphys,
733     ,     plev,play,phi,
734     ,     u, v,temp, q, 
735     ,     w,
736C - outputs
737     s     du, dv, dtemp, dq,dpsurf)
738!      write(*,*) "testphys1d apres q", q(1,:)
739
740
741c       wind increment : specific for 1D
742c       --------------------------------
743 
744c       The physics compute the tendencies on u and v,
745c       here we just add Coriolos effect
746c
747c       DO ilayer=1,nlayer
748c          du(ilayer)=du(ilayer)+ptif*(v(ilayer)-grv)
749c          dv(ilayer)=dv(ilayer)+ptif*(-u(ilayer)+gru)
750c       ENDDO
751
752c       For some tests : No coriolis force at equator
753c       if(latitude(1).eq.0.) then
754          DO ilayer=1,nlayer
755             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
756             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
757          ENDDO
758c       end if
759c     
760c
761c       Compute time for next time step
762c       ---------------------------------------
763        firstcall=.false.
764        time=time+dtphys/daysec
765        IF (time.gt.1.E+0) then
766            time=time-1.E+0
767            day=day+1
768        ENDIF
769
770c       compute winds and temperature for next time step
771c       ----------------------------------------------------------
772
773        DO ilayer=1,nlayer
774           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
775           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
776           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
777        ENDDO
778
779c       compute pressure for next time step
780c       ----------------------------------------------------------
781
782           psurf=psurf+dtphys*dpsurf(1)   ! surface pressure change
783           DO ilevel=1,nlevel
784             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
785           ENDDO
786           DO ilayer=1,nlayer
787             play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
788           ENDDO
789
790!       increment tracers
791        DO iq = 1, nq
792          DO ilayer=1,nlayer
793             q(ilayer,iq)=q(ilayer,iq)+dtphys*dq(ilayer,iq)
794          ENDDO
795        ENDDO
796
797      ENDDO   ! of idt=1,ndt ! end of time stepping loop
798
799c    ========================================================
800c    OUTPUTS
801c    ========================================================
802
803c    finalize and close grads files "g1d.dat" and "g1d.ctl"
804
805c        CALL endg1d(1,nlayer,zphi/(g*1000.),ndt)
806        CALL endg1d(1,nlayer,zlay/1000.,ndt)
807
808      write(*,*) "testphys1d: Everything is cool."
809
810      END
811 
812c***********************************************************************
813c***********************************************************************
814c     Dummy subroutines used only in 3D, but required to
815c     compile testphys1d (to cleanly use writediagfi)
816
817      subroutine gr_fi_dyn(nfield,ngrid,im,jm,pfi,pdyn)
818
819      IMPLICIT NONE
820
821      INTEGER im,jm,ngrid,nfield
822      REAL pdyn(im,jm,nfield)
823      REAL pfi(ngrid,nfield)
824     
825      if (ngrid.ne.1) then
826        write(*,*) "gr_fi_dyn error: in 1D ngrid should be 1!!!"
827        stop
828      endif
829     
830      pdyn(1,1,1:nfield)=pfi(1,1:nfield)
831     
832      end
833 
834c***********************************************************************
835c***********************************************************************
836
Note: See TracBrowser for help on using the repository browser.