source: LMDZ5/trunk/libf/dyn3d/gcm.F @ 2026

Last change on this file since 2026 was 2021, checked in by lguez, 10 years ago

Removed unused variables pks, pk, pkf from main program unit gcm.

Encapsulated procedures exner_hyb, exner_hyb_p, exner_hyb_loc,
exner_milieu, exner_milieu_p and exner_milieu_loc into
modules. (Compulsory to allow optional arguments.)

In the procedures exner_hyb, exner_hyb_p, exner_hyb_loc, donwgraded
arguments alpha and beta to local variables. In exner_milieu,
exner_milieu_p and exner_milieu_loc, removed beta altogether. In the
six procedures exner_*, made pkf an optional argument. Made some
cosmetic modifications in order to keep the six procedures exner_* as
close as possible.

In the six procedures exner_*, removed the averaging of pks at the
poles: this is not useful because ps is already the same at all
longitudes at the poles. This modification changes the results of the
program. Motivation: had to do this for exner_hyb because we call it
from test_disvert with a few surface pressure values.

In all the procedures calling exner_*, removed the variables alpha and
beta. Also removed variables alpha and beta from module leapfrog_mod
and from module call_calfis_mod.

Removed actual argument pkf in call to exner_hyb* and exner_milieu*
from guide_interp, guide_main, iniacademic and iniacademic_loc (pkf
was not used in those procedures).

Argument workvar of startget_dyn is used only if varname is tpot or

  1. When varname is tpot or q, the actual argument associated to

workvar in etat0_netcdf is not y. So y in etat0_netcdf is a
place-holder, never used. So we remove optional argument y in the
calls to exner_hyb and exner_milieu from etat0_netcdf.

Created procedure test_disvert, called only by etat0_netcdf. This
procedure tests the order of pressure values at half-levels and full
levels.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 16.2 KB
RevLine 
[1279]1!
2! $Id: gcm.F 2021 2014-04-25 10:20:14Z fhourdin $
3!
[524]4c
5c
6      PROGRAM gcm
7
8#ifdef CPP_IOIPSL
9      USE IOIPSL
[1146]10#else
11! if not using IOIPSL, we still need to use (a local version of) getin
12      USE ioipsl_getincom
[524]13#endif
[956]14
[1825]15
16#ifdef CPP_XIOS
17    ! ug Pour les sorties XIOS
18        USE wxios
19#endif
20
[1146]21      USE filtreg_mod
22      USE infotrac
[1403]23      USE control_mod
[1146]24
[1785]25#ifdef INCA
26! Only INCA needs these informations (from the Earth's physics)
27      USE indice_sol_mod
28#endif
29
[956]30!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
31! FH 2008/05/09 On elimine toutes les clefs physiques dans la dynamique
32! A nettoyer. On ne veut qu'une ou deux routines d'interface
33! dynamique -> physique pour l'initialisation
[1615]34#ifdef CPP_PHYS
[762]35      USE dimphy
36      USE comgeomphy
[863]37      USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb
[956]38#endif
39!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
40
[524]41      IMPLICIT NONE
42
43c      ......   Version  du 10/01/98    ..........
44
45c             avec  coordonnees  verticales hybrides
46c   avec nouveaux operat. dissipation * ( gradiv2,divgrad2,nxgraro2 )
47
48c=======================================================================
49c
50c   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
51c   -------
52c
53c   Objet:
54c   ------
55c
56c   GCM LMD nouvelle grille
57c
58c=======================================================================
59c
60c  ... Dans inigeom , nouveaux calculs pour les elongations  cu , cv
61c      et possibilite d'appeler une fonction f(y)  a derivee tangente
62c      hyperbolique a la  place de la fonction a derivee sinusoidale.
63c  ... Possibilite de choisir le schema pour l'advection de
64c        q  , en modifiant iadv dans traceur.def  (MAF,10/02) .
65c
66c      Pour Van-Leer + Vapeur d'eau saturee, iadv(1)=4. (F.Codron,10/99)
67c      Pour Van-Leer iadv=10
68c
69c-----------------------------------------------------------------------
70c   Declarations:
71c   -------------
72
73#include "dimensions.h"
74#include "paramet.h"
75#include "comconst.h"
76#include "comdissnew.h"
77#include "comvert.h"
78#include "comgeom.h"
79#include "logic.h"
80#include "temps.h"
[1403]81!!!!!!!!!!!#include "control.h"
[524]82#include "ener.h"
83#include "description.h"
84#include "serre.h"
[1403]85!#include "com_io_dyn.h"
[524]86#include "iniprint.h"
[541]87#include "tracstoke.h"
[1403]88#ifdef INCA
89! Only INCA needs these informations (from the Earth's physics)
[1785]90!#include "indicesol.h"
[1403]91#endif
[524]92      INTEGER         longcles
93      PARAMETER     ( longcles = 20 )
94      REAL  clesphy0( longcles )
95      SAVE  clesphy0
96
97
98
99      REAL zdtvr
100
101c   variables dynamiques
102      REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm) ! vents covariants
103      REAL teta(ip1jmp1,llm)                 ! temperature potentielle
[1146]104      REAL, ALLOCATABLE, DIMENSION(:,:,:):: q! champs advectes
[524]105      REAL ps(ip1jmp1)                       ! pression  au sol
106      REAL p (ip1jmp1,llmp1  )               ! pression aux interfac.des couches
107      REAL masse(ip1jmp1,llm)                ! masse d'air
108      REAL phis(ip1jmp1)                     ! geopotentiel au sol
109      REAL phi(ip1jmp1,llm)                  ! geopotentiel
110      REAL w(ip1jmp1,llm)                    ! vitesse verticale
111
112c variables dynamiques intermediaire pour le transport
113
114c   variables pour le fichier histoire
115      REAL dtav      ! intervalle de temps elementaire
116
117      REAL time_0
118
119      LOGICAL lafin
120      INTEGER ij,iq,l,i,j
121
122
123      real time_step, t_wrt, t_ops
124
125      LOGICAL first
126
127      LOGICAL call_iniphys
128      data call_iniphys/.true./
129
130c+jld variables test conservation energie
[962]131c      REAL ecin(ip1jmp1,llm),ecin0(ip1jmp1,llm)
[524]132C     Tendance de la temp. potentiel d (theta)/ d t due a la
133C     tansformation d'energie cinetique en energie thermique
134C     cree par la dissipation
135      REAL dhecdt(ip1jmp1,llm)
[962]136c      REAL vcont(ip1jm,llm),ucont(ip1jmp1,llm)
137c      REAL      d_h_vcol, d_qt, d_qw, d_ql, d_ec
138      CHARACTER (len=15) :: ztit
[524]139c-jld
140
141
[962]142      character (len=80) :: dynhist_file, dynhistave_file
143      character (len=20) :: modname
144      character (len=80) :: abort_message
[1279]145! locales pour gestion du temps
146      INTEGER :: an, mois, jour
147      REAL :: heure
[524]148
149
150c-----------------------------------------------------------------------
151c    variables pour l'initialisation de la physique :
152c    ------------------------------------------------
[1146]153      INTEGER ngridmx
[524]154      PARAMETER( ngridmx = 2+(jjm-1)*iim - 1/jjm   )
155      REAL zcufi(ngridmx),zcvfi(ngridmx)
156      REAL latfi(ngridmx),lonfi(ngridmx)
157      REAL airefi(ngridmx)
158      SAVE latfi, lonfi, airefi
159
160c-----------------------------------------------------------------------
161c   Initialisations:
162c   ----------------
163
164      abort_message = 'last timestep reached'
165      modname = 'gcm'
166      descript = 'Run GCM LMDZ'
167      lafin    = .FALSE.
168      dynhist_file = 'dyn_hist.nc'
169      dynhistave_file = 'dyn_hist_ave.nc'
170
[762]171
172
[524]173c----------------------------------------------------------------------
174c  lecture des fichiers gcm.def ou run.def
175c  ---------------------------------------
176c
[1146]177! Ehouarn: dump possibility of using defrun
178!#ifdef CPP_IOIPSL
[524]179      CALL conf_gcm( 99, .TRUE. , clesphy0 )
[1146]180!#else
181!      CALL defrun( 99, .TRUE. , clesphy0 )
182!#endif
[762]183
[956]184!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[1825]185! Initialisation de XIOS
186!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
187
188#ifdef CPP_XIOS
189        CALL wxios_init("LMDZ")
190#endif
191
192
193!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[956]194! FH 2008/05/02
195! A nettoyer. On ne veut qu'une ou deux routines d'interface
196! dynamique -> physique pour l'initialisation
[1615]197#ifdef CPP_PHYS
[1403]198      CALL Init_Phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
[762]199      call InitComgeomphy
[956]200#endif
201!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[1279]202c-----------------------------------------------------------------------
203c   Choix du calendrier
204c   -------------------
[762]205
[1279]206c      calend = 'earth_365d'
207
208#ifdef CPP_IOIPSL
209      if (calend == 'earth_360d') then
210        call ioconf_calendar('360d')
211        write(lunout,*)'CALENDRIER CHOISI: Terrestre a 360 jours/an'
212      else if (calend == 'earth_365d') then
213        call ioconf_calendar('noleap')
214        write(lunout,*)'CALENDRIER CHOISI: Terrestre a 365 jours/an'
215      else if (calend == 'earth_366d') then
216        call ioconf_calendar('gregorian')
217        write(lunout,*)'CALENDRIER CHOISI: Terrestre bissextile'
218      else
219        abort_message = 'Mauvais choix de calendrier'
220        call abort_gcm(modname,abort_message,1)
221      endif
222#endif
223c-----------------------------------------------------------------------
224
[1563]225      IF (type_trac == 'inca') THEN
[762]226#ifdef INCA
[1315]227      call init_const_lmdz(nbtr,anneeref,dayref,iphysiq,day_step,nday,
228     $        nbsrf, is_oce,is_sic,is_ter,is_lic)
[863]229      call init_inca_para(iim,jjm+1,klon,1,klon_mpi_para_nb,0)
[762]230#endif
[960]231      END IF
[524]232c
233c
[762]234c------------------------------------
235c   Initialisation partie parallele
236c------------------------------------
237
238c
239c
[524]240c-----------------------------------------------------------------------
241c   Initialisation des traceurs
242c   ---------------------------
[1146]243c  Choix du nombre de traceurs et du schema pour l'advection
244c  dans fichier traceur.def, par default ou via INCA
245      call infotrac_init
[524]246
[1146]247c Allocation de la tableau q : champs advectes   
248      allocate(q(ip1jmp1,llm,nqtot))
249
[524]250c-----------------------------------------------------------------------
251c   Lecture de l'etat initial :
252c   ---------------------------
253
254c  lecture du fichier start.nc
255      if (read_start) then
[1146]256      ! we still need to run iniacademic to initialize some
[1403]257      ! constants & fields, if we run the 'newtonian' or 'SW' cases:
258        if (iflag_phys.ne.1) then
[1146]259          CALL iniacademic(vcov,ucov,teta,q,masse,ps,phis,time_0)
260        endif
[1403]261
[1454]262!        if (planet_type.eq."earth") then
[1146]263! Load an Earth-format start file
264         CALL dynetat0("start.nc",vcov,ucov,
[1403]265     &              teta,q,masse,ps,phis, time_0)
[1454]266!        endif ! of if (planet_type.eq."earth")
[1403]267       
[524]268c       write(73,*) 'ucov',ucov
269c       write(74,*) 'vcov',vcov
270c       write(75,*) 'teta',teta
271c       write(76,*) 'ps',ps
272c       write(77,*) 'q',q
273
[1146]274      endif ! of if (read_start)
[524]275
[1563]276      IF (type_trac == 'inca') THEN
[762]277#ifdef INCA
[960]278         call init_inca_dim(klon,llm,iim,jjm,
279     $        rlonu,rlatu,rlonv,rlatv)
[762]280#endif
[960]281      END IF
[524]282
283
284c le cas echeant, creation d un etat initial
285      IF (prt_level > 9) WRITE(lunout,*)
[1146]286     .              'GCM: AVANT iniacademic AVANT AVANT AVANT AVANT'
[524]287      if (.not.read_start) then
[1146]288         CALL iniacademic(vcov,ucov,teta,q,masse,ps,phis,time_0)
[524]289      endif
290
291
292c-----------------------------------------------------------------------
293c   Lecture des parametres de controle pour la simulation :
294c   -------------------------------------------------------
295c  on recalcule eventuellement le pas de temps
296
297      IF(MOD(day_step,iperiod).NE.0) THEN
298        abort_message =
299     .  'Il faut choisir un nb de pas par jour multiple de iperiod'
300        call abort_gcm(modname,abort_message,1)
301      ENDIF
302
303      IF(MOD(day_step,iphysiq).NE.0) THEN
304        abort_message =
305     * 'Il faut choisir un nb de pas par jour multiple de iphysiq'
306        call abort_gcm(modname,abort_message,1)
307      ENDIF
308
[1403]309      zdtvr    = daysec/REAL(day_step)
[524]310        IF(dtvr.NE.zdtvr) THEN
311         WRITE(lunout,*)
312     .    'WARNING!!! changement de pas de temps',dtvr,'>',zdtvr
313        ENDIF
314
315C
316C on remet le calendrier à zero si demande
317c
[1577]318      IF (start_time /= starttime) then
319        WRITE(lunout,*)' GCM: Attention l''heure de depart lue dans le'
320     &,' fichier restart ne correspond pas à celle lue dans le run.def'
321        IF (raz_date == 1) then
322          WRITE(lunout,*)'Je prends l''heure lue dans run.def'
323          start_time = starttime
324        ELSE
[1930]325          call abort_gcm("gcm", "'Je m''arrete'", 1)
[1577]326        ENDIF
327      ENDIF
[1403]328      IF (raz_date == 1) THEN
329        annee_ref = anneeref
330        day_ref = dayref
331        day_ini = dayref
332        itau_dyn = 0
333        itau_phy = 0
334        time_0 = 0.
[524]335        write(lunout,*)
[1403]336     .   'GCM: On reinitialise a la date lue dans gcm.def'
337      ELSE IF (annee_ref .ne. anneeref .or. day_ref .ne. dayref) THEN
338        write(lunout,*)
[1146]339     .  'GCM: Attention les dates initiales lues dans le fichier'
[524]340        write(lunout,*)
341     .  ' restart ne correspondent pas a celles lues dans '
342        write(lunout,*)' gcm.def'
[1403]343        write(lunout,*)' annee_ref=',annee_ref," anneeref=",anneeref
344        write(lunout,*)' day_ref=',day_ref," dayref=",dayref
345        write(lunout,*)' Pas de remise a zero'
346      ENDIF
[1279]347
[1403]348c      if (annee_ref .ne. anneeref .or. day_ref .ne. dayref) then
349c        write(lunout,*)
350c     .  'GCM: Attention les dates initiales lues dans le fichier'
351c        write(lunout,*)
352c     .  ' restart ne correspondent pas a celles lues dans '
353c        write(lunout,*)' gcm.def'
354c        write(lunout,*)' annee_ref=',annee_ref," anneeref=",anneeref
355c        write(lunout,*)' day_ref=',day_ref," dayref=",dayref
356c        if (raz_date .ne. 1) then
357c          write(lunout,*)
358c     .    'GCM: On garde les dates du fichier restart'
359c        else
360c          annee_ref = anneeref
361c          day_ref = dayref
362c          day_ini = dayref
363c          itau_dyn = 0
364c          itau_phy = 0
365c          time_0 = 0.
366c          write(lunout,*)
367c     .   'GCM: On reinitialise a la date lue dans gcm.def'
368c        endif
369c      ELSE
370c        raz_date = 0
371c      endif
372
[1147]373#ifdef CPP_IOIPSL
[1279]374      mois = 1
375      heure = 0.
376      call ymds2ju(annee_ref, mois, day_ref, heure, jD_ref)
377      jH_ref = jD_ref - int(jD_ref)
378      jD_ref = int(jD_ref)
379
380      call ioconf_startdate(INT(jD_ref), jH_ref)
381
382      write(lunout,*)'DEBUG'
383      write(lunout,*)'annee_ref, mois, day_ref, heure, jD_ref'
384      write(lunout,*)annee_ref, mois, day_ref, heure, jD_ref
385      call ju2ymds(jD_ref+jH_ref,an, mois, jour, heure)
386      write(lunout,*)'jD_ref+jH_ref,an, mois, jour, heure'
387      write(lunout,*)jD_ref+jH_ref,an, mois, jour, heure
388#else
389! Ehouarn: we still need to define JD_ref and JH_ref
390! and since we don't know how many days there are in a year
391! we set JD_ref to 0 (this should be improved ...)
392      jD_ref=0
393      jH_ref=0
[1147]394#endif
[524]395
396
[1403]397      if (iflag_phys.eq.1) then
398      ! these initialisations have already been done (via iniacademic)
399      ! if running in SW or Newtonian mode
[524]400c-----------------------------------------------------------------------
401c   Initialisation des constantes dynamiques :
402c   ------------------------------------------
[1403]403        dtvr = zdtvr
404        CALL iniconst
[524]405
406c-----------------------------------------------------------------------
407c   Initialisation de la geometrie :
408c   --------------------------------
[1403]409        CALL inigeom
[524]410
411c-----------------------------------------------------------------------
412c   Initialisation du filtre :
413c   --------------------------
[1403]414        CALL inifilr
415      endif ! of if (iflag_phys.eq.1)
[524]416c
417c-----------------------------------------------------------------------
418c   Initialisation de la dissipation :
419c   ----------------------------------
420
421      CALL inidissip( lstardis, nitergdiv, nitergrot, niterh   ,
[1697]422     *                tetagdiv, tetagrot , tetatemp, vert_prof_dissip)
[524]423
424c-----------------------------------------------------------------------
425c   Initialisation de la physique :
426c   -------------------------------
[1146]427
[1529]428      IF (call_iniphys.and.(iflag_phys==1.or.iflag_phys>=100)) THEN
[524]429         latfi(1)=rlatu(1)
430         lonfi(1)=0.
431         zcufi(1) = cu(1)
432         zcvfi(1) = cv(1)
433         DO j=2,jjm
434            DO i=1,iim
435               latfi((j-2)*iim+1+i)= rlatu(j)
436               lonfi((j-2)*iim+1+i)= rlonv(i)
437               zcufi((j-2)*iim+1+i) = cu((j-1)*iip1+i)
438               zcvfi((j-2)*iim+1+i) = cv((j-1)*iip1+i)
439            ENDDO
440         ENDDO
441         latfi(ngridmx)= rlatu(jjp1)
442         lonfi(ngridmx)= 0.
443         zcufi(ngridmx) = cu(ip1jm+1)
444         zcvfi(ngridmx) = cv(ip1jm-iim)
445         CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,aire,airefi)
446         WRITE(lunout,*)
[1146]447     .       'GCM: WARNING!!! vitesse verticale nulle dans la physique'
[1615]448! Physics:
449#ifdef CPP_PHYS
[1671]450         CALL iniphysiq(ngridmx,llm,daysec,day_ini,dtphys/nsplit_phys,
451     &                latfi,lonfi,airefi,zcufi,zcvfi,rad,g,r,cpp,
452     &                iflag_phys)
[1146]453#endif
[524]454         call_iniphys=.false.
[1146]455      ENDIF ! of IF (call_iniphys.and.(iflag_phys.eq.1))
[524]456
457c  numero de stockage pour les fichiers de redemarrage:
458
459c-----------------------------------------------------------------------
460c   Initialisation des I/O :
461c   ------------------------
462
463
464      day_end = day_ini + nday
465      WRITE(lunout,300)day_ini,day_end
[1146]466 300  FORMAT('1'/,15x,'run du jour',i7,2x,'au jour',i7//)
[524]467
[1279]468#ifdef CPP_IOIPSL
469      call ju2ymds(jD_ref + day_ini - day_ref, an, mois, jour, heure)
470      write (lunout,301)jour, mois, an
471      call ju2ymds(jD_ref + day_end - day_ref, an, mois, jour, heure)
472      write (lunout,302)jour, mois, an
473 301  FORMAT('1'/,15x,'run du ', i2,'/',i2,'/',i4)
474 302  FORMAT('1'/,15x,'    au ', i2,'/',i2,'/',i4)
475#endif
476
[1454]477!      if (planet_type.eq."earth") then
478! Write an Earth-format restart file
[1529]479
[1279]480        CALL dynredem0("restart.nc", day_end, phis)
[1454]481!      endif
[524]482
483      ecripar = .TRUE.
484
[1146]485#ifdef CPP_IOIPSL
[524]486      time_step = zdtvr
[1403]487      if (ok_dyn_ins) then
488        ! initialize output file for instantaneous outputs
489        ! t_ops = iecri * daysec ! do operations every t_ops
490        t_ops =((1.0*iecri)/day_step) * daysec 
491        t_wrt = daysec ! iecri * daysec ! write output every t_wrt
492        CALL inithist(day_ref,annee_ref,time_step,
493     &              t_ops,t_wrt)
494      endif
[524]495
[1403]496      IF (ok_dyn_ave) THEN
497        ! initialize output file for averaged outputs
498        t_ops = iperiod * time_step ! do operations every t_ops
499        t_wrt = periodav * daysec   ! write output every t_wrt
500        CALL initdynav(day_ref,annee_ref,time_step,
501     &       t_ops,t_wrt)
502      END IF
[524]503      dtav = iperiod*dtvr/daysec
504#endif
[1146]505! #endif of #ifdef CPP_IOIPSL
[524]506
[541]507c  Choix des frequences de stokage pour le offline
508c      istdyn=day_step/4     ! stockage toutes les 6h=1jour/4
509c      istdyn=day_step/12     ! stockage toutes les 2h=1jour/12
510      istdyn=day_step/4     ! stockage toutes les 6h=1jour/12
511      istphy=istdyn/iphysiq     
512
513
[524]514c
515c-----------------------------------------------------------------------
516c   Integration temporelle du modele :
517c   ----------------------------------
518
519c       write(78,*) 'ucov',ucov
520c       write(78,*) 'vcov',vcov
521c       write(78,*) 'teta',teta
522c       write(78,*) 'ps',ps
523c       write(78,*) 'q',q
524
525
[1146]526      CALL leapfrog(ucov,vcov,teta,ps,masse,phis,q,clesphy0,
[550]527     .              time_0)
[524]528
529      END
530
Note: See TracBrowser for help on using the repository browser.