source: LMDZ6/trunk/libf/dynphy_lonlat/calfis.F90 @ 5248

Last change on this file since 5248 was 5246, checked in by abarral, 23 hours ago

Convert fixed-form to free-form sources .F -> .{f,F}90
(WIP: some .F remain, will be handled in subsequent commits)

  • 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: 19.7 KB
RevLine 
[524]1!
[1279]2! $Id: calfis.F90 5246 2024-10-21 12:58:45Z abarral $
[524]3!
[5246]4!
5!
6SUBROUTINE calfis(lafin, &
7        jD_cur, jH_cur, &
8        pucov, &
9        pvcov, &
10        pteta, &
11        pq, &
12        pmasse, &
13        pps, &
14        pp, &
15        ppk, &
16        pphis, &
17        pphi, &
18        pducov, &
19        pdvcov, &
20        pdteta, &
21        pdq, &
22        flxw, &
23        pdufi, &
24        pdvfi, &
25        pdhfi, &
26        pdqfi, &
27        pdpsfi)
28  !
29  !    Auteur :  P. Le Van, F. Hourdin
30  !   .........
31  USE infotrac, ONLY: nqtot, tracers
32  USE control_mod, ONLY: planet_type, nsplit_phys
[2418]33#ifdef CPP_PHYS
[5246]34  USE callphysiq_mod, ONLY: call_physiq
35#endif
36  USE comconst_mod, ONLY: cpp, daysec, dtphys, dtvr, kappa, pi
37  USE comvert_mod, ONLY: preff, presnivs
[524]38
[5246]39  IMPLICIT NONE
40  !=======================================================================
41  !
42  !   1. rearrangement des tableaux et transformation
43  !  variables dynamiques  >  variables physiques
44  !   2. calcul des termes physiques
45  !   3. retransformation des tendances physiques en tendances dynamiques
46  !
47  !   remarques:
48  !   ----------
49  !
50  !    - les vents sont donnes dans la physique par leurs composantes
51  !  naturelles.
52  !    - la variable thermodynamique de la physique est une variable
53  !  intensive :   T
54  !  pour la dynamique on prend    T * ( preff / p(l) ) **kappa
55  !    - les deux seules variables dependant de la geometrie necessaires
56  !  pour la physique sont la latitude pour le rayonnement et
57  !  l'aire de la maille quand on veut integrer une grandeur
58  !  horizontalement.
59  !    - les points de la physique sont les points scalaires de la
60  !  la dynamique; numerotation:
61  !      1 pour le pole nord
62  !      (jjm-1)*iim pour l'interieur du domaine
63  !      ngridmx pour le pole sud
64  !  ---> ngridmx=2+(jjm-1)*iim
65  !
66  ! Input :
67  ! -------
68  !   pucov           covariant zonal velocity
69  !   pvcov           covariant meridional velocity
70  !   pteta           potential temperature
71  !   pps             surface pressure
72  !   pmasse          masse d'air dans chaque maille
73  !   pts             surface temperature  (K)
74  !   callrad         clef d'appel au rayonnement
75  !
76  !    Output :
77  !    --------
78  !    pdufi          tendency for the natural zonal velocity (ms-1)
79  !    pdvfi          tendency for the natural meridional velocity
80  !    pdhfi          tendency for the potential temperature
81  !    pdtsfi         tendency for the surface temperature
82  !
83  !    pdtrad         radiative tendencies  \  both input
84  !    pfluxrad       radiative fluxes      /  and output
85  !
86  !=======================================================================
87  !
88  !-----------------------------------------------------------------------
89  !
90  !    0.  Declarations :
91  !    ------------------
[524]92
[5246]93  include "dimensions.h"
94  include "paramet.h"
[524]95
[5246]96  INTEGER :: ngridmx
97  PARAMETER( ngridmx = 2+(jjm-1)*iim - 1/jjm   )
[524]98
[5246]99  include "comgeom2.h"
100  include "iniprint.h"
[524]101
[5246]102  !    Arguments :
103  !    -----------
104  LOGICAL,INTENT(IN) ::  lafin ! .true. for the very last call to physics
105  REAL,INTENT(IN):: jD_cur, jH_cur
106  REAL,INTENT(IN) :: pvcov(iip1,jjm,llm) ! covariant meridional velocity
107  REAL,INTENT(IN) :: pucov(iip1,jjp1,llm) ! covariant zonal velocity
108  REAL,INTENT(IN) :: pteta(iip1,jjp1,llm) ! potential temperature
109  REAL,INTENT(IN) :: pmasse(iip1,jjp1,llm) ! mass in each cell ! not used
110  REAL,INTENT(IN) :: pq(iip1,jjp1,llm,nqtot) ! tracers
111  REAL,INTENT(IN) :: pphis(iip1,jjp1) ! surface geopotential
112  REAL,INTENT(IN) :: pphi(iip1,jjp1,llm) ! geopotential
[1279]113
[5246]114  REAL,INTENT(IN) :: pdvcov(iip1,jjm,llm) ! dynamical tendency on vcov
115  REAL,INTENT(IN) :: pducov(iip1,jjp1,llm) ! dynamical tendency on ucov
116  REAL,INTENT(IN) :: pdteta(iip1,jjp1,llm) ! dynamical tendency on teta
117  ! ! NB: pdteta is used only to compute pcvgt which is in fact not used...
118  REAL,INTENT(IN) :: pdq(iip1,jjp1,llm,nqtot) ! dynamical tendency on tracers
119  ! ! NB: pdq is only used to compute pcvgq which is in fact not used...
[524]120
[5246]121  REAL,INTENT(IN) :: pps(iip1,jjp1) ! surface pressure (Pa)
122  REAL,INTENT(IN) :: pp(iip1,jjp1,llmp1) ! pressure at mesh interfaces (Pa)
123  REAL,INTENT(IN) :: ppk(iip1,jjp1,llm) ! Exner at mid-layer
124  REAL,INTENT(IN) :: flxw(iip1,jjp1,llm) ! Vertical mass flux on lower mesh interfaces (kg/s) (on llm because flxw(:,:,llm+1)=0)
[524]125
[5246]126  ! ! tendencies (in */s) from the physics
127  REAL,INTENT(OUT) :: pdvfi(iip1,jjm,llm) ! tendency on covariant meridional wind
128  REAL,INTENT(OUT) :: pdufi(iip1,jjp1,llm) ! tendency on covariant zonal wind
129  REAL,INTENT(OUT) :: pdhfi(iip1,jjp1,llm) ! tendency on potential temperature (K/s)
130  REAL,INTENT(OUT) :: pdqfi(iip1,jjp1,llm,nqtot) ! tendency on tracers
131  REAL,INTENT(OUT) :: pdpsfi(iip1,jjp1) ! tendency on surface pressure (Pa/s)
[524]132
133
[5246]134  !    Local variables :
135  !    -----------------
[1403]136
[5246]137  INTEGER :: i,j,l,ig0,ig,iq,itr
138  REAL :: zpsrf(ngridmx)
139  REAL :: zplev(ngridmx,llm+1),zplay(ngridmx,llm)
140  REAL :: zphi(ngridmx,llm),zphis(ngridmx)
141  !
142  REAL :: zrot(iip1,jjm,llm) ! AdlC May 2014
143  REAL :: zufi(ngridmx,llm), zvfi(ngridmx,llm)
144  REAL :: zrfi(ngridmx,llm) ! relative wind vorticity
145  REAL :: ztfi(ngridmx,llm),zqfi(ngridmx,llm,nqtot)
146  REAL :: zpk(ngridmx,llm)
147  !
148  REAL :: pcvgu(ngridmx,llm), pcvgv(ngridmx,llm)
149  REAL :: pcvgt(ngridmx,llm), pcvgq(ngridmx,llm,2)
150  !
151  REAL :: zdufi(ngridmx,llm),zdvfi(ngridmx,llm)
152  REAL :: zdtfi(ngridmx,llm),zdqfi(ngridmx,llm,nqtot)
153  REAL :: zdpsrf(ngridmx)
154  !
155  REAL :: zdufic(ngridmx,llm),zdvfic(ngridmx,llm)
156  REAL :: zdtfic(ngridmx,llm),zdqfic(ngridmx,llm,nqtot)
157  REAL :: jH_cur_split,zdt_split
158  LOGICAL :: debut_split,lafin_split
159  INTEGER :: isplit
[524]160
[5246]161  REAL :: zsin(iim),zcos(iim),z1(iim)
162  REAL :: zsinbis(iim),zcosbis(iim),z1bis(iim)
163  REAL :: unskap, pksurcp
164  !
165  REAL :: flxwfi(ngridmx,llm)  ! Flux de masse verticale sur la grille physiq
166  !
[1654]167
[5246]168  REAL :: SSUM
[524]169
[5246]170  LOGICAL,SAVE :: firstcal=.true., debut=.true.
171   ! REAL rdayvrai
[524]172
[5246]173  !
174  !-----------------------------------------------------------------------
175  !
176  !    1. Initialisations :
177  !    --------------------
178  !
179  !
180  IF ( firstcal )  THEN
181    debut = .TRUE.
182    IF (ngridmx.NE.2+(jjm-1)*iim) THEN
183     write(lunout,*) 'STOP dans calfis'
184     write(lunout,*) &
185           'La dimension ngridmx doit etre egale a 2 + (jjm-1)*iim'
186     write(lunout,*) '  ngridmx  jjm   iim   '
187     write(lunout,*) ngridmx,jjm,iim
188     call abort_gcm("calfis", "", 1)
189    ENDIF
190  ELSE
191    debut = .FALSE.
192  ENDIF ! of IF (firstcal)
[524]193
[5246]194  !
195  !
196  !-----------------------------------------------------------------------
197  !   40. transformation des variables dynamiques en variables physiques:
198  !   ---------------------------------------------------------------
[524]199
[5246]200  !   41. pressions au sol (en Pascals)
201  !   ----------------------------------
[524]202
203
[5246]204  zpsrf(1) = pps(1,1)
[524]205
[5246]206  ig0  = 2
207  DO j = 2,jjm
208     CALL SCOPY( iim,pps(1,j),1,zpsrf(ig0), 1 )
209     ig0 = ig0+iim
210  ENDDO
[524]211
[5246]212  zpsrf(ngridmx) = pps(1,jjp1)
213
214
215  !   42. pression intercouches et fonction d'Exner:
216  !
217  !   -----------------------------------------------------------------
218  ! .... zplev  definis aux (llm +1) interfaces des couches  ....
219  ! .... zplay  definis aux (  llm )    milieux des couches  ....
220  !   -----------------------------------------------------------------
221
222  !    ...    Exner = cp * ( p(l) / preff ) ** kappa     ....
223  !
224   unskap   = 1./ kappa
225  !
226  DO l = 1, llm
227    zpk(   1,l ) = ppk(1,1,l)
228    zplev( 1,l ) = pp(1,1,l)
229    ig0 = 2
230      DO j = 2, jjm
231         DO i =1, iim
232          zpk(   ig0,l ) = ppk(i,j,l)
233          zplev( ig0,l ) = pp(i,j,l)
234          ig0 = ig0 +1
235         ENDDO
[524]236      ENDDO
[5246]237    zpk(   ngridmx,l ) = ppk(1,jjp1,l)
238    zplev( ngridmx,l ) = pp(1,jjp1,l)
239  ENDDO
240    zplev( 1,llmp1 ) = pp(1,1,llmp1)
241    ig0 = 2
242      DO j = 2, jjm
243         DO i =1, iim
244          zplev( ig0,llmp1 ) = pp(i,j,llmp1)
245          ig0 = ig0 +1
246         ENDDO
247      ENDDO
248    zplev( ngridmx,llmp1 ) = pp(1,jjp1,llmp1)
249  !
250  !
[524]251
[5246]252  !   43. temperature naturelle (en K) et pressions milieux couches .
253  !   ---------------------------------------------------------------
[524]254
[5246]255  DO l=1,llm
[524]256
[5246]257     pksurcp     =  ppk(1,1,l) / cpp
258     zplay(1,l)  =  preff * pksurcp ** unskap
259     ztfi(1,l)   =  pteta(1,1,l) *  pksurcp
260     pcvgt(1,l)  =  pdteta(1,1,l) * pksurcp / pmasse(1,1,l)
261     ig0         = 2
[524]262
[5246]263     DO j = 2, jjm
264        DO i = 1, iim
265          pksurcp        = ppk(i,j,l) / cpp
266          zplay(ig0,l)   = preff * pksurcp ** unskap
267          ztfi(ig0,l)    = pteta(i,j,l)  * pksurcp
268          pcvgt(ig0,l)   = pdteta(i,j,l) * pksurcp / pmasse(i,j,l)
269          ig0            = ig0 + 1
270        ENDDO
271     ENDDO
[524]272
[5246]273     pksurcp       = ppk(1,jjp1,l) / cpp
274     zplay(ig0,l)  = preff * pksurcp ** unskap
275     ztfi (ig0,l)  = pteta(1,jjp1,l)  * pksurcp
276     pcvgt(ig0,l)  = pdteta(1,jjp1,l) * pksurcp/ pmasse(1,jjp1,l)
[524]277
[5246]278  ENDDO
[524]279
[5246]280  !   43.bis traceurs
281  !   ---------------
282  !
283  itr=0
284  DO iq=1,nqtot
285     IF(.NOT.tracers(iq)%isAdvected) CYCLE
286     itr = itr + 1
287     DO l=1,llm
288        zqfi(1,l,itr) = pq(1,1,l,iq)
289        ig0           = 2
290        DO j=2,jjm
291           DO i = 1, iim
292              zqfi(ig0,l,itr) = pq(i,j,l,iq)
293              ig0             = ig0 + 1
294           ENDDO
295        ENDDO
296        zqfi(ig0,l,itr) = pq(1,jjp1,l,iq)
297     ENDDO
298  ENDDO
[524]299
[5246]300  !   convergence dynamique pour les traceurs "EAU"
301  ! Earth-specific treatment of first 2 tracers (water)
302   if (planet_type=="earth") then
303    DO iq=1,2
304     DO l=1,llm
305        pcvgq(1,l,iq)= pdq(1,1,l,iq) / pmasse(1,1,l)
306        ig0          = 2
307        DO j=2,jjm
308           DO i = 1, iim
309              pcvgq(ig0,l,iq) = pdq(i,j,l,iq) / pmasse(i,j,l)
310              ig0             = ig0 + 1
311           ENDDO
[1279]312        ENDDO
[5246]313        pcvgq(ig0,l,iq)= pdq(1,jjp1,l,iq) / pmasse(1,jjp1,l)
314     ENDDO
315    ENDDO
316   endif ! of if (planet_type=="earth")
[524]317
318
[5246]319  !   Geopotentiel calcule par rapport a la surface locale:
320  !   -----------------------------------------------------
[524]321
[5246]322  CALL gr_dyn_fi(llm,iip1,jjp1,ngridmx,pphi,zphi)
323  CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,pphis,zphis)
324  DO l=1,llm
325     DO ig=1,ngridmx
326       zphi(ig,l)=zphi(ig,l)-zphis(ig)
327     ENDDO
328  ENDDO
[524]329
[5246]330  !   ....  Calcul de la vitesse  verticale  ( en Pa*m*s  ou Kg/s )  ....
331  ! JG : ancien calcule de omega utilise dans physiq.F. Maintenant le flux
332  !    de masse est calclue dans advtrac.F
333   ! DO l=1,llm
334   !   pvervel(1,l)=pw(1,1,l) * g /apoln
335   !   ig0=2
336   !  DO j=2,jjm
337   !      DO i = 1, iim
338   !         pvervel(ig0,l) = pw(i,j,l) * g * unsaire(i,j)
339   !         ig0 = ig0 + 1
340   !      ENDDO
341  !       ENDDO
342   !   pvervel(ig0,l)=pw(1,jjp1,l) * g /apols
343   ! ENDDO
[524]344
[5246]345  !
346  !   45. champ u:
347  !   ------------
[524]348
[5246]349  DO l=1,llm
[524]350
[5246]351     DO j=2,jjm
352        ig0 = 1+(j-2)*iim
353        zufi(ig0+1,l)= 0.5 * &
354              ( pucov(iim,j,l)/cu(iim,j) + pucov(1,j,l)/cu(1,j) )
355        pcvgu(ig0+1,l)= 0.5 * &
356              ( pducov(iim,j,l)/cu(iim,j) + pducov(1,j,l)/cu(1,j) )
357        DO i=2,iim
358           zufi(ig0+i,l)= 0.5 * &
359                 ( pucov(i-1,j,l)/cu(i-1,j) + pucov(i,j,l)/cu(i,j) )
360           pcvgu(ig0+i,l)= 0.5 * &
361                 ( pducov(i-1,j,l)/cu(i-1,j) + pducov(i,j,l)/cu(i,j) )
362        END DO
363     END DO
[524]364
[5246]365  END DO
[524]366
367
[5246]368  !  Alvaro de la Camara (May 2014)
369  !  46.1 Calcul de la vorticite et passage sur la grille physique
370  !  --------------------------------------------------------------
371  DO l=1,llm
372    do i=1,iim
373      do j=1,jjm
374        zrot(i,j,l) = (pvcov(i+1,j,l) - pvcov(i,j,l) &
375              + pucov(i,j+1,l) - pucov(i,j,l)) &
376              / (cu(i,j)+cu(i,j+1)) &
377              / (cv(i+1,j)+cv(i,j)) *4
378      enddo
379    enddo
380  ENDDO
[2333]381
[5246]382  !   46.champ v:
383  !   -----------
[524]384
[5246]385  DO l=1,llm
386     DO j=2,jjm
387        ig0=1+(j-2)*iim
388        DO i=1,iim
389           zvfi(ig0+i,l)= 0.5 * &
390                 ( pvcov(i,j-1,l)/cv(i,j-1) + pvcov(i,j,l)/cv(i,j) )
391           pcvgv(ig0+i,l)= 0.5 * &
392                 ( pdvcov(i,j-1,l)/cv(i,j-1) + pdvcov(i,j,l)/cv(i,j) )
393        ENDDO
394           zrfi(ig0 + 1,l)= 0.25 *(zrot(iim,j-1,l)+zrot(iim,j,l) &
395                 +zrot(1,j-1,l)+zrot(1,j,l))
396        DO i=2,iim
397           zrfi(ig0 + i,l)= 0.25 *(zrot(i-1,j-1,l)+zrot(i-1,j,l) &
398                 +zrot(i,j-1,l)+zrot(i,j,l))   !  AdlC MAY 2014
399        ENDDO
400     ENDDO
401  ENDDO
[524]402
403
[5246]404  !   47. champs de vents aux pole nord
405  !   ------------------------------
406     ! U = 1 / pi  *  integrale [ v * cos(long) * d long ]
407     ! V = 1 / pi  *  integrale [ v * sin(long) * d long ]
[524]408
[5246]409  DO l=1,llm
[524]410
[5246]411     z1(1)   =(rlonu(1)-rlonu(iim)+2.*pi)*pvcov(1,1,l)/cv(1,1)
412     z1bis(1)=(rlonu(1)-rlonu(iim)+2.*pi)*pdvcov(1,1,l)/cv(1,1)
413     DO i=2,iim
414        z1(i)   =(rlonu(i)-rlonu(i-1))*pvcov(i,1,l)/cv(i,1)
415        z1bis(i)=(rlonu(i)-rlonu(i-1))*pdvcov(i,1,l)/cv(i,1)
416     ENDDO
[524]417
[5246]418     DO i=1,iim
419        zcos(i)   = COS(rlonv(i))*z1(i)
420        zcosbis(i)= COS(rlonv(i))*z1bis(i)
421        zsin(i)   = SIN(rlonv(i))*z1(i)
422        zsinbis(i)= SIN(rlonv(i))*z1bis(i)
423     ENDDO
[524]424
[5246]425     zufi(1,l)  = SSUM(iim,zcos,1)/pi
426     pcvgu(1,l) = SSUM(iim,zcosbis,1)/pi
427     zvfi(1,l)  = SSUM(iim,zsin,1)/pi
428     pcvgv(1,l) = SSUM(iim,zsinbis,1)/pi
429     zrfi(1, l) = 0.
430  ENDDO
[524]431
432
[5246]433  !   48. champs de vents aux pole sud:
434  !   ---------------------------------
435     ! U = 1 / pi  *  integrale [ v * cos(long) * d long ]
436     ! V = 1 / pi  *  integrale [ v * sin(long) * d long ]
[524]437
[5246]438  DO l=1,llm
[524]439
[5246]440     z1(1)   =(rlonu(1)-rlonu(iim)+2.*pi)*pvcov(1,jjm,l)/cv(1,jjm)
441     z1bis(1)=(rlonu(1)-rlonu(iim)+2.*pi)*pdvcov(1,jjm,l)/cv(1,jjm)
442     DO i=2,iim
443        z1(i)   =(rlonu(i)-rlonu(i-1))*pvcov(i,jjm,l)/cv(i,jjm)
444        z1bis(i)=(rlonu(i)-rlonu(i-1))*pdvcov(i,jjm,l)/cv(i,jjm)
445     ENDDO
[524]446
[5246]447     DO i=1,iim
448        zcos(i)    = COS(rlonv(i))*z1(i)
449        zcosbis(i) = COS(rlonv(i))*z1bis(i)
450        zsin(i)    = SIN(rlonv(i))*z1(i)
451        zsinbis(i) = SIN(rlonv(i))*z1bis(i)
452     ENDDO
[524]453
[5246]454     zufi(ngridmx,l)  = SSUM(iim,zcos,1)/pi
455     pcvgu(ngridmx,l) = SSUM(iim,zcosbis,1)/pi
456     zvfi(ngridmx,l)  = SSUM(iim,zsin,1)/pi
457     pcvgv(ngridmx,l) = SSUM(iim,zsinbis,1)/pi
458     zrfi(ngridmx, l) = 0.
459  ENDDO
460  !
461  ! On change de grille, dynamique vers physiq, pour le flux de masse verticale
462  CALL gr_dyn_fi(llm,iip1,jjp1,ngridmx,flxw,flxwfi)
[524]463
[5246]464  !-----------------------------------------------------------------------
465  !   Appel de la physique:
466  !   ---------------------
[524]467
468
[1403]469
[5246]470   ! write(lunout,*) 'PHYSIQUE AVEC NSPLIT_PHYS=',nsplit_phys
471  zdt_split=dtphys/nsplit_phys
472  zdufic(:,:)=0.
473  zdvfic(:,:)=0.
474  zdtfic(:,:)=0.
475  zdqfic(:,:,:)=0.
[1403]476
[1615]477#ifdef CPP_PHYS
[1403]478
[5246]479   do isplit=1,nsplit_phys
[1615]480
[5246]481     jH_cur_split=jH_cur+(isplit-1) * dtvr / (daysec *nsplit_phys)
482     debut_split=debut.and.isplit==1
483     lafin_split=lafin.and.isplit==nsplit_phys
[1403]484
[5246]485   ! if (planet_type=="earth") then
486    CALL call_physiq(ngridmx,llm,nqtot,tracers(:)%name, &
487          debut_split,lafin_split, &
488          jD_cur,jH_cur_split,zdt_split, &
489          zplev,zplay, &
490          zpk,zphi,zphis, &
491          presnivs, &
492          zufi,zvfi,zrfi,ztfi,zqfi, &
493          flxwfi,pducov, &
494          zdufi,zdvfi,zdtfi,zdqfi,zdpsrf)
495  !
496  !  else if ( planet_type=="generic" ) then
497  !
498  !     CALL physiq (ngridmx,     !! ngrid
499  ! .             llm,            !! nlayer
500  ! .             nqtot,          !! nq
501  ! .             tracers(:)%name,!! tracer names from dynamical core (given in infotrac)
502  ! .             debut_split,    !! firstcall
503  ! .             lafin_split,    !! lastcall
504  ! .             jD_cur,         !! pday. see leapfrog
505  ! .             jH_cur_split,   !! ptime "fraction of day"
506  ! .             zdt_split,      !! ptimestep
507  ! .             zplev,          !! pplev
508  ! .             zplay,          !! pplay
509  ! .             zphi,           !! pphi
510  ! .             zufi,           !! pu
511  ! .             zvfi,           !! pv
512  ! .             ztfi,           !! pt
513  ! .             zqfi,           !! pq
514  ! .             flxwfi,         !! pw !! or 0. anyway this is for diagnostic. not used in physiq.
515  ! .             zdufi,          !! pdu
516  ! .             zdvfi,          !! pdv
517  ! .             zdtfi,          !! pdt
518  ! .             zdqfi,          !! pdq
519  ! .             zdpsrf,         !! pdpsrf
520  ! .             tracerdyn)      !! tracerdyn <-- utilite ???
521  !
522  !  endif ! of if (planet_type=="earth")
[1654]523
[5246]524     zufi(:,:)=zufi(:,:)+zdufi(:,:)*zdt_split
525     zvfi(:,:)=zvfi(:,:)+zdvfi(:,:)*zdt_split
526     ztfi(:,:)=ztfi(:,:)+zdtfi(:,:)*zdt_split
527     zqfi(:,:,:)=zqfi(:,:,:)+zdqfi(:,:,:)*zdt_split
[1403]528
[5246]529     zdufic(:,:)=zdufic(:,:)+zdufi(:,:)
530     zdvfic(:,:)=zdvfic(:,:)+zdvfi(:,:)
531     zdtfic(:,:)=zdtfic(:,:)+zdtfi(:,:)
532     zdqfic(:,:,:)=zdqfic(:,:,:)+zdqfi(:,:,:)
[1403]533
[5246]534   enddo ! of do isplit=1,nsplit_phys
[1615]535
536#endif
[5246]537  ! of #ifdef CPP_PHYS
[1615]538
[5246]539  zdufi(:,:)=zdufic(:,:)/nsplit_phys
540  zdvfi(:,:)=zdvfic(:,:)/nsplit_phys
541  zdtfi(:,:)=zdtfic(:,:)/nsplit_phys
542  zdqfi(:,:,:)=zdqfic(:,:,:)/nsplit_phys
[1403]543
[524]544
545500   CONTINUE
546
[5246]547  !-----------------------------------------------------------------------
548  !   transformation des tendances physiques en tendances dynamiques:
549  !   ---------------------------------------------------------------
[524]550
[5246]551  !  tendance sur la pression :
552  !  -----------------------------------
[524]553
[5246]554  CALL gr_fi_dyn(1,ngridmx,iip1,jjp1,zdpsrf,pdpsfi)
555  !
556  !   62. enthalpie potentielle
557  !   ---------------------
[524]558
[5246]559  DO l=1,llm
[524]560
[5246]561     DO i=1,iip1
562      pdhfi(i,1,l)    = cpp *  zdtfi(1,l)      / ppk(i, 1  ,l)
563      pdhfi(i,jjp1,l) = cpp *  zdtfi(ngridmx,l)/ ppk(i,jjp1,l)
564     ENDDO
[524]565
[5246]566     DO j=2,jjm
567        ig0=1+(j-2)*iim
568        DO i=1,iim
569           pdhfi(i,j,l) = cpp * zdtfi(ig0+i,l) / ppk(i,j,l)
570        ENDDO
571           pdhfi(iip1,j,l) =  pdhfi(1,j,l)
572     ENDDO
[524]573
[5246]574  ENDDO
[524]575
576
[5246]577  !   62. humidite specifique
578  !   ---------------------
579  ! Ehouarn: removed this useless bit: was overwritten at step 63 anyways
580   ! DO iq=1,nqtot
581   !    DO l=1,llm
582   !       DO i=1,iip1
583   !          pdqfi(i,1,l,iq)    = zdqfi(1,l,iq)
584   !          pdqfi(i,jjp1,l,iq) = zdqfi(ngridmx,l,iq)
585   !       ENDDO
586   !       DO j=2,jjm
587   !          ig0=1+(j-2)*iim
588   !          DO i=1,iim
589   !             pdqfi(i,j,l,iq) = zdqfi(ig0+i,l,iq)
590   !          ENDDO
591   !          pdqfi(iip1,j,l,iq) = pdqfi(1,j,l,iq)
592   !       ENDDO
593   !    ENDDO
594   ! ENDDO
[524]595
[5246]596  !   63. traceurs
597  !   ------------
598  ! initialisation des tendances
599  pdqfi(:,:,:,:)=0.
600  !
601  itr = 0
602  DO iq=1,nqtot
603     IF(.NOT.tracers(iq)%isAdvected) CYCLE
604     itr = itr + 1
605     DO l=1,llm
606        DO i=1,iip1
607           pdqfi(i,1,l,iq)    = zdqfi(1,l,itr)
608           pdqfi(i,jjp1,l,iq) = zdqfi(ngridmx,l,itr)
609        ENDDO
610        DO j=2,jjm
611           ig0=1+(j-2)*iim
612           DO i=1,iim
613              pdqfi(i,j,l,iq) = zdqfi(ig0+i,l,itr)
614           ENDDO
615           pdqfi(iip1,j,l,iq) = pdqfi(1,j,l,itr)
616        ENDDO
617     ENDDO
618  ENDDO
[524]619
[5246]620  !   65. champ u:
621  !   ------------
[524]622
[5246]623  DO l=1,llm
[524]624
[5246]625     DO i=1,iip1
626        pdufi(i,1,l)    = 0.
627        pdufi(i,jjp1,l) = 0.
628     ENDDO
[524]629
[5246]630     DO j=2,jjm
631        ig0=1+(j-2)*iim
632        DO i=1,iim-1
633           pdufi(i,j,l)= &
634                 0.5*(zdufi(ig0+i,l)+zdufi(ig0+i+1,l))*cu(i,j)
635        ENDDO
636        pdufi(iim,j,l)= &
637              0.5*(zdufi(ig0+1,l)+zdufi(ig0+iim,l))*cu(iim,j)
638        pdufi(iip1,j,l)=pdufi(1,j,l)
639     ENDDO
[524]640
[5246]641  ENDDO
[524]642
643
[5246]644  !   67. champ v:
645  !   ------------
[524]646
[5246]647  DO l=1,llm
[524]648
[5246]649     DO j=2,jjm-1
650        ig0=1+(j-2)*iim
651        DO i=1,iim
652           pdvfi(i,j,l)= &
653                 0.5*(zdvfi(ig0+i,l)+zdvfi(ig0+i+iim,l))*cv(i,j)
654        ENDDO
655        pdvfi(iip1,j,l) = pdvfi(1,j,l)
656     ENDDO
657  ENDDO
[524]658
659
[5246]660  !   68. champ v pres des poles:
661  !   ---------------------------
662   ! v = U * cos(long) + V * SIN(long)
[524]663
[5246]664  DO l=1,llm
[524]665
[5246]666     DO i=1,iim
667        pdvfi(i,1,l)= &
668              zdufi(1,l)*COS(rlonv(i))+zdvfi(1,l)*SIN(rlonv(i))
669        pdvfi(i,jjm,l)=zdufi(ngridmx,l)*COS(rlonv(i)) &
670              +zdvfi(ngridmx,l)*SIN(rlonv(i))
671        pdvfi(i,1,l)= &
672              0.5*(pdvfi(i,1,l)+zdvfi(i+1,l))*cv(i,1)
673        pdvfi(i,jjm,l)= &
674              0.5*(pdvfi(i,jjm,l)+zdvfi(ngridmx-iip1+i,l))*cv(i,jjm)
675      ENDDO
[524]676
[5246]677     pdvfi(iip1,1,l)  = pdvfi(1,1,l)
678     pdvfi(iip1,jjm,l)= pdvfi(1,jjm,l)
[524]679
[5246]680  ENDDO
[524]681
[5246]682  !-----------------------------------------------------------------------
[524]683
684700   CONTINUE
685
[5246]686  firstcal = .FALSE.
687
688  RETURN
689END SUBROUTINE calfis
Note: See TracBrowser for help on using the repository browser.