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

Last change on this file since 2511 was 2410, checked in by emillour, 5 years ago

Mars GCM:
Follow-up of previous commit (adaptations for 1D and other main programs).
EM

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