source: trunk/mars/libf/phymars/testphys1d.F @ 131

Last change on this file since 131 was 86, checked in by aslmd, 14 years ago

*
mars + LMD_MM_MARS
* Precompilation flag MESOSCALE for better transparency

* in shared phymars between GCM and mesoscale model

*

M 85 mars/libf/phymars/meso_physiq.F
M 85 mars/libf/phymars/meso_inifis.F
Added a pre-compilation flag MESOSCALE so that the LMDZ.MARS GCM
will compile without stating errors because of mesoscale routines.

M 85 mars/libf/phymars/newcondens.F
M 85 mars/libf/phymars/testphys1d.F
M 85 mars/libf/phymars/dustlift.F
D 85 mars/libf/phymars/meso_testphys1d.F
D 85 mars/libf/phymars/meso_dustlift.F
D 85 mars/libf/phymars/meso_newcondens.F
Now, this MESOSCALE precompilation flag can be used to lower
the number of meso_* routines when adaptations for mesoscale
applications are not very extended.
--> Three meso_* routines were deleted and changes are
now impacted under the MESOSCALE flag in the original GCM routines
--> Completely transparent for GCM compilation since it is devoid
of the -DMESOSCALE option
--> Very good for syncing because changes in dustlift, newcondens
will be directly available in the mesoscale model

M 84 mesoscale/LMD_MM_MARS/makemeso
Changed meso_testphys1d in testphys1d

M 84 mesoscale/LMD_MM_MARS/SRC/WRFV2/mars_lmd_new/makegcm_pgf
M 84 mesoscale/LMD_MM_MARS/SRC/WRFV2/mars_lmd_new/makegcm_mpifort
M 84 mesoscale/LMD_MM_MARS/SRC/WRFV2/mars_lmd_new/makegcm_ifort
M 84 mesoscale/LMD_MM_MARS/SRC/WRFV2/mars_lmd_new/makegcm_g95
M 84 mesoscale/LMD_MM_MARS/SRC/WRFV2/mars_lmd_new/makegcm_mpi
Added the option -DMESOSCALE in these scripts

*
LMD_MM_MARS
* Various minor changes related to water cycle and plotting routines

* Also included the GW test case

*

A 0 mesoscale/LMDZ.MARS.new/myGCM/DEFS_JB/callphys.def.orig
M 84 mesoscale/NOTES.txt
D 84 mesoscale/LMD_MM_MARS/SRC/ARWpost/idl
M 84 mesoscale/LMD_MM_MARS/SRC/WRFV2/Registry/Registry.EM
M 84 mesoscale/LMD_MM_MARS/SIMU/gnome_launch.meso
M 85 mesoscale/PLOT/MINIMAL/map_latlon.pro
D 85 mesoscale/PLOT/SPEC/LES/getget.pro
M 85 mesoscale/PLOT/SPEC/MAP/map_uvt.pro
A + - mesoscale/PLOT/SPEC/getget.pro
A 0 mesoscale/PLOT/RESERVE/obsolete
A 0 mesoscale/TESTS/TESTGW.tar.gz
M 84 000-USERS

File size: 22.6 KB
Line 
1
2      PROGRAM testphys1d
3! to use  'getin'
4      USE ioipsl_getincom
5      IMPLICIT NONE
6
7c=======================================================================
8c   subject:
9c   --------
10c   PROGRAM useful to run physical part of the martian GCM in a 1D column
11c       
12c Can be compiled with a command like (e.g. for 25 layers)
13c  "makegcm -p mars -d 25 testphys1d"
14c It requires the files "testphys1d.def" "callphys.def"
15c   and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line)
16c      and a file describing the sigma layers (e.g. "z2sig.def")
17c
18c   author: Frederic Hourdin, R.Fournier,F.Forget
19c   -------
20c   
21c   update: 12/06/2003 including chemistry (S. Lebonnois)
22c                            and water ice (F. Montmessin)
23c
24c=======================================================================
25
26#include "dimensions.h"
27#include "dimphys.h"
28#include "dimradmars.h"
29#include "comgeomfi.h"
30#include "surfdat.h"
31#include "comsoil.h"
32#include "comdiurn.h"
33#include "callkeys.h"
34#include "comcstfi.h"
35#include "planete.h"
36#include "comsaison.h"
37#include "yomaer.h"
38#include "control.h"
39#include "comvert.h"
40#include "netcdf.inc"
41#include "comg1d.h"
42#include "watercap.h"
43#include "logic.h"
44#include "advtrac.h"
45
46c --------------------------------------------------------------
47c  Declarations
48c --------------------------------------------------------------
49c
50      INTEGER unitstart      ! unite d'ecriture de "startfi"
51      INTEGER nlayer,nlevel,nsoil,ndt
52      INTEGER ilayer,ilevel,isoil,idt,iq
53      LOGICAl firstcall,lastcall
54c
55      INTEGER day0          ! date initial (sol ; =0 a Ls=0)
56      REAL day           ! date durant le run
57      REAL time             ! time (0<time<1 ; time=0.5 a midi)
58      REAL play(nlayermx)   ! Pressure at the middle of the layers (Pa)
59      REAL plev(nlayermx+1) ! intermediate pressure levels (pa)
60      REAL psurf,tsurf     
61      REAL u(nlayermx),v(nlayermx)  ! zonal, meridional wind
62      REAL gru,grv   ! prescribed "geostrophic" background wind
63      REAL temp(nlayermx)   ! temperature at the middle of the layers
64      REAL q(nlayermx,nqmx) ! tracer mixing ratio (e.g. kg/kg)
65      REAL qsurf(nqmx)      ! tracer surface budget (e.g. kg.m-2)
66      REAL tsoil(nsoilmx)   ! subsurface soik temperature (K)
67      REAL co2ice           ! co2ice layer (kg.m-2)
68      REAL emis             ! surface layer
69      REAL q2(nlayermx+1)   ! Turbulent Kinetic Energy
70      REAL zlay(nlayermx)   ! altitude estimee dans les couches (km)
71
72c    Physical and dynamical tandencies (e.g.  m.s-2, K/s, Pa/s)
73      REAL du(nlayermx),dv(nlayermx),dtemp(nlayermx)
74      REAL dudyn(nlayermx),dvdyn(nlayermx),dtempdyn(nlayermx)
75      REAL dpsurf   
76      REAL dq(nlayermx,nqmx)
77      REAL dqdyn(nlayermx,nqmx)
78
79c   Various intermediate variables
80      INTEGER thermo
81      REAL zls
82      REAL phi(nlayermx),h(nlayermx),s(nlayermx)
83      REAL pks, ptif, w(nlayermx)
84      REAL qtotinit, mqtot(nqmx),qtot
85      INTEGER ierr, aslun
86      REAL tmp1(0:nlayermx),tmp2(0:nlayermx)
87      Logical  tracerdyn
88      integer :: nq=1 ! number of tracers
89
90      character*2 str2
91      character (len=7) :: str7
92      character(len=44) :: txt
93
94c=======================================================================
95
96c=======================================================================
97c INITIALISATION
98c=======================================================================
99
100c ------------------------------------------------------
101c  Constantes prescrites ICI
102c ------------------------------------------------------
103
104      pi=2.E+0*asin(1.E+0)
105
106c     Constante de la Planete Mars
107c     ----------------------------
108      rad=3397200.               ! rayon de mars (m)  ~3397200 m
109      daysec=88775.              ! duree du sol (s)  ~88775 s
110      omeg=4.*asin(1.)/(daysec)  ! vitesse de rotation (rad.s-1)
111      g=3.72                     ! gravite (m.s-2) ~3.72 
112      mugaz=43.49                ! Masse molaire de l'atm (g.mol-1) ~43.49
113      rcp=.256793                ! = r/cp  ~0.256793
114      r= 8.314511E+0 *1000.E+0/mugaz
115      cpp= r/rcp
116      year_day = 669           ! duree de l'annee (sols) ~668.6
117      periheli = 206.66          ! dist.min. soleil-mars (Mkm) ~206.66
118      aphelie = 249.22           ! dist.max. soleil-mars (Mkm) ~249.22
119      peri_day =  485.           ! date du perihelie (sols depuis printemps)
120      obliquit = 25.2           ! Obliquite de la planete (deg) ~25.2         
121 
122c     Parametres Couche limite et Turbulence
123c     --------------------------------------
124      z0 =  1.e-2                ! surface roughness (m) ~0.01
125      emin_turb = 1.e-6          ! energie minimale ~1.e-8
126      lmixmin = 30               ! longueur de melange ~100
127 
128c     propriete optiques des calottes et emissivite du sol
129c     ----------------------------------------------------
130      emissiv= 0.95              ! Emissivite du sol martien ~.95
131      emisice(1)=0.95            ! Emissivite calotte nord
132      emisice(2)=0.95            ! Emissivite calotte sud 
133      albedice(1)=0.5           ! Albedo calotte nord
134      albedice(2)=0.5            ! Albedo calotte sud
135      iceradius(1) = 100.e-6     ! mean scat radius of CO2 snow (north)
136      iceradius(2) = 100.e-6     ! mean scat radius of CO2 snow (south)
137      dtemisice(1) = 2.          ! time scale for snow metamorphism (north)
138      dtemisice(2) = 2.          ! time scale for snow metamorphism (south
139      hybrid=.false.
140
141c ------------------------------------------------------
142c  Lecture des parametres dans "run.def"
143c ------------------------------------------------------
144
145
146! check if 'run.def' file is around (otherwise reading parameters
147! from callphys.def via getin() routine won't work.
148      open(99,file='run.def',status='old',form='formatted',
149     &     iostat=ierr)
150      if (ierr.ne.0) then
151        write(*,*) 'Cannot find required file "run.def"'
152        write(*,*) '  (which should contain some input parameters'
153        write(*,*) '   along with the following line:'
154        write(*,*) '   INCLUDEDEF=callphys.def'
155        write(*,*) '   )'
156        write(*,*) ' ... might as well stop here ...'
157        stop
158      else
159        close(99)
160      endif
161
162! check if we are going to run with or without tracers
163      write(*,*) "Run with or without tracer transport ?"
164      tracer=.false. ! default value
165      call getin("tracer",tracer)
166      write(*,*) " tracer = ",tracer
167
168! while we're at it, check if there is a 'traceur.def' file
169! and preocess it, if necessary. Otherwise initialize tracer names
170      if (tracer) then
171      ! load tracer names from file 'traceur.def'
172        open(90,file='traceur.def',status='old',form='formatted',
173     &       iostat=ierr)
174        if (ierr.ne.0) then
175          write(*,*) 'Cannot find required file "traceur.def"'
176          write(*,*) ' If you want to run with tracers, I need it'
177          write(*,*) ' ... might as well stop here ...'
178          stop
179        else
180          write(*,*) "testphys1d: Reading file traceur.def"
181          ! read number of tracers:
182          read(90,*,iostat=ierr) nq
183          if (ierr.ne.0) then
184            write(*,*) "testphys1d: error reading number of tracers"
185            write(*,*) "   (first line of traceur.def) "
186            stop
187          else
188            ! check that the number of tracers is indeed nqmx
189            if (nq.ne.nqmx) then
190              write(*,*) "testphys1d: error, wrong number of tracers:"
191              write(*,*) "nq=",nq," whereas nqmx=",nqmx
192              stop
193            endif
194          endif
195        endif
196        ! read tracer names from file traceur.def
197        do iq=1,nqmx
198          read(90,*,iostat=ierr) tnom(iq)
199          if (ierr.ne.0) then
200            write(*,*) 'testphys1d: error reading tracer names...'
201            stop
202          endif
203        enddo
204        close(90)
205
206        ! initialize tracers here:
207        write(*,*) "testphys1d: initializing tracers"
208        q(:,:)=0 ! default, set everything to zero
209        qsurf(:)=0
210        ! "smarter" initialization of some tracers
211        ! (get values from "profile_*" files, if these are available)
212        do iq=1,nqmx
213          txt=""
214          write(txt,"(a)") tnom(iq)
215          write(*,*)"  tracer:",trim(txt)
216          ! CO2
217          if (txt.eq."co2") then
218            q(:,iq)=0.95   ! kg /kg of atmosphere
219            qsurf(iq)=0. ! kg/m2 (not used for CO2)
220            ! even better, look for a "profile_co2" input file
221            open(91,file='profile_co2',status='old',
222     &       form='formatted',iostat=ierr)
223            if (ierr.eq.0) then
224              read(91,*) qsurf(iq)
225              do ilayer=1,nlayermx
226                read(91,*) q(ilayer,iq)
227              enddo
228            endif
229            close(91)
230          endif ! of if (txt.eq."co2")
231          ! WATER VAPOUR
232          if (txt.eq."h2o_vap") then
233            !look for a "profile_h2o_vap" input file
234            open(91,file='profile_h2o_vap',status='old',
235     &       form='formatted',iostat=ierr)
236            if (ierr.eq.0) then
237              read(91,*) qsurf(iq)
238              do ilayer=1,nlayermx
239                read(91,*) q(ilayer,iq)
240              enddo
241            else
242              write(*,*) "No profile_h2o_vap file!"
243            endif
244            close(91)
245          endif ! of if (txt.eq."h2o_ice")
246          ! WATER ICE
247          if (txt.eq."h2o_ice") then
248            !look for a "profile_h2o_vap" input file
249            open(91,file='profile_h2o_ice',status='old',
250     &       form='formatted',iostat=ierr)
251            if (ierr.eq.0) then
252              read(91,*) qsurf(iq)
253              do ilayer=1,nlayermx
254                read(91,*) q(ilayer,iq)
255              enddo
256            else
257              write(*,*) "No profile_h2o_ice file!"
258            endif
259            close(91)
260          endif ! of if (txt.eq."h2o_ice")
261          ! DUST
262          !if (txt(1:4).eq."dust") then
263          !  q(:,iq)=0.4    ! kg/kg of atmosphere
264          !  qsurf(iq)=100 ! kg/m2
265          !endif
266          ! DUST MMR
267          if (txt.eq."dust_mass") then
268            !look for a "profile_dust_mass" input file
269            open(91,file='profile_dust_mass',status='old',
270     &       form='formatted',iostat=ierr)
271            if (ierr.eq.0) then
272              read(91,*) qsurf(iq)
273              do ilayer=1,nlayermx
274                read(91,*) q(ilayer,iq)
275!                write(*,*) "l=",ilayer," q(ilayer,iq)=",q(ilayer,iq)
276              enddo
277            else
278              write(*,*) "No profile_dust_mass file!"
279            endif
280            close(91)
281          endif ! of if (txt.eq."dust_mass")
282          ! DUST NUMBER
283          if (txt.eq."dust_number") then
284            !look for a "profile_dust_number" input file
285            open(91,file='profile_dust_number',status='old',
286     &       form='formatted',iostat=ierr)
287            if (ierr.eq.0) then
288              read(91,*) qsurf(iq)
289              do ilayer=1,nlayermx
290                read(91,*) q(ilayer,iq)
291              enddo
292            else
293              write(*,*) "No profile_dust_number file!"
294            endif
295            close(91)
296          endif ! of if (txt.eq."dust_number")
297        enddo ! of do iq=1,nqmx
298        ! NB: some more initializations (chemistry) is done later
299
300      else
301      ! we still need to set (dummy) tracer names for physdem1
302        nq=nqmx
303        do iq=1,nq
304          write(str7,'(a1,i2.2)')'q',iq
305          tnom(iq)=str7
306        enddo
307      endif ! of if (tracer)
308
309c  Date et heure locale du debut du run
310c  ------------------------------------
311c    Date (en sols depuis le solstice de printemps) du debut du run
312      day0 = 0 ! default value for day0
313      write(*,*) 'Initial date (in martian sols ; =0 at Ls=0)?'
314      call getin("day0",day0)
315      day=float(day0)
316      write(*,*) " day0 = ",day0
317c  Heure de demarrage
318      time=0 ! default value for time
319      write(*,*)'Initial local time (in hours, between 0 and 24)?'
320      call getin("time",time)
321      write(*,*)" time = ",time
322      time=time/24.E+0 ! convert time (hours) to fraction of sol
323
324c  Discretisation (Definition de la grille et des pas de temps)
325c  --------------
326c
327      nlayer=nlayermx
328      nlevel=nlayer+1
329      nsoil=nsoilmx
330
331      day_step=48 ! default value for day_step
332      PRINT *,'Number of time steps per sol ?'
333      call getin("day_step",day_step)
334      write(*,*) " day_step = ",day_step
335
336      ndt=10 ! default value for ndt
337      PRINT *,'Number of sols to run ?'
338      call getin("ndt",ndt)
339      write(*,*) " ndt = ",ndt
340
341      ndt=ndt*day_step     
342      dtphys=daysec/day_step 
343c Pression de surface sur la planete
344c ------------------------------------
345c
346      psurf=610. ! default value for psurf
347      PRINT *,'Surface pressure (Pa) ?'
348      call getin("psurf",psurf)
349      write(*,*) " psurf = ",psurf
350c Pression de reference
351      pa=20.
352      preff=610.     
353 
354c Proprietes des poussiere aerosol
355c --------------------------------
356      tauvis=0.2 ! default value for tauvis
357      print *,'Reference dust opacity at 700 Pa ?'
358      call getin("tauvis",tauvis)
359      write(*,*) " tauvis = ",tauvis
360 
361c  latitude/longitude
362c  ------------------
363      lati(1)=0 ! default value for lati(1)
364      PRINT *,'latitude (in degrees) ?'
365      call getin("latitude",lati(1))
366      write(*,*) " latitude = ",lati(1)
367      lati(1)=lati(1)*pi/180.E+0
368      long(1)=0.E+0
369      long(1)=long(1)*pi/180.E+0
370
371c  Initialisation albedo / inertie du sol
372c  --------------------------------------
373c
374      albedodat(1)=0.2 ! default value for albedodat
375      PRINT *,'Albedo of bare ground ?'
376      call getin("albedo",albedodat(1))
377      write(*,*) " albedo = ",albedodat(1)
378
379      inertiedat(1,1)=400 ! default value for inertiedat
380      PRINT *,'Soil thermal inertia (SI) ?'
381      call getin("inertia",inertiedat(1,1))
382      write(*,*) " inertia = ",inertiedat(1,1)
383c
384c  pour le schema d'ondes de gravite
385c  ---------------------------------
386c
387      zmea(1)=0.E+0
388      zstd(1)=0.E+0
389      zsig(1)=0.E+0
390      zgam(1)=0.E+0
391      zthe(1)=0.E+0
392
393
394
395c   Initialisation speciales "physiq"
396c   ---------------------------------
397c   la surface de chaque maille est inutile en 1D --->
398      area(1)=1.E+0
399
400c   le geopotentiel au sol est inutile en 1D car tout est controle
401c   par la pression de surface --->
402      phisfi(1)=0.E+0
403
404c  "inifis" reproduit un certain nombre d'initialisations deja faites
405c  + lecture des clefs de callphys.def
406c  + calcul de la frequence d'appel au rayonnement
407c  + calcul des sinus et cosinus des longitude latitude
408
409!Mars possible matter with dtphys in input and include!!!
410#ifdef MESOSCALE
411      CALL meso_inifis(1,llm,day0,daysec,dtphys,
412#else
413      CALL inifis(1,llm,day0,daysec,dtphys,
414#endif
415     .            lati,long,area,rad,g,r,cpp)
416c   Initialisation pour prendre en compte les vents en 1-D
417c   ------------------------------------------------------
418      ptif=2.E+0*omeg*sinlat(1)
419 
420c    vent geostrophique
421      gru=10. ! default value for gru
422      PRINT *,'zonal eastward component of the geostrophic wind (m/s) ?'
423      call getin("u",gru)
424      write(*,*) " u = ",gru
425      grv=0. !default value for grv
426      PRINT *,'meridional northward component of the geostrophic',
427     &' wind (m/s) ?'
428      call getin("v",grv)
429      write(*,*) " v = ",grv
430
431c     Initialisation des vents  au premier pas de temps
432      DO ilayer=1,nlayer
433         u(ilayer)=gru
434         v(ilayer)=grv
435      ENDDO
436
437c     energie cinetique turbulente
438      DO ilevel=1,nlevel
439         q2(ilevel)=0.E+0
440      ENDDO
441
442c  glace de CO2 au sol
443c  -------------------
444      co2ice=0.E+0 ! default value for co2ice
445      PRINT *,'Initial CO2 ice on the surface (kg.m-2)'
446      call getin("co2ice",co2ice)
447      write(*,*) " co2ice = ",co2ice
448
449c
450c  emissivite
451c  ----------
452      emis=emissiv
453      IF (co2ice.eq.1.E+0) THEN
454         emis=emisice(1) ! northern hemisphere
455         IF(lati(1).LT.0) emis=emisice(2) ! southern hemisphere
456      ENDIF
457
458 
459
460c  calcul des pressions et altitudes en utilisant les niveaux sigma
461c  ----------------------------------------------------------------
462
463c    Vertical Coordinates
464c    """"""""""""""""""""
465      hybrid=.true.
466      PRINT *,'Hybrid coordinates ?'
467      call getin("hybrid",hybrid)
468      write(*,*) " hybrid = ", hybrid
469
470      CALL  disvert
471
472      DO ilevel=1,nlevel
473        plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
474      ENDDO
475
476      DO ilayer=1,nlayer
477        play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
478      ENDDO
479
480      DO ilayer=1,nlayer
481        zlay(ilayer)=-200.E+0 *r*log(play(ilayer)/plev(1))
482     &   /g
483      ENDDO
484
485
486c  profil de temperature au premier appel
487c  --------------------------------------
488      pks=psurf**rcp
489
490c altitude en km dans profile: on divise zlay par 1000
491      tmp1(0)=0.E+0
492      DO ilayer=1,nlayer
493        tmp1(ilayer)=zlay(ilayer)/1000.E+0
494      ENDDO
495
496      call profile(nlayer+1,tmp1,tmp2)
497
498      tsurf=tmp2(0)
499      DO ilayer=1,nlayer
500        temp(ilayer)=tmp2(ilayer)
501      ENDDO
502     
503
504
505! Initialize soil properties and temperature
506! ------------------------------------------
507      volcapa=1.e6 ! volumetric heat capacity
508      DO isoil=1,nsoil
509         inertiedat(1,isoil)=inertiedat(1,1) ! soil thermal inertia
510         tsoil(isoil)=tsurf  ! soil temperature
511      ENDDO
512
513! Initialize depths
514! -----------------
515      do isoil=0,nsoil-1
516        mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) ! mid-layer depth
517      enddo
518      do isoil=1,nsoil
519        layer(isoil)=2.e-4*(2.**(isoil-1)) ! layer depth
520      enddo
521
522c    Initialisation des traceurs
523c    ---------------------------
524
525      if (photochem.or.callthermos) then
526         write(*,*) 'Especes chimiques initialisees'
527         ! thermo=0: initialisation pour toutes les couches
528         thermo=0
529         call inichim_newstart(q,psurf,sig,nqmx,lati,long,area,
530     $   thermo,qsurf)
531      endif
532      watercaptag(ngridmx)=.false.
533     
534
535c    Initialisation pour les sorties GRADS dans "g1d.dat" et "g1d.ctl"
536c    ----------------------------------------------------------------
537c    (effectuee avec "writeg1d" a partir de "physiq.F"
538c    ou tout autre subroutine physique, y compris celle ci).
539
540        g1d_nlayer=nlayer
541        g1d_nomfich='g1d.dat'
542        g1d_unitfich=40
543        g1d_nomctl='g1d.ctl'
544        g1d_unitctl=41
545        g1d_premier=.true.
546        g2d_premier=.true.
547
548c  Ecriture de "startfi"
549c  --------------------
550c  (Ce fichier sera aussitot relu au premier
551c   appel de "physiq", mais il est necessaire pour passer
552c   les variables purement physiques a "physiq"...
553
554      call physdem1("startfi.nc",long,lati,nsoilmx,nqmx,
555     .              dtphys,float(day0),
556     .              time,tsurf,tsoil,co2ice,emis,q2,qsurf,
557     .              area,albedodat,inertiedat,zmea,zstd,zsig,zgam,zthe)
558c=======================================================================
559c  BOUCLE TEMPORELLE DU MODELE 1D
560c=======================================================================
561c
562      firstcall=.true.
563      lastcall=.false.
564
565      DO idt=1,ndt
566c        IF (idt.eq.ndt) lastcall=.true.
567        IF (idt.eq.ndt-day_step-1) then       !test
568         lastcall=.true.
569         call solarlong(day*1.0,zls)
570         write(103,*) 'Ls=',zls*180./pi
571         write(103,*) 'Lat=', lati(1)*180./pi
572         write(103,*) 'Tau=', tauvis/700*psurf
573         write(103,*) 'RunEnd - Atmos. Temp. File'
574         write(103,*) 'RunEnd - Atmos. Temp. File'
575         write(104,*) 'Ls=',zls*180./pi
576         write(104,*) 'Lat=', lati(1)
577         write(104,*) 'Tau=', tauvis/700*psurf
578         write(104,*) 'RunEnd - Atmos. Temp. File'
579        ENDIF
580
581c    calcul du geopotentiel
582c     ~~~~~~~~~~~~~~~~~~~~~
583      DO ilayer=1,nlayer
584        s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp
585        h(ilayer)=cpp*temp(ilayer)/(pks*s(ilayer))
586      ENDDO
587      phi(1)=pks*h(1)*(1.E+0-s(1))
588      DO ilayer=2,nlayer
589         phi(ilayer)=phi(ilayer-1)+
590     &               pks*(h(ilayer-1)+h(ilayer))*.5E+0
591     &                  *(s(ilayer-1)-s(ilayer))
592
593      ENDDO
594
595c       appel de la physique
596c       --------------------
597#ifdef MESOSCALE
598      CALL meso_physiq (1,llm,nqmx,
599#else
600      CALL physiq (1,llm,nqmx,
601#endif
602     ,     firstcall,lastcall,
603     ,     day,time,dtphys,
604     ,     plev,play,phi,
605     ,     u, v,temp, q, 
606     ,     w,
607C - sorties
608     s     du, dv, dtemp, dq,dpsurf,tracerdyn)
609
610c       evolution du vent : modele 1D
611c       -----------------------------
612 
613c       la physique calcule les derivees temporelles de u et v.
614c       on y rajoute betement un effet Coriolis.
615c
616c       DO ilayer=1,nlayer
617c          du(ilayer)=du(ilayer)+ptif*(v(ilayer)-grv)
618c          dv(ilayer)=dv(ilayer)+ptif*(-u(ilayer)+gru)
619c       ENDDO
620
621c       Pour certain test : pas de coriolis a l'equateur
622c       if(lati(1).eq.0.) then
623          DO ilayer=1,nlayer
624             du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4
625             dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4
626          ENDDO
627c       end if
628c     
629c
630c       Calcul du temps au pas de temps suivant
631c       ---------------------------------------
632        firstcall=.false.
633        time=time+dtphys/daysec
634        IF (time.gt.1.E+0) then
635            time=time-1.E+0
636            day=day+1
637        ENDIF
638
639c       calcul des vitesses et temperature au pas de temps suivant
640c       ----------------------------------------------------------
641
642        DO ilayer=1,nlayer
643           u(ilayer)=u(ilayer)+dtphys*du(ilayer)
644           v(ilayer)=v(ilayer)+dtphys*dv(ilayer)
645           temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer)
646        ENDDO
647
648c       calcul des pressions au pas de temps suivant
649c       ----------------------------------------------------------
650
651           psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
652           DO ilevel=1,nlevel
653             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
654           ENDDO
655           DO ilayer=1,nlayer
656             play(ilayer)=aps(ilayer)+psurf*bps(ilayer)
657           ENDDO
658
659!       increment tracer values
660        DO iq = 1, nqmx
661          DO ilayer=1,nlayer
662             q(ilayer,iq)=q(ilayer,iq)+dtphys*dq(ilayer,iq)
663          ENDDO
664        ENDDO
665
666      ENDDO   ! fin de la boucle temporelle
667
668c    ========================================================
669c    GESTION DES SORTIE
670c    ========================================================
671
672c    fermeture pour conclure les sorties format grads dans "g1d.dat"
673c    et "g1d.ctl" (effectuee avec "writeg1d" a partir de "physiq.F"
674c    ou tout autre subroutine physique, y compris celle ci).
675
676c        CALL endg1d(1,nlayer,zphi/(g*1000.),ndt)
677        CALL endg1d(1,nlayer,zlay/1000.,ndt)
678
679c    ========================================================
680      END
681 
682c***********************************************************************
683c***********************************************************************
684c     Subroutines Bidons utilise seulement en 3D, mais
685c     necessaire a la compilation de testphys1d en 1-D
686
687      subroutine gr_fi_dyn
688      RETURN
689      END
690 
691c***********************************************************************
692c***********************************************************************
693
694#include "../dyn3d/disvert.F"
Note: See TracBrowser for help on using the repository browser.