/[lmdze]/trunk/Sources/dyn3d/leapfrog.f
ViewVC logotype

Diff of /trunk/Sources/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 18 by guez, Thu Aug 7 12:29:13 2008 UTC revision 29 by guez, Tue Mar 30 10:44:42 2010 UTC
# Line 1  Line 1 
1  module leapfrog_m  module leapfrog_m
2    
   ! This module is clean: no C preprocessor directive, no include line.  
   
3    IMPLICIT NONE    IMPLICIT NONE
4    
5  contains  contains
6    
7    SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, nq, q, time_0)    SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, q, time_0)
   
     ! From dyn3d/leapfrog.F, version 1.6 2005/04/13 08:58:34  
   
     ! Version du 10/01/98, avec coordonnees verticales hybrides, avec  
     ! nouveaux operat. dissipation * (gradiv2, divgrad2, nxgraro2)  
8    
9      ! Auteur: P. Le Van /L. Fairhead/F.Hourdin      ! From dyn3d/leapfrog.F, version 1.6, 2005/04/13 08:58:34
10      ! Objet:      ! Authors: P. Le Van, L. Fairhead, F. Hourdin
     ! GCM LMD nouvelle grille  
   
     ! ... Dans inigeom, nouveaux calculs pour les elongations cu, cv  
     ! et possibilite d'appeler une fonction f(y) a derivee tangente  
     ! hyperbolique a la place de la fonction a derivee sinusoidale.  
   
     ! ... Possibilité de choisir le schéma pour l'advection de  
     ! q, en modifiant iadv dans "traceur.def" (10/02) .  
   
     ! Pour Van-Leer + Vapeur d'eau saturee, iadv(1)=4. (F.Codron, 10/99)  
     ! Pour Van-Leer iadv=10  
   
     use dimens_m, only: iim, jjm, llm, nqmx  
     use paramet_m, only: ip1jmp1, ip1jm, ijmllm, ijp1llm, jjp1, iip1, iip2  
     use comconst, only: dtvr, daysec, dtphys  
     use comvert, only: ap, bp  
     use conf_gcm_m, only: day_step, iconser, idissip, iphysiq, iperiod, nday, &  
          offline, periodav  
     use logic, only: ok_guide, iflag_phys  
     use comgeom  
     use serre  
     use temps, only: itaufin, day_ini, dt  
     use iniprint, only: prt_level  
     use com_io_dyn  
     use ener  
     use calfis_m, only: calfis  
     use exner_hyb_m, only: exner_hyb  
     use guide_m, only: guide  
     use pression_m, only: pression  
     use pressure_var, only: p3d  
11    
12      integer, intent(in):: nq      USE calfis_m, ONLY: calfis
13        USE com_io_dyn, ONLY: histaveid
14        USE comconst, ONLY: daysec, dtphys, dtvr
15        USE comgeom, ONLY: aire_2d, apoln, apols
16        USE comvert, ONLY: ap, bp
17        USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &
18             periodav
19        USE dimens_m, ONLY: iim, jjm, llm, nqmx
20        USE dynetat0_m, ONLY: day_ini
21        use dynredem1_m, only: dynredem1
22        USE exner_hyb_m, ONLY: exner_hyb
23        use filtreg_m, only: filtreg
24        USE guide_m, ONLY: guide
25        use inidissip_m, only: idissip
26        USE logic, ONLY: iflag_phys, ok_guide
27        USE paramet_m, ONLY: ip1jmp1
28        USE pression_m, ONLY: pression
29        USE pressure_var, ONLY: p3d
30        USE temps, ONLY: itau_dyn
31    
32      ! Variables dynamiques:      ! Variables dynamiques:
33      REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm) ! vents covariants      REAL vcov((iim + 1) * jjm, llm), ucov(ip1jmp1, llm) ! vents covariants
34      REAL teta(ip1jmp1, llm) ! temperature potentielle      REAL, intent(inout):: teta(iim + 1, jjm + 1, llm) ! potential temperature
35      REAL q(ip1jmp1, llm, nqmx) ! mass fractions of advected fields      REAL ps(iim + 1, jjm + 1) ! pression au sol, en Pa
36      REAL ps(ip1jmp1) ! pression au sol, en Pa  
37      REAL masse(ip1jmp1, llm) ! masse d'air      REAL masse(ip1jmp1, llm) ! masse d'air
38      REAL phis(ip1jmp1) ! geopotentiel au sol      REAL phis(ip1jmp1) ! geopotentiel au sol
39        REAL q(ip1jmp1, llm, nqmx) ! mass fractions of advected fields
40      REAL time_0      REAL, intent(in):: time_0
41    
42      ! Variables local to the procedure:      ! Variables local to the procedure:
43    
44      ! Variables dynamiques:      ! Variables dynamiques:
45    
46      REAL pks(ip1jmp1) ! exner au sol      REAL pks(ip1jmp1) ! exner au sol
47      REAL pk(ip1jmp1, llm) ! exner au milieu des couches      REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
48      REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches      REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches
49      REAL phi(ip1jmp1, llm) ! geopotential      REAL phi(ip1jmp1, llm) ! geopotential
50      REAL w(ip1jmp1, llm) ! vitesse verticale      REAL w(ip1jmp1, llm) ! vitesse verticale
51    
52      ! variables dynamiques intermediaire pour le transport      ! variables dynamiques intermediaire pour le transport
53      REAL pbaru(ip1jmp1, llm), pbarv(ip1jm, llm) !flux de masse      REAL pbaru(ip1jmp1, llm), pbarv((iim + 1) * jjm, llm) !flux de masse
54    
55      ! variables dynamiques au pas - 1      ! variables dynamiques au pas - 1
56      REAL vcovm1(ip1jm, llm), ucovm1(ip1jmp1, llm)      REAL vcovm1((iim + 1) * jjm, llm), ucovm1(ip1jmp1, llm)
57      REAL tetam1(ip1jmp1, llm), psm1(ip1jmp1)      REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
58      REAL massem1(ip1jmp1, llm)      REAL massem1(ip1jmp1, llm)
59    
60      ! tendances dynamiques      ! tendances dynamiques
61      REAL dv(ip1jm, llm), du(ip1jmp1, llm)      REAL dv((iim + 1) * jjm, llm), du(ip1jmp1, llm)
62      REAL dteta(ip1jmp1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)      REAL dteta(ip1jmp1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)
63    
64      ! tendances de la dissipation      ! tendances de la dissipation
65      REAL dvdis(ip1jm, llm), dudis(ip1jmp1, llm)      REAL dvdis((iim + 1) * jjm, llm), dudis(ip1jmp1, llm)
66      REAL dtetadis(ip1jmp1, llm)      REAL dtetadis(iim + 1, jjm + 1, llm)
67    
68      ! tendances physiques      ! tendances physiques
69      REAL dvfi(ip1jm, llm), dufi(ip1jmp1, llm)      REAL dvfi((iim + 1) * jjm, llm), dufi(ip1jmp1, llm)
70      REAL dtetafi(ip1jmp1, llm), dqfi(ip1jmp1, llm, nqmx), dpfi(ip1jmp1)      REAL dtetafi(ip1jmp1, llm), dqfi(ip1jmp1, llm, nqmx), dpfi(ip1jmp1)
71    
72      ! variables pour le fichier histoire      ! variables pour le fichier histoire
73    
74      REAL tppn(iim), tpps(iim), tpn, tps      INTEGER itau ! index of the time step of the dynamics, starts at 0
75        INTEGER itaufin
     INTEGER itau, itaufinp1  
76      INTEGER iday ! jour julien      INTEGER iday ! jour julien
77      REAL time ! Heure de la journee en fraction d'1 jour      REAL time ! time of day, as a fraction of day length
   
     REAL SSUM  
78      real finvmaold(ip1jmp1, llm)      real finvmaold(ip1jmp1, llm)
79        LOGICAL:: lafin=.false.
80      LOGICAL :: lafin=.false.      INTEGER l
     INTEGER ij, l  
81    
82      REAL rdayvrai, rdaym_ini      REAL rdayvrai, rdaym_ini
     LOGICAL:: callinigrads = .true.  
83    
84      !+jld variables test conservation energie      ! Variables test conservation energie
85      REAL ecin(ip1jmp1, llm), ecin0(ip1jmp1, llm)      REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
86      ! Tendance de la temp. potentiel d (theta) / d t due a la      ! Tendance de la temp. potentiel d (theta) / d t due a la
87      ! tansformation d'energie cinetique en energie thermique      ! tansformation d'energie cinetique en energie thermique
88      ! cree par la dissipation      ! cree par la dissipation
89      REAL dtetaecdt(ip1jmp1, llm)      REAL dtetaecdt(iim + 1, jjm + 1, llm)
90      REAL vcont(ip1jm, llm), ucont(ip1jmp1, llm)      REAL vcont((iim + 1) * jjm, llm), ucont(ip1jmp1, llm)
91      CHARACTER*15 ztit      logical forward, leapf
92      INTEGER:: ip_ebil_dyn = 0 ! PRINT level for energy conserv. diag.      REAL dt
   
     logical:: dissip_conservative = .true.  
     LOGICAL:: prem = .true.  
     logical forward, leapf, apphys, conser, apdiss  
93    
94      !---------------------------------------------------      !---------------------------------------------------
95    
96      print *, "Call sequence information: leapfrog"      print *, "Call sequence information: leapfrog"
97    
98      itaufin = nday * day_step      itaufin = nday * day_step
     itaufinp1 = itaufin + 1  
   
99      itau = 0      itau = 0
100      iday = day_ini      iday = day_ini
101      time = time_0      time = time_0
102      IF (time > 1.) THEN      dq = 0.
        time = time - 1.  
        iday = iday + 1  
     ENDIF  
   
103      ! On initialise la pression et la fonction d'Exner :      ! On initialise la pression et la fonction d'Exner :
     dq=0.  
104      CALL pression(ip1jmp1, ap, bp, ps, p3d)      CALL pression(ip1jmp1, ap, bp, ps, p3d)
105      CALL exner_hyb(ps, p3d, pks, pk, pkf)      CALL exner_hyb(ps, p3d, pks, pk, pkf)
106    
107      ! Debut de l'integration temporelle:      ! Début de l'integration temporelle :
108      outer_loop:do      outer_loop:do
109         if (ok_guide.and.(itaufin - itau - 1) * dtvr > 21600) then         if (ok_guide .and. (itaufin - itau - 1) * dtvr > 21600.) &
110            call guide(itau, ucov, vcov, teta, q, masse, ps)              call guide(itau, ucov, vcov, teta, q, masse, ps)
111         else         vcovm1 = vcov
112            IF (prt_level > 9) print *, &         ucovm1 = ucov
113                 'Attention : on ne guide pas les 6 dernieres heures.'         tetam1 = teta
114         endif         massem1 = masse
115           psm1 = ps
        CALL SCOPY(ijmllm, vcov, 1, vcovm1, 1)  
        CALL SCOPY(ijp1llm, ucov, 1, ucovm1, 1)  
        CALL SCOPY(ijp1llm, teta, 1, tetam1, 1)  
        CALL SCOPY(ijp1llm, masse, 1, massem1, 1)  
        CALL SCOPY(ip1jmp1, ps, 1, psm1, 1)  
   
116         forward = .TRUE.         forward = .TRUE.
117         leapf = .FALSE.         leapf = .FALSE.
118         dt = dtvr         dt = dtvr
119           finvmaold = masse
120         CALL SCOPY(ijp1llm, masse, 1, finvmaold, 1)         CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)
        CALL filtreg(finvmaold, jjp1, llm, - 2, 2, .TRUE., 1)  
121    
122         do         do
123            ! gestion des appels de la physique et des dissipations:            ! Calcul des tendances dynamiques:
   
           apphys = .FALSE.  
           conser = .FALSE.  
           apdiss = .FALSE.  
   
           IF (MOD(itau, iconser) == 0) conser = .TRUE.  
           IF (MOD(itau + 1, idissip) == 0) apdiss = .TRUE.  
           IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) apphys=.TRUE.  
   
           ! calcul des tendances dynamiques:  
   
124            CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)            CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
   
125            CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &            CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
126                 conser, du, dv, dteta, dp, w, pbaru, pbarv, &                 MOD(itau, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &
127                 time + iday - day_ini)                 time + iday - day_ini)
128    
           ! calcul des tendances advection des traceurs (dont l'humidite)  
   
129            IF (forward .OR. leapf) THEN            IF (forward .OR. leapf) THEN
130                 ! Calcul des tendances advection des traceurs (dont l'humidité)
131               CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)               CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)
132               IF (offline) THEN               IF (offline) THEN
133                  !maf stokage du flux de masse pour traceurs OFF-LINE                  ! Stokage du flux de masse pour traceurs off-line
134                  CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, dtvr, &                  CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, dtvr, &
135                       itau)                       itau)
136               ENDIF               ENDIF
# Line 194  contains Line 139  contains
139            ! integrations dynamique et traceurs:            ! integrations dynamique et traceurs:
140            CALL integrd(2, vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, &            CALL integrd(2, vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, &
141                 dteta, dq, dp, vcov, ucov, teta, q, ps, masse, phis, &                 dteta, dq, dp, vcov, ucov, teta, q, ps, masse, phis, &
142                 finvmaold, leapf)                 finvmaold, leapf, dt)
   
           ! calcul des tendances physiques:  
143    
144            IF (apphys) THEN            IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
145                 ! calcul des tendances physiques:
146               IF (itau + 1 == itaufin) lafin = .TRUE.               IF (itau + 1 == itaufin) lafin = .TRUE.
147    
148               CALL pression(ip1jmp1, ap, bp, ps, p3d)               CALL pression(ip1jmp1, ap, bp, ps, p3d)
# Line 207  contains Line 151  contains
151               rdaym_ini = itau * dtvr / daysec               rdaym_ini = itau * dtvr / daysec
152               rdayvrai = rdaym_ini + day_ini               rdayvrai = rdaym_ini + day_ini
153    
154               ! Interface avec les routines de phylmd (phymars ...)               CALL calfis(nqmx, lafin, rdayvrai, time, ucov, vcov, teta, q, &
   
              ! Diagnostique de conservation de l'énergie : initialisation  
              IF (ip_ebil_dyn >= 1) THEN  
                 ztit='bil dyn'  
                 CALL diagedyn(ztit, 2, 1, 1, dtphys, ucov, vcov, ps, p3d, pk, &  
                      teta, q(:, :, 1), q(:, :, 2))  
              ENDIF  
   
              CALL calfis(nq, lafin, rdayvrai, time, ucov, vcov, teta, q, &  
155                    masse, ps, pk, phis, phi, du, dv, dteta, dq, w, &                    masse, ps, pk, phis, phi, du, dv, dteta, dq, w, &
156                    dufi, dvfi, dtetafi, dqfi, dpfi)                    dufi, dvfi, dtetafi, dqfi, dpfi)
157    
158               ! ajout des tendances physiques:               ! ajout des tendances physiques:
159               CALL addfi(nqmx, dtphys, &               CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &
160                    ucov, vcov, teta, q, ps, &                    dtetafi, dqfi, dpfi)
                   dufi, dvfi, dtetafi, dqfi, dpfi)  
   
              ! Diagnostique de conservation de l'énergie : difference  
              IF (ip_ebil_dyn >= 1) THEN  
                 ztit = 'bil phys'  
                 CALL diagedyn(ztit, 2, 1, 1, dtphys, ucov, vcov, ps, p3d, pk, &  
                      teta, q(:, :, 1), q(:, :, 2))  
              ENDIF  
161            ENDIF            ENDIF
162    
163            CALL pression(ip1jmp1, ap, bp, ps, p3d)            CALL pression(ip1jmp1, ap, bp, ps, p3d)
164            CALL exner_hyb(ps, p3d, pks, pk, pkf)            CALL exner_hyb(ps, p3d, pks, pk, pkf)
165    
166            ! dissipation horizontale et verticale des petites echelles:            IF (MOD(itau + 1, idissip) == 0) THEN
167                 ! dissipation horizontale et verticale des petites echelles:
168    
           IF (apdiss) THEN  
169               ! calcul de l'energie cinetique avant dissipation               ! calcul de l'energie cinetique avant dissipation
170               call covcont(llm, ucov, vcov, ucont, vcont)               call covcont(llm, ucov, vcov, ucont, vcont)
171               call enercin(vcov, ucov, vcont, ucont, ecin0)               call enercin(vcov, ucov, vcont, ucont, ecin0)
# Line 248  contains Line 175  contains
175               ucov=ucov + dudis               ucov=ucov + dudis
176               vcov=vcov + dvdis               vcov=vcov + dvdis
177    
178               if (dissip_conservative) then               ! On rajoute la tendance due à la transformation Ec -> E
179                  ! On rajoute la tendance due a la transform. Ec -> E               ! thermique créée lors de la dissipation
180                  ! therm. cree lors de la dissipation               call covcont(llm, ucov, vcov, ucont, vcont)
181                  call covcont(llm, ucov, vcov, ucont, vcont)               call enercin(vcov, ucov, vcont, ucont, ecin)
182                  call enercin(vcov, ucov, vcont, ucont, ecin)               dtetaecdt= (ecin0 - ecin) / pk
183                  dtetaecdt= (ecin0 - ecin) / pk               dtetadis=dtetadis + dtetaecdt
                 dtetadis=dtetadis + dtetaecdt  
              endif  
184               teta=teta + dtetadis               teta=teta + dtetadis
185    
186               ! Calcul de la valeur moyenne, unique de h aux poles .....               ! Calcul de la valeur moyenne unique de h aux pôles
187                 forall (l = 1: llm)
188               DO l = 1, llm                  teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
189                  DO ij = 1, iim                       / apoln
190                     tppn(ij) = aire(ij) * teta(ij, l)                  teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
191                     tpps(ij) = aire(ij + ip1jm) * teta(ij + ip1jm, l)                       * teta(:iim, jjm + 1, l)) / apols
192                  ENDDO               END forall
193                  tpn = SSUM(iim, tppn, 1) / apoln  
194                  tps = SSUM(iim, tpps, 1) / apols               ps(:, 1) = SUM(aire_2d(:iim, 1) * ps(:iim, 1)) / apoln
195                 ps(:, jjm + 1) = SUM(aire_2d(:iim, jjm+1) * ps(:iim, jjm + 1)) &
196                  DO ij = 1, iip1                    / apols
                    teta(ij, l) = tpn  
                    teta(ij + ip1jm, l) = tps  
                 ENDDO  
              ENDDO  
   
              DO ij = 1, iim  
                 tppn(ij) = aire(ij) * ps(ij)  
                 tpps(ij) = aire(ij + ip1jm) * ps(ij + ip1jm)  
              ENDDO  
              tpn = SSUM(iim, tppn, 1) / apoln  
              tps = SSUM(iim, tpps, 1) / apols  
   
              DO ij = 1, iip1  
                 ps(ij) = tpn  
                 ps(ij + ip1jm) = tps  
              ENDDO  
   
197            END IF            END IF
198    
199            ! fin de l'intégration dynamique et physique pour le pas "itau"            ! fin de l'intégration dynamique et physique pour le pas "itau"
# Line 303  contains Line 211  contains
211               ENDIF               ENDIF
212            ENDIF            ENDIF
213    
214            IF (itau == itaufinp1) exit outer_loop            IF (itau == itaufin + 1) exit outer_loop
   
           ! ecriture du fichier histoire moyenne:  
215    
           ! Comment out the following calls when you do not want the output  
           ! files "dyn_hist_ave.nc" and "dynzon.nc"  
216            IF (MOD(itau, iperiod) == 0 .OR. itau == itaufin) THEN            IF (MOD(itau, iperiod) == 0 .OR. itau == itaufin) THEN
217                 ! ecriture du fichier histoire moyenne:
218               CALL writedynav(histaveid, nqmx, itau, vcov, &               CALL writedynav(histaveid, nqmx, itau, vcov, &
219                    ucov, teta, pk, phi, q, masse, ps, phis)                    ucov, teta, pk, phi, q, masse, ps, phis)
220               call bilan_dyn(2, dtvr * iperiod, dtvr * day_step * periodav, &               call bilan_dyn(2, dtvr * iperiod, dtvr * day_step * periodav, &
221                    ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, q)                    ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, q)
222            ENDIF            ENDIF
223    
224            IF (itau == itaufin) THEN            IF (itau == itaufin) CALL dynredem1("restart.nc", vcov, ucov, teta, &
225               CALL dynredem1("restart.nc", 0., vcov, ucov, teta, q, masse, ps)                 q, masse, ps, itau=itau_dyn+itaufin)
              CLOSE(99)  
           ENDIF  
226    
227            ! gestion de l'integration temporelle:            ! gestion de l'integration temporelle:
   
228            IF (MOD(itau, iperiod) == 0) exit            IF (MOD(itau, iperiod) == 0) exit
229            IF (MOD(itau - 1, iperiod) == 0) THEN            IF (MOD(itau - 1, iperiod) == 0) THEN
230               IF (forward) THEN               IF (forward) THEN
# Line 335  contains Line 237  contains
237                  dt = 2. * dtvr                  dt = 2. * dtvr
238               END IF               END IF
239            ELSE            ELSE
240               ! ...... pas leapfrog .....               ! pas leapfrog
241               leapf = .TRUE.               leapf = .TRUE.
242               dt = 2. * dtvr               dt = 2. * dtvr
243            END IF            END IF

Legend:
Removed from v.18  
changed lines
  Added in v.29

  ViewVC Help
Powered by ViewVC 1.1.21