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

Last change on this file since 2316 was 2312, checked in by lrossi, 5 years ago

MARS GCM:
Implementation of HDO in the physics
New tracers hdo_vap and hdo_ice are added. Cf. traceur.def.isotopes in deftank for exemple of traceur.def.
All HDO related computations are activated by flag hdo=.true. in callphys.def. (see callphys.def.hdo in deftank for an example)
Additional option hdofrac (true by default) allows for turning on/off fractionation (for tests).
For now, only functional with simplified cloud scheme (so microphys=.false. and activice=.false. also recommended).
Initialisation of start files can be done with option 'inihdo' in newstart.
LR

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