source: trunk/LMDZ.TITAN/libf/phytitan/physiq.F @ 251

Last change on this file since 251 was 175, checked in by slebonnois, 14 years ago

S.LEBONNOIS:

  • Revision majeure de la physique Titan => ajout des nuages version 10 bins (Jeremie Burgalat) Cette version reste a tester mais avec clouds=0, on reste sur l'ancienne.
  • Quelques ajouts dans la doc.
File size: 51.6 KB
Line 
1!
2! $Header: /home/cvsroot/LMDZ4/libf/phylmd/physiq.F,v 1.8 2005/02/24 09:58:18 fairhead Exp $
3!
4c
5      SUBROUTINE physiq (nlon,nlev,nqmax,
6     .            debut,lafin,rjourvrai,gmtime,pdtphys,
7     .            paprs,pplay,ppk,pphi,pphis,presnivs,
8     .            u,v,t,qx,
9     .            omega,
10     .            d_u, d_v, d_t, d_qx, d_ps)
11
12c======================================================================
13c
14c Modifications pour la physique de Titan
15c  adaptation a partir de celle de Venus
16c     S. Lebonnois (LMD/CNRS) Mai 2008
17c
18c ---------------------------------------------------------------------
19c Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
20c
21c Objet: Moniteur general de la physique du modele
22cAA      Modifications quant aux traceurs :
23cAA                  -  uniformisation des parametrisations ds phytrac
24cAA                  -  stockage des moyennes des champs necessaires
25cAA                     en mode traceur off-line
26c    modif   ( P. Le Van ,  12/10/98 )
27c
28c  Arguments:
29c
30c nlon------input-I-nombre de points horizontaux
31c nlev------input-I-nombre de couches verticales
32c nqmax-----input-I-nombre de traceurs
33c debut-----input-L-variable logique indiquant le premier passage
34c lafin-----input-L-variable logique indiquant le dernier passage
35c rjourvrai-input-R-NBjours
36c gmtime----input-R-temps universel dans la journee (fraction de jour)
37c pdtphys---input-R-pas d'integration pour la physique (seconde)
38c paprs-----input-R-pression pour chaque inter-couche (en Pa)
39c pplay-----input-R-pression pour le mileu de chaque couche (en Pa)
40c ppk  -----input-R-fonction d'Exner au milieu de couche
41c pphi------input-R-geopotentiel de chaque couche (g z) (reference sol)
42c pphis-----input-R-geopotentiel du sol
43c presnivs--input_R_pressions approximat. des milieux couches ( en PA)
44c u---------input-R-vitesse dans la direction X (de O a E) en m/s
45c v---------input-R-vitesse Y (de S a N) en m/s
46c t---------input-R-temperature (K)
47c qx--------input-R-mass mixing ratio traceurs (kg/kg)
48c d_t_dyn---input-R-tendance dynamique pour "t" (K/s)
49c omega-----input-R-vitesse verticale en Pa/s
50c
51c d_u-----output-R-tendance physique de "u" (m/s/s)
52c d_v-----output-R-tendance physique de "v" (m/s/s)
53c d_t-----output-R-tendance physique de "t" (K/s)
54c d_qx----output-R-tendance physique de "qx" (kg/kg/s)
55c d_ps----output-R-tendance physique de la pression au sol
56c======================================================================
57      USE ioipsl
58      USE histcom
59      USE infotrac
60      USE control_mod
61      use dimphy
62      USE comgeomphy
63      IMPLICIT none
64c======================================================================
65c   CLEFS CPP POUR LES IO
66c   =====================
67#define histmth
68#define histday
69#define histins
70c======================================================================
71#include "dimensions.h"
72      integer jjmp1
73      parameter (jjmp1=jjm+1-1/jjm)
74#include "dimsoil.h"
75#include "clesphys.h"
76#include "temps.h"
77#include "iniprint.h"
78#include "timerad.h"
79#include "logic.h"
80#include "comorbit.h"
81#include "microtab.h"
82#include "diagmuphy.h"
83#include "itemps.h"
84c======================================================================
85      LOGICAL ok_mensuel ! sortir le fichier mensuel
86      save ok_mensuel
87c      PARAMETER (ok_mensuel=.true.)
88c
89      LOGICAL ok_journe ! sortir le fichier journalier
90      save ok_journe
91c      PARAMETER (ok_journe=.true.)
92c
93      LOGICAL ok_instan ! sortir le fichier instantane
94      save ok_instan
95c      PARAMETER (ok_instan=.true.)
96c
97c======================================================================
98c
99c Variables argument:
100c
101      INTEGER nlon
102      INTEGER nlev
103      INTEGER nqmax
104      REAL rjourvrai
105      REAL gmtime
106      REAL pdtphys
107      LOGICAL debut, lafin
108      REAL paprs(klon,klev+1)
109      REAL pplay(klon,klev)
110      REAL pphi(klon,klev)
111      REAL pphis(klon)
112      REAL presnivs(klev)
113
114! ADAPTATION GCM POUR CP(T)
115      REAL ppk(klon,klev)
116
117      REAL u(klon,klev)
118      REAL v(klon,klev)
119      REAL t(klon,klev)
120      REAL qx(klon,klev,nqmax)
121
122      REAL,save,allocatable :: t_ancien(:,:)
123      REAL,save,allocatable :: u_ancien(:,:)
124      LOGICAL ancien_ok
125      SAVE ancien_ok
126
127      REAL d_u_dyn(klon,klev)
128      REAL d_t_dyn(klon,klev)
129
130      REAL omega(klon,klev)
131
132      REAL d_u(klon,klev)
133      REAL d_v(klon,klev)
134      REAL d_t(klon,klev)
135      REAL d_qx(klon,klev,nqmax)
136      REAL d_ps(klon)
137
138      REAL,save,allocatable :: swnet(:,:)
139      REAL,save,allocatable :: lwnet(:,:)
140c
141c Variables propres a la physique
142c
143      REAL dtime
144      SAVE dtime                  ! pas temporel de la physique
145c
146      INTEGER radpas
147      SAVE radpas                 ! frequence d'appel rayonnement
148c
149      INTEGER chimpas
150      SAVE chimpas                 ! frequence d'appel chimie
151c
152      REAL,save,allocatable :: radsol(:) ! bilan radiatif au sol calcule par code radiatif
153      REAL,save,allocatable :: rlev(:,:) ! altitude a chaque niveau (interface inferieure de la couche)
154      INTEGER,save :: itap        ! compteur pour la physique
155      REAL,save,allocatable :: ftsol(:)    ! temperature du sol
156      REAL,save,allocatable :: ftsoil(:,:) ! temperature dans le sol
157      REAL,save,allocatable :: falbe(:)    ! albedo
158      REAL delp(klon,klev)        ! epaisseur d'une couche
159     
160CMODDEB FLOTT
161c
162c  Parametres de l'Orographie a l'Echelle Sous-Maille (OESM):
163c
164      REAL,save,allocatable :: zmea(:)   ! orographie moyenne
165      REAL,save,allocatable :: zstd(:)   ! deviation standard de l'OESM
166      REAL,save,allocatable :: zsig(:)   ! pente de l'OESM
167      REAL,save,allocatable :: zgam(:)   ! anisotropie de l'OESM
168      REAL,save,allocatable :: zthe(:)   ! orientation de l'OESM
169      REAL,save,allocatable :: zpic(:)   ! Maximum de l'OESM
170      REAL,save,allocatable :: zval(:)   ! Minimum de l'OESM
171      REAL,save,allocatable :: rugoro(:) ! longueur de rugosite de l'OESM
172
173      INTEGER igwd,idx(klon),itest(klon)
174c
175c  Diagnostiques 2D de drag_noro, lift_noro et gw_nonoro
176
177      REAL zulow(klon),zvlow(klon)
178      REAL zustrdr(klon), zvstrdr(klon)
179      REAL zustrli(klon), zvstrli(klon)
180      REAL zustrhi(klon), zvstrhi(klon)
181
182c Pour calcul GW drag oro et nonoro: CALCUL de N2:
183      real zdzlev(klon,klev)
184      real ztlev(klon,klev),zpklev(klon,klev)
185      real ztetalay(klon,klev),ztetalev(klon,klev)
186      real zdtetalev(klon,klev)
187      real zn2(klon,klev) ! BV^2 at plev
188
189c Pour les bilans de moment angulaire,
190      integer bilansmc
191c Pour le transport de ballons
192      integer ballons
193c j'ai aussi besoin
194c du stress de couche limite a la surface:
195
196      REAL zustrcl(klon),zvstrcl(klon)
197
198c et du stress total c de la physique:
199
200      REAL zustrph(klon),zvstrph(klon)
201c
202      REAL,save,allocatable :: zuthe(:),zvthe(:)
203
204c Variables locales:
205c
206      REAL cdragh(klon) ! drag coefficient pour T and Q
207      REAL cdragm(klon) ! drag coefficient pour vent
208c
209cAA  Pour TRACEURS
210cAA
211      REAL,save,allocatable :: source(:,:)
212      integer nmicro
213      save    nmicro
214
215      REAL ycoefh(klon,klev)    ! coef d'echange pour phytrac
216      REAL yu1(klon)            ! vents dans la premiere couche U
217      REAL yv1(klon)            ! vents dans la premiere couche V
218      character*8 nom
219      REAL qaer(klon,klev,nqmax)
220
221      REAL sens(klon), dsens(klon) ! chaleur sensible et sa derivee
222      REAL,save,allocatable :: dlw(:)  ! derivee infra rouge
223      REAL,save,allocatable :: fder(:) ! Derive de flux (sensible et latente)
224      REAL ve(klon) ! integr. verticale du transport meri. de l'energie
225      REAL vq(klon) ! integr. verticale du transport meri. de l'eau
226      REAL ue(klon) ! integr. verticale du transport zonal de l'energie
227      REAL uq(klon) ! integr. verticale du transport zonal de l'eau
228c
229
230c======================================================================
231c
232c Declaration des procedures appelees
233c
234      EXTERNAL ajsec     ! ajustement sec
235      EXTERNAL clmain    ! couche limite
236      EXTERNAL hgardfou  ! verifier les temperatures
237      EXTERNAL orbite    ! calculer l'orbite
238      EXTERNAL phyetat0  ! lire l'etat initial de la physique
239      EXTERNAL phyredem  ! ecrire l'etat de redemarrage de la physique
240      EXTERNAL radlwsw   ! rayonnements solaire et infrarouge
241      EXTERNAL suphec    ! initialiser certaines constantes
242c      EXTERNAL transp    ! transport total de l'eau et de l'energie
243      EXTERNAL abort_gcm
244      EXTERNAL printflag
245      EXTERNAL zenang
246      EXTERNAL diagetpq
247      EXTERNAL conf_phys
248      EXTERNAL diagphy
249      EXTERNAL mucorr
250      EXTERNAL phytrac
251c
252c Variables locales
253c
254CXXX PB
255      REAL fluxt(klon,klev)   ! flux turbulent de chaleur
256      REAL fluxu(klon,klev)   ! flux turbulent de vitesse u
257      REAL fluxv(klon,klev)   ! flux turbulent de vitesse v
258c
259      REAL flux_dyn(klon,klev)  ! flux de chaleur produit par la dynamique
260      REAL flux_ajs(klon,klev)  ! flux de chaleur ajustement sec
261      REAL flux_ec(klon,klev)   ! flux de chaleur Ec
262c
263c Le rayonnement n'est pas calcule tous les pas, il faut donc
264c                      sauvegarder les sorties du rayonnement
265      REAL,save,allocatable :: heat(:,:)    ! chauffage solaire
266      REAL,save,allocatable :: cool(:,:)    ! refroidissement infrarouge
267      REAL,save,allocatable :: dtrad(:,:)   ! K s-1
268      REAL,save,allocatable :: topsw(:), toplw(:)
269      REAL,save,allocatable :: solsw(:), sollw(:)
270      REAL,save,allocatable :: sollwdown(:) ! downward LW flux at surface
271      REAL tmpout(klon,klev)  ! K s-1
272
273      INTEGER itaprad
274      SAVE itaprad
275      REAL zdtime
276c
277
278c CHIMIE
279
280      REAL    dtimechim
281      INTEGER itapchim,appel_chim
282      SAVE itapchim,dtimechim
283
284c ORBITE
285
286      REAL dist, rmu0(klon), fract(klon), pdecli
287      REAL zday
288      REAL zls
289c
290      INTEGER i, k, iq, ig, j, ll, l
291c
292      REAL zphi(klon,klev)
293      REAL zzlev(klon,klev+1),zzlay(klon,klev),z1,z2
294c
295c Variables du changement
296c
297c ajs: ajustement sec
298c vdf: couche limite (Vertical DiFfusion)
299c mph: microphysique
300c kim: chimie
301      REAL d_t_ajs(klon,klev), d_tr_ajs(klon,klev,nqmax)
302      REAL d_u_ajs(klon,klev), d_v_ajs(klon,klev)
303c
304      REAL d_ts(klon)
305c
306      REAL d_u_vdf(klon,klev), d_v_vdf(klon,klev)
307      REAL d_t_vdf(klon,klev), d_tr_vdf(klon,klev,nqmax)
308c
309      REAL d_tr_mph(klon,klev,nqmax),d_tr_kim(klon,klev,nqmax)
310
311CMOD LOTT: Tendances Orography Sous-maille
312      REAL d_u_oro(klon,klev), d_v_oro(klon,klev)
313      REAL d_t_oro(klon,klev)
314      REAL d_u_lif(klon,klev), d_v_lif(klon,klev)
315      REAL d_t_lif(klon,klev)
316C          Tendances Ondes de G non oro (runs strato).
317      REAL d_u_hin(klon,klev), d_v_hin(klon,klev)
318      REAL d_t_hin(klon,klev)
319
320c
321c Variables liees a l'ecriture de la bande histoire physique
322c
323      INTEGER ecrit_mth
324      SAVE ecrit_mth   ! frequence d'ecriture (fichier mensuel)
325c
326      INTEGER ecrit_day
327      SAVE ecrit_day   ! frequence d'ecriture (fichier journalier)
328c
329      INTEGER ecrit_ins
330      SAVE ecrit_ins   ! frequence d'ecriture (fichier instantane)
331c
332      integer itau_w   ! pas de temps ecriture = itap + itau_phy
333
334c Variables locales pour effectuer les appels en serie
335c
336      REAL t_seri(klon,klev)
337      REAL u_seri(klon,klev), v_seri(klon,klev)
338c
339      REAL tr_seri(klon,klev,nqmax)
340
341      INTEGER        length
342      PARAMETER    ( length = 100 )
343      REAL tabcntr0( length       )
344c
345c pour ioipsl
346      INTEGER ndex2d(iim*jjmp1),ndex3d(iim*jjmp1*klev)
347      REAL zx_tmp_fi2d(klon)      ! variable temporaire grille physique
348      REAL zx_tmp_fi3d(klon,klev) ! variable temporaire pour champs 3D
349      REAL zx_tmp_2d(iim,jjmp1),zx_tmp_3d(iim,jjmp1,klev)
350      REAL zx_lon(iim,jjmp1), zx_lat(iim,jjmp1)
351
352      INTEGER nid_day, nid_mth, nid_ins
353      SAVE nid_day, nid_mth, nid_ins
354c
355      INTEGER nhori, nvert, idayref
356      REAL zsto, zout, zsto1, zsto2, zero
357      parameter (zero=0.0e0)
358      real zjulian
359      save zjulian
360
361      CHARACTER*2  str2
362      character*20 modname
363      character*80 abort_message
364      logical ok_sync
365
366      character*30 nom_fichier
367      character*10 varname
368      character*40 vartitle
369      character*20 varunits
370C     Variables liees au bilan d'energie et d'enthalpi
371      REAL ztsol(klon)
372      REAL      h_vcol_tot, h_dair_tot, h_qw_tot, h_ql_tot
373     $        , h_qs_tot, qw_tot, ql_tot, qs_tot , ec_tot
374      SAVE      h_vcol_tot, h_dair_tot, h_qw_tot, h_ql_tot
375     $        , h_qs_tot, qw_tot, ql_tot, qs_tot , ec_tot
376      REAL      d_h_vcol, d_h_dair, d_qt, d_qw, d_ql, d_qs, d_ec
377      REAL      d_h_vcol_phy
378      REAL      fs_bound, fq_bound
379      SAVE      d_h_vcol_phy
380      REAL      zero_v(klon),zero_v2(klon,klev)
381      CHARACTER*15 ztit
382      INTEGER   ip_ebil  ! PRINT level for energy conserv. diag.
383      SAVE      ip_ebil
384      DATA      ip_ebil/2/
385      INTEGER   if_ebil ! level for energy conserv. dignostics
386      SAVE      if_ebil
387c+jld ec_conser
388      REAL d_t_ec(klon,klev)    ! tendance du a la conversion Ec -> E thermique
389c-jld ec_conser
390
391c TEST VENUS...
392      REAL mang(klon,klev)    ! moment cinetique
393      REAL mangtot            ! moment cinetique total
394
395      CHARACTER*2 str1
396
397c Temporaire avant de trouver mieux :
398c Recuperation des TAU du TR
399      REAL t_tauhvd(klon,klev),t_khvd(klon,klev)
400      REAL t_tcld(klon,klev),t_kcld(klon,klev)
401      REAL t_kcvd(klon,klev)
402c  ASTUCE POUR EVITER klon... EN ATTENDANT MIEUX
403      INTEGER   ngrid
404      PARAMETER (ngrid=(jjm-1)*iim+2)  ! = klon
405      INTEGER NSPECV,NSPECI,NLAYER
406      PARAMETER (NSPECV=24,NSPECI=46,NLAYER=llm)
407      REAL TAUHID(ngrid,NLAYER,NSPECI)
408     &               ,TAUCID(ngrid,NLAYER,NSPECI)
409     &               ,TAUGID(ngrid,NLAYER,NSPECI)
410     &               ,TAUHVD(ngrid,NLAYER,NSPECV)
411     &               ,TAUCVD(ngrid,NLAYER,NSPECV)
412     &               ,TAUGVD(ngrid,NLAYER,NSPECV)
413
414      COMMON /TAUD/   TAUHID,TAUCID,TAUGID,TAUHVD,TAUCVD,TAUGVD
415* common relatifs au nuages
416      real rmcbar(ngrid,NLAYER),xfbar(ngrid,NLAYER,4)
417      integer ncount(ngrid,NLAYER)
418      COMMON/rnuabar/ncount,rmcbar,xfbar
419
420       REAL ch4(klon,jjm+1),dch4(jjm+1)
421       INTEGER ig0
422       integer ich4
423       common/ch4ind/ich4
424
425c     flux de chaleur latente d'evaporation CH4
426      REAL fclat(klon)
427c     reservoir de surface
428      REAL,save,allocatable :: reservoir(:)
429
430c Declaration des constantes et des fonctions thermodynamiques
431c
432#include "YOMCST.h"
433
434c======================================================================
435c INITIALISATIONS
436c======================================================================
437
438      modname = 'physiq'
439      ok_sync=.TRUE.
440
441      bilansmc = 0
442      ballons  = 0
443
444      IF (if_ebil.ge.1) THEN
445        DO i=1,klon
446          zero_v(i)=0.
447        END DO
448        DO i=1,klon
449         DO j=1,klev
450          zero_v2(i,j)=0.
451         END DO
452        END DO
453      END IF
454     
455c PREMIER APPEL SEULEMENT
456c========================
457      IF (debut) THEN
458         allocate(t_ancien(klon,klev),u_ancien(klon,klev))
459         allocate(swnet(klon,klevp1),lwnet(klon,klevp1))
460         allocate(radsol(klon),ftsol(klon),falbe(klon))
461         allocate(rlev(klon,klevp1),ftsoil(klon,nsoilmx))
462         allocate(zmea(klon),zstd(klon),zsig(klon),zgam(klon))
463         allocate(zthe(klon),zpic(klon),zval(klon),rugoro(klon))
464         allocate(zuthe(klon),zvthe(klon),dlw(klon),fder(klon))
465         allocate(heat(klon,klev),cool(klon,klev))
466         allocate(dtrad(klon,klev),topsw(klon),toplw(klon))
467         allocate(solsw(klon),sollw(klon),sollwdown(klon))
468         allocate(source(klon,nqmax))
469         allocate(reservoir(klon))
470
471         CALL suphec ! initialiser constantes et parametres phys.
472
473         IF (if_ebil.ge.1) d_h_vcol_phy=0.
474c
475c appel a la lecture du physiq.def
476c
477         call conf_phys(ok_mensuel,ok_journe,ok_instan,if_ebil)
478
479c
480c Initialiser les compteurs:
481c
482         itap        = 0
483         itaprad     = 0
484         itapchim    = 0
485         ncount(:,:) = 0
486
487c         
488c Lecture startphy.nc :
489c
490c REMETTRE TOUS LES PARAMETRES POUR OROGW...  A FAIRE POUR TITAN
491         CALL phyetat0 ("startphy.nc",dtime,
492     .       rlatd,rlond,ftsol,ftsoil,
493     .       falbe, solsw, sollw,
494     .       dlw,radsol,reservoir,
495c     .      zmea,zstd,zsig,zgam,zthe,zpic,zval,
496     .       tabcntr0,
497     .       t_ancien, ancien_ok)
498
499c dtime est lu dans startphy
500c pdtphys est calcule a partir des nouvelles conditions:
501c Reinitialisation du pas de temps physique quand changement
502         IF (ABS(dtime-pdtphys).GT.0.001) THEN
503            WRITE(lunout,*) 'Pas physique a change',dtime,
504     .                        pdtphys
505c           abort_message='Pas physique n est pas correct '
506c           call abort_gcm(modname,abort_message,1)
507c----------------
508c pour initialiser convenablement le time_counter, il faut tenir compte
509c du changement de dtime en changeant itau_phy (point de depart)
510            itau_phy = NINT(itau_phy*dtime/pdtphys)
511c----------------
512            dtime=pdtphys
513         ENDIF
514
515         radpas  = NINT( RDAY/pdtphys/nbapp_rad)
516         chimpas =   radpas*nbapp_rad/nbapp_chim
517
518         CALL printflag( tabcntr0,radpas,chimpas,ok_mensuel,
519     .                             ok_journe, ok_instan )
520
521c
522c Initialiser les pas de temps:
523c
524      dtimerad = dtime*FLOAT(radpas)  ! pas de temps du rayonnement (s)
525c      PRINT*,'dtimerad,dtime,radpas',dtimerad,dtime,radpas
526           
527      dtimechim = dtime*FLOAT(chimpas)  ! pas de temps de la chimie (s)
528c      PRINT*,'dtimechim,dtime,chimpas',dtimechim,dtime,chimpas
529
530
531c INITIALISATION ORBITE
532
533         CALL iniorbit(aphelie,periheli,year_day,peri_day,obliquit)
534
535c---------
536c FLOTT
537c       IF (ok_orodr) THEN
538c         DO i=1,klon
539c         rugoro(i) = MAX(1.0e-05, zstd(i)*zsig(i)/2.0)
540c         ENDDO
541c         CALL SUGWD(klon,klev,paprs,pplay)
542c         DO i=1,klon
543c         zuthe(i)=0.
544c         zvthe(i)=0.
545c         if(zstd(i).gt.10.)then
546c           zuthe(i)=(1.-zgam(i))*cos(zthe(i))
547c           zvthe(i)=(1.-zgam(i))*sin(zthe(i))
548c         endif
549c         ENDDO
550c       ENDIF
551
552      if (bilansmc.eq.1) then
553C  OUVERTURE D'UN FICHIER FORMATTE POUR STOCKER LES COMPOSANTES
554C  DU BILAN DE MOMENT ANGULAIRE.
555      open(27,file='aaam_bud.out',form='formatted')
556      open(28,file='fields_2d.out',form='formatted')
557      write(*,*)'Ouverture de aaam_bud.out (FL Vous parle)'
558      write(*,*)'Ouverture de fields_2d.out (FL Vous parle)'
559      endif !bilansmc
560
561c--------------SLEBONNOIS
562C  OUVERTURE DES FICHIERS FORMATTES CONTENANT LES POSITIONS ET VITESSES
563C  DES BALLONS
564      if (ballons.eq.1) then
565      open(30,file='ballons-lat.out',form='formatted')
566      open(31,file='ballons-lon.out',form='formatted')
567      open(32,file='ballons-u.out',form='formatted')
568      open(33,file='ballons-v.out',form='formatted')
569      open(34,file='ballons-alt.out',form='formatted')
570      write(*,*)'Ouverture des ballons*.out'
571      endif !ballons
572c-------------
573
574c---------
575C TRACEURS
576C source dans couche limite
577          source = 0.0 ! pas de source, pour l'instant
578C
579c Si microphysique offline, pas besoin d'avoir de traceurs microphysiques
580c car on lit les profils verticaux des qaer dans une look-up table pour
581c le rayonnement.
582
583c  calcul de nmicro
584c !!!! Les traceurs microphysiques doivent etre toujours en premiers!!
585
586      nmicro = 0
587      do iq=1,nqmax
588         nom = tname(iq)
589c        print*,iq,"nom=",nom,"tname=",tname(iq)
590         print*,iq,"nom=",nom
591         if (nom(1:1).eq."q") then
592           nmicro = nmicro+1
593         endif
594      enddo
595      print*,"nmicro=",nmicro
596
597c
598c Verifications:
599c
600         IF ((nmicro.eq.0).and.(microfi.eq.1)) THEN
601           abort_message="MICROPHYSIQUE ONLINE, MAIS NMICRO=0..."
602           call abort_gcm(modname,abort_message,1)
603         ENDIF
604         IF (microfi.lt.1.and.clouds.eq.1) THEN
605          write(lunout,*)"microfi.lt.1.and.clouds.eq.1"
606          abort_message =
607     &    "Impossible de faire des nuages sans microphysique..."
608          call abort_gcm(modname,abort_message,1)
609         ENDIF
610         IF (nlon .NE. klon) THEN
611            WRITE(lunout,*)'nlon et klon ne sont pas coherents', nlon,
612     .                      klon
613            abort_message='nlon et klon ne sont pas coherents'
614            call abort_gcm(modname,abort_message,1)
615         ENDIF
616         IF (nlev .NE. klev) THEN
617            WRITE(lunout,*)'nlev et klev ne sont pas coherents', nlev,
618     .                       klev
619            abort_message='nlev et klev ne sont pas coherents'
620            call abort_gcm(modname,abort_message,1)
621         ENDIF
622c
623         IF (dtime*FLOAT(radpas).GT.(RDAY*0.25).AND.cycle_diurne)
624     $    THEN
625           WRITE(lunout,*)'Nbre d appels au rayonnement insuffisant'
626           WRITE(lunout,*)"Au minimum 4 appels par jour si cycle diurne"
627           abort_message='Nbre d appels au rayonnement insuffisant'
628           call abort_gcm(modname,abort_message,1)
629         ENDIF
630c
631         WRITE(lunout,*)"Clef pour la convection seche, iflag_ajs=",
632     .                   iflag_ajs
633c
634         ecrit_mth = NINT(RDAY/dtime) *nday  ! tous les nday jours
635         IF (ok_mensuel) THEN
636         WRITE(lunout,*)'La frequence de sortie mensuelle est de ',
637     .                   ecrit_mth
638         ENDIF
639
640         ecrit_day = NINT(RDAY/dtime *1.0)  ! tous les jours
641         IF (ok_journe) THEN
642         WRITE(lunout,*)'La frequence de sortie journaliere est de ',
643     .                   ecrit_day
644         ENDIF
645
646         ecrit_ins = NINT(RDAY/dtime*ecritphy)  ! Fraction de jour reglable
647         IF (ok_instan) THEN
648         WRITE(lunout,*)'La frequence de sortie instant. est de ',
649     .                   ecrit_ins
650         ENDIF
651
652c Initialisation des sorties
653c===========================
654
655#ifdef CPP_IOIPSL
656
657#ifdef histmth
658#include "ini_histmth.h"
659#endif
660
661#ifdef histday
662#include "ini_histday.h"
663#endif
664
665#ifdef histins
666#include "ini_histins.h"
667#endif
668
669#endif
670
671c
672c Initialiser les valeurs de u pour calculs tendances
673c (pour T, c'est fait dans phyetat0)
674c
675      DO k = 1, klev
676      DO i = 1, klon
677         u_ancien(i,k) = u(i,k)
678      ENDDO
679      ENDDO
680
681      rmcbar  = 0.
682      xfbar   = 0.
683         
684      ENDIF ! debut
685c====================================================================
686c======================================================================
687
688c   Creer un reservoir de surface infini
689c
690      reservoir(:) = 2.
691
692c Mettre a zero des variables de sortie (pour securite)
693c
694      DO i = 1, klon
695         d_ps(i) = 0.0
696      ENDDO
697      DO k = 1, klev
698      DO i = 1, klon
699         d_t(i,k) = 0.0
700         d_u(i,k) = 0.0
701         d_v(i,k) = 0.0
702      ENDDO
703      ENDDO
704      DO iq = 1, nqmax
705      DO k = 1, klev
706      DO i = 1, klon
707         d_qx(i,k,iq) = 0.0
708      ENDDO
709      ENDDO
710      ENDDO
711c
712c Ne pas affecter les valeurs entrees de u, v, h, et q
713c
714      DO k = 1, klev
715      DO i = 1, klon
716         t_seri(i,k)  = t(i,k)
717         u_seri(i,k)  = u(i,k)
718         v_seri(i,k)  = v(i,k)
719      ENDDO
720      ENDDO
721      DO iq = 1, nqmax
722      DO  k = 1, klev
723      DO  i = 1, klon
724         tr_seri(i,k,iq) = qx(i,k,iq)
725      ENDDO
726      ENDDO
727      ENDDO
728C
729      DO i = 1, klon
730          ztsol(i) = ftsol(i)
731      ENDDO
732C
733      IF (if_ebil.ge.1) THEN
734        ztit='after dynamic'
735        CALL diagetpq(airephy,ztit,ip_ebil,1,1,dtime
736     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
737     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
738C     Comme les tendances de la physique sont ajoute dans la dynamique,
739C     on devrait avoir que la variation d'entalpie par la dynamique
740C     est egale a la variation de la physique au pas de temps precedent.
741C     Donc la somme de ces 2 variations devrait etre nulle.
742        call diagphy(airephy,ztit,ip_ebil
743     e      , zero_v, zero_v, zero_v, zero_v, zero_v
744     e      , zero_v, zero_v, zero_v, ztsol
745     e      , d_h_vcol+d_h_vcol_phy, d_qt, 0.
746     s      , fs_bound, fq_bound )
747      END IF
748
749c====================================================================
750c Diagnostiquer la tendance dynamique
751c
752      IF (ancien_ok) THEN
753         DO k = 1, klev
754         DO i = 1, klon
755            d_u_dyn(i,k) = (u_seri(i,k)-u_ancien(i,k))/dtime
756            d_t_dyn(i,k) = (t_seri(i,k)-t_ancien(i,k))/dtime
757         ENDDO
758         ENDDO
759
760! ADAPTATION GCM POUR CP(T)
761         do i=1,klon
762          flux_dyn(i,1) = 0.0
763          do j=2,klev
764            flux_dyn(i,j) = flux_dyn(i,j-1)
765     . +cpdet(t_seri(i,j-1))/RG*d_t_dyn(i,j-1)*(paprs(i,j-1)-paprs(i,j))
766          enddo
767         enddo
768         
769      ELSE
770         DO k = 1, klev
771         DO i = 1, klon
772            d_u_dyn(i,k) = 0.0
773            d_t_dyn(i,k) = 0.0
774         ENDDO
775         ENDDO
776         ancien_ok = .TRUE.
777      ENDIF
778c====================================================================
779c
780c Ajouter le geopotentiel du sol:
781c
782      DO k = 1, klev
783      DO i = 1, klon
784         zphi(i,k) = pphi(i,k) + pphis(i)
785      ENDDO
786      ENDDO
787
788c   calcul du geopotentiel aux niveaux intercouches
789c   ponderation des altitudes au niveau des couches en dp/p
790
791      DO l=1,klev
792         DO i=1,klon
793c           zzlay(i,l)=zphi(i,l)/RG
794c SI ON TIENT COMPTE DE LA VARIATION DE G AVEC L'ALTITUDE:
795            zzlay(i,l)=RG*RA*RA/(RG*RA-zphi(i,l))-RA
796         ENDDO
797      ENDDO
798      DO i=1,klon
799c        zzlev(i,1)=0.
800c CORRECTION 13/01/2011 
801c (correspond a la position de la surface en ce point vs RA)
802         zzlev(i,1)=pphis(i)/RG
803      ENDDO
804      DO l=2,klev
805         DO i=1,klon
806            z1=(pplay(i,l-1)+paprs(i,l))/(pplay(i,l-1)-paprs(i,l))
807            z2=(paprs(i,l)+pplay(i,l))/(paprs(i,l)-pplay(i,l))
808            zzlev(i,l)=(z1*zzlay(i,l-1)+z2*zzlay(i,l))/(z1+z2)
809         ENDDO
810      ENDDO
811      DO i=1,klon
812         zzlev(i,klev+1)=zzlay(i,klev)+(zzlay(i,klev)-zzlev(i,klev))
813      ENDDO
814
815c- - - - - - - - - - - - - - - -
816c DIAGNOSTIQUE GRILLE VERTICALE
817c- - - - - - - - - - - - - - - -
818c     print*,"DIAGNOSTIQUE GRILLE VERTICALE"
819c     i=klon/2
820c     print*,"Niveau  Pression  Altitude    (lev puis lay)"
821c     do l=1,klev
822c      print*,l,paprs(i,l),zzlev(i,l)
823c      print*,l,pplay(i,l),zzlay(i,l)
824c     enddo
825c     print*,klev+1,paprs(i,klev+1),zzlev(i,klev+1)
826c     stop
827
828c====================================================================
829c
830c Verifier les temperatures
831c
832      CALL hgardfou(t_seri,ftsol,'debutphy')
833c====================================================================
834c
835c Incrementer le compteur de la physique
836c
837      itap   = itap + 1
838
839c====================================================================
840c
841c Epaisseurs couches
842
843      DO k = 1, klev
844      DO i = 1, klon
845         delp(i,k) = paprs(i,k)-paprs(i,k+1)
846      ENDDO
847      ENDDO
848
849
850
851
852c====================================================================
853c ORBITE ET ECLAIREMENT
854c====================================================================
855
856
857c Pour TITAN:
858c  calcul de la longitude solaire
859          CALL solarlong(rjourvrai+gmtime,zls)
860          print*,'Ls',zls*180./RPI      ! zls est en radians !!
861      CALL orbite(zls,dist,pdecli)
862
863c dans zenang, Ls en degres ; dans mucorr, Ls en radians
864      IF (cycle_diurne) THEN
865        zdtime=dtime*FLOAT(radpas) ! pas de temps du rayonnement (s)
866        CALL zenang(zls*180./RPI,gmtime,zdtime,rlatd,rlond,rmu0,fract)
867      ELSE
868        call mucorr(klon,zls,rlatd,rmu0,fract)
869      ENDIF
870
871c====================================================================
872c COUCHE LIMITE
873c====================================================================
874
875c-------------------------------
876c TEST: on ne tient pas compte des calculs de clmain mais on force
877c l'equilibre radiatif du sol
878      if (1.eq.0) then
879              if (debut) then
880                print*,"ATTENTION, CLMAIN SHUNTEE..."
881              endif
882
883      DO i = 1, klon
884         sens(i) = 0.0e0 ! flux de chaleur sensible au sol
885         fder(i) = 0.0e0
886         dlw(i)  = 0.0e0
887      ENDDO
888
889c Incrementer la temperature du sol
890c
891      DO i = 1, klon
892         d_ts(i)  = dtime * radsol(i)/22000. !valeur calculee par GCM pour I=200
893         ftsol(i) = ftsol(i) + d_ts(i)
894         do j=1,nsoilmx
895           ftsoil(i,j)=ftsol(i)
896         enddo
897      ENDDO
898
899c-------------------------------
900      else
901c-------------------------------
902
903      fder = dlw
904
905c     print*,"radsol avant clmain=",radsol(klon/2)
906c     print*,"solsw avant clmain=",solsw(klon/2)
907c     print*,"sollw avant clmain=",sollw(klon/2)
908
909c  CLMAIN
910
911! ADAPTATION GCM POUR CP(T)
912      CALL clmain(dtime,itap,
913     e            t_seri,u_seri,v_seri,
914     e            rmu0,
915     e            ftsol,
916     $            ftsoil,
917     $            paprs,pplay,ppk,radsol,falbe,
918     e            solsw, sollw, sollwdown, fder,
919     e            rlond, rlatd, cuphy, cvphy,
920     e            debut, lafin,
921     s            d_t_vdf,d_u_vdf,d_v_vdf,d_ts,
922     s            fluxt,fluxu,fluxv,cdragh,cdragm,
923     s            dsens,
924     s            ycoefh,yu1,yv1)
925
926c     print*,"radsol apres clmain=",radsol(klon/2)
927c     print*,"solsw apres clmain=",solsw(klon/2)
928c     print*,"sollw apres clmain=",sollw(klon/2)
929
930CXXX Incrementation des flux
931      DO i = 1, klon
932         sens(i) = - fluxt(i,1) ! flux de chaleur sensible au sol
933         fder(i) = dlw(i) + dsens(i)
934      ENDDO
935CXXX
936
937      DO k = 1, klev
938      DO i = 1, klon
939         t_seri(i,k) = t_seri(i,k) + d_t_vdf(i,k)
940         d_t_vdf(i,k)= d_t_vdf(i,k)/dtime          ! K/s
941         u_seri(i,k) = u_seri(i,k) + d_u_vdf(i,k)
942         d_u_vdf(i,k)= d_u_vdf(i,k)/dtime          ! (m/s)/s
943         v_seri(i,k) = v_seri(i,k) + d_v_vdf(i,k)
944         d_v_vdf(i,k)= d_v_vdf(i,k)/dtime          ! (m/s)/s
945      ENDDO
946      ENDDO
947
948c        print*,"d_t_vdf1=",d_t_vdf(1,:)*dtime
949c        print*,"d_t_vdf2=",d_t_vdf(klon/2,:)*dtime
950c        print*,"d_t_vdf3=",d_t_vdf(klon,:)*dtime
951c        print*,"d_u_vdf=",d_u_vdf(klon/2,:)*dtime
952c        print*,"d_v_vdf=",d_v_vdf(klon/2,:)*dtime
953
954C TRACEURS
955
956      d_tr_vdf = 0.
957      if (iflag_trac.eq.1) then
958      DO iq=1, nqmax
959          CALL cltrac(dtime,ycoefh,t_seri,
960     s               tr_seri(1,1,iq), source,
961     e               paprs, pplay, delp,
962     s               d_tr_vdf(1,1,iq))
963
964          tr_seri(:,:,iq) = tr_seri(:,:,iq) + d_tr_vdf(:,:,iq)
965          d_tr_vdf(:,:,iq)= d_tr_vdf(:,:,iq)/dtime          ! /s
966      ENDDO
967      endif
968
969      IF (if_ebil.ge.2) THEN
970        ztit='after clmain'
971        CALL diagetpq(airephy,ztit,ip_ebil,2,1,dtime
972     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
973     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
974         call diagphy(airephy,ztit,ip_ebil
975     e      , zero_v, zero_v, zero_v, zero_v, sens
976     e      , zero_v, zero_v, zero_v, ztsol
977     e      , d_h_vcol, d_qt, d_ec
978     s      , fs_bound, fq_bound )
979      END IF
980C
981c
982c Incrementer la temperature du sol
983c
984c     print*,'Tsol avant clmain:',ftsol(klon/2)
985      DO i = 1, klon
986         ftsol(i) = ftsol(i) + d_ts(i)
987      ENDDO
988c     print*,'DTsol apres clmain:',d_ts(klon/2)
989c     print*,'Tsol apres clmain:',ftsol(klon/2)
990
991c Calculer la derive du flux infrarouge
992c
993      DO i = 1, klon
994            dlw(i) = - 4.0*RSIGMA*ftsol(i)**3
995      ENDDO
996
997c-------------------------------
998      endif  ! fin du TEST
999
1000c
1001c Appeler l'ajustement sec
1002c
1003c===================================================================
1004c Convection seche
1005c===================================================================
1006c
1007      d_t_ajs(:,:)=0.
1008      d_u_ajs(:,:)=0.
1009      d_v_ajs(:,:)=0.
1010      d_tr_ajs(:,:,:)=0.
1011c
1012      IF(prt_level>9)WRITE(lunout,*)
1013     .    'AVANT LA CONVECTION SECHE , iflag_ajs='
1014     s   ,iflag_ajs
1015
1016      if(iflag_ajs.eq.0) then
1017c  Rien
1018c  ====
1019         IF(prt_level>9)WRITE(lunout,*)'pas de convection'
1020
1021      else if(iflag_ajs.eq.1) then
1022
1023c  Ajustement sec
1024c  ==============
1025         IF(prt_level>9)WRITE(lunout,*)'ajsec'
1026
1027! ADAPTATION GCM POUR CP(T)
1028         CALL ajsec(paprs, pplay, ppk, t_seri, u_seri, v_seri, nqmax,
1029     .              tr_seri, d_t_ajs, d_u_ajs, d_v_ajs, d_tr_ajs)
1030
1031! ADAPTATION GCM POUR CP(T)
1032         do i=1,klon
1033          flux_ajs(i,1) = 0.0
1034          do j=2,klev
1035            flux_ajs(i,j) = flux_ajs(i,j-1)
1036     .        + cpdet(t_seri(i,j-1))/RG*d_t_ajs(i,j-1)/dtime
1037     .                                 *delp(i,j-1)
1038          enddo
1039         enddo
1040         
1041         t_seri(:,:) = t_seri(:,:) + d_t_ajs(:,:)
1042         d_t_ajs(:,:)= d_t_ajs(:,:)/dtime          ! K/s
1043         u_seri(:,:) = u_seri(:,:) + d_u_ajs(:,:)
1044         d_u_ajs(:,:)= d_u_ajs(:,:)/dtime          ! (m/s)/s
1045         v_seri(:,:) = v_seri(:,:) + d_v_ajs(:,:)
1046         d_v_ajs(:,:)= d_v_ajs(:,:)/dtime          ! (m/s)/s
1047      if (iflag_trac.eq.1) then
1048         tr_seri(:,:,:) = tr_seri(:,:,:) + d_tr_ajs(:,:,:)
1049         d_tr_ajs(:,:,:)= d_tr_ajs(:,:,:)/dtime          ! /s
1050      endif
1051     
1052c        print*,"d_t_ajs1=",d_t_ajs(1,:)*dtime
1053c        print*,"d_t_ajs2=",d_t_ajs(klon/2,:)*dtime
1054c        print*,"d_t_ajs3=",d_t_ajs(klon,:)*dtime
1055c        print*,"d_u_ajs=",d_u_ajs(klon/2,:)*dtime
1056c        print*,"d_v_ajs=",d_v_ajs(klon/2,:)*dtime
1057
1058      endif
1059c
1060      IF (if_ebil.ge.2) THEN
1061        ztit='after dry_adjust'
1062        CALL diagetpq(airephy,ztit,ip_ebil,2,2,dtime
1063     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
1064     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1065        call diagphy(airephy,ztit,ip_ebil
1066     e      , zero_v, zero_v, zero_v, zero_v, sens
1067     e      , zero_v, zero_v, zero_v, ztsol
1068     e      , d_h_vcol, d_qt, d_ec
1069     s      , fs_bound, fq_bound )
1070      END IF
1071
1072c====================================================================
1073c   MICROPHYSIQUE ET CHIMIE
1074c====================================================================
1075
1076      d_tr_mph(:,:,:)=0.
1077      d_tr_kim(:,:,:)=0.
1078     
1079c on recupere tr_seri inchange, d_tr_micro, d_tr_chim, tous les trois sur nqmax
1080c on recupere aussi qaer pour le mettre dans les sorties
1081c  si microfi=1, sortie de qaer(1:nmicro)
1082c  si nmicro != nqmax et si chimi, sortie de tr_seri(nmicro+1:nqmax)
1083
1084c faire un test comme pour rayonnement, avec chimi en plus comme flag,
1085c pour voir si chimie appelee -> bouleen, qui passe dans phytrac.
1086c faut aussi le pas de temps chimique: dtimechim, a passer..
1087
1088      appel_chim = 0
1089      IF (MOD(itapchim,chimpas).EQ.0) THEN
1090c             print*,'CHIMIE ',
1091c    $             ' (itapchim=',itapchim,'/chimpas=',chimpas,')'
1092       appel_chim = 1
1093       itapchim = 0
1094      ENDIF
1095      itapchim = itapchim + 1
1096
1097      if (iflag_trac.eq.1) then
1098c         call begintime(tt0)
1099         call phytrac (debut,lafin,
1100     .                   nqmax,nmicro,dtime,appel_chim,dtimechim,
1101     .                   paprs,pplay,delp,t,rmu0,fract,pdecli,zls,
1102     .                   yu1,yv1,zzlev,zzlay,ftsol,
1103     .                   tr_seri,qaer,d_tr_mph,d_tr_kim,
1104     .                   fclat,reservoir)
1105
1106c         call endtime(tt0,tt1)
1107c         ttphytra=ttphytra+tt1
1108
1109c ----- ICI on ajuste radsol en tenant compte du flux de chaleur latente
1110c       d'evaporation du reservoir.
1111c       NOTE : c'est pas tres elegant mais ca permet d'eviter d'aller
1112c              toucher a clmain.
1113        if (clouds.eq.1) then
1114          radsol(:) = radsol(:)+fclat(:)    !test pas de flx de chaleur latente
1115        endif
1116
1117        if (microfi.ge.1) then
1118         tr_seri(:,:,1:nmicro) = tr_seri(:,:,1:nmicro)
1119     .                        + d_tr_mph(:,:,1:nmicro)*dtime
1120        endif
1121c       PAS ELEGANT mais je n'ai pas trouve d'autres solutions :
1122c       Il semblerait qu'il y ait un probleme lorsque les tendances de traceurs
1123c       retourne des traceurs nuls et il y a parfois des valeurs negatives qui trainent.
1124c       Pour ne diffuser le probleme, on force les valeurs negatives a ZERO.
1125        DO iq=1,nmicro
1126          DO i=1,klon
1127            DO l=1,klev
1128              if (tr_seri(i,l,iq).lt.0.) then
1129                 tr_seri(i,l,iq) = 0.
1130              endif
1131            ENDDO
1132          ENDDO
1133        ENDDO
1134
1135c condensation:
1136c       NE PAS OUBLIER LA CONDENSATION DES NUAGES !!!!
1137        if ((clouds.eq.1.or.(chimi)).and.nqmax.gt.nmicro) then
1138          tr_seri(:,:,nmicro+1:nqmax) = tr_seri(:,:,nmicro+1:nqmax)
1139     .                         + d_tr_mph(:,:,nmicro+1:nqmax)*dtime
1140        endif
1141        if ((chimi).and.(nqmax.gt.nmicro)) then
1142c chimie:
1143         tr_seri(:,:,:) = tr_seri(:,:,:) + d_tr_kim(:,:,:)*dtime
1144        endif
1145
1146      endif
1147
1148c       ch4=0.
1149c       do l=1,llm
1150c         ch4(1,l) = tr_seri(1,l,ich4)
1151c         do j=2,jjm
1152c           ig0=1+(j-2)*iim
1153c           do i=1,iim
1154c             ch4(j,l)= ch4(j,l)  + tr_seri(ig0+i,l,ich4)/iim
1155c           enddo
1156c         enddo
1157c         ch4(jjm+1,l) = tr_seri(klon,l,ich4)
1158c       enddo
1159c       do j=1,jjm+1
1160c         write(501,*) j,ch4(j,1)
1161c       enddo
1162c       do l=1,llm
1163c         write(502,'(I3,49(ES24.17,1X))') l,
1164c     &   (ch4(j,l),j=1,jjm+1)
1165c       enddo
1166c       write(501,*) ""
1167c       write(502,*) ""
1168
1169c------------------
1170c test condensation
1171c     do i=1,nqmax
1172c       if(tname(i).eq."HCN") then
1173c          print*,"HCN="
1174c          do k=1,klev
1175c           print*,k,tr_seri(klon/2,k,i),d_tr_mph(klon/2,k,i)*dtime
1176c    v      ,d_tr_kim(klon/2,k,i)*dtime
1177c          enddo
1178c          stop
1179c       endif
1180c     enddo
1181c------------------
1182
1183c====================================================================
1184c RAYONNEMENT
1185c====================================================================
1186
1187      IF (MOD(itaprad,radpas).EQ.0) THEN
1188c             print*,'RAYONNEMENT ',
1189c    $             ' (itaprad=',itaprad,'/radpas=',radpas,')'
1190
1191c ATTENTION, (klon/2) ne marche pas toujours............
1192c     print*,"radsol avant radlwsw=",radsol(klon/2)
1193c     print*,"solsw avant radlwsw=",solsw(klon/2)
1194c     print*,"sollw avant radlwsw=",sollw(klon/2)
1195c     print*,"avant radlwsw"
1196
1197c   ----------------
1198c   Calcul du rayon moyen des gouttes et des fractions volumique pour le TR
1199c  ----------------
1200      IF (clouds.eq.1) THEN
1201        DO i=1,klon
1202          DO j=1,klev
1203            rmcbar(i,j)=rmcbar(i,j)/MAX(FLOAT(ncount(i,j)),1.)
1204            xfbar(i,j,:)=xfbar(i,j,:)/MAX(FLOAT(ncount(i,j)),1.)
1205          ENDDO
1206        ENDDO
1207      ENDIF
1208     
1209c      call begintime(tt0)
1210      CALL radlwsw
1211     e            (dist, rmu0, fract, dtimerad, zzlev,
1212     e             paprs, pplay,ftsol, t_seri, nqmax, nmicro,
1213     c             tr_seri, qaer,
1214     s             heat,cool,radsol,
1215     s             topsw,toplw,solsw,sollw,
1216     s             sollwdown,
1217     s             lwnet, swnet)
1218c      call endtime(tt0,tt1)
1219c      ttrad=ttrad+tt1
1220
1221c     print*,"apres radlwsw"
1222c     mise a zero du rayon moyen des gouttes et des fractions volumique
1223      IF (clouds.eq.1) THEN
1224        rmcbar(:,:)  = 0.
1225        xfbar(:,:,:) = 0.
1226        ncount(:,:)  = 0
1227      ENDIF
1228
1229c     print*,"radsol apres radlwsw=",radsol(klon/2)
1230c     print*,"solsw apres radlwsw=",solsw(klon/2)
1231c     print*,"sollw apres radlwsw=",sollw(klon/2)
1232      itaprad = 0
1233      DO k = 1, klev
1234       DO i = 1, klon
1235         dtrad(i,k) = heat(i,k)-cool(i,k)     !K/s
1236       ENDDO
1237      ENDDO
1238c       print*,"heat (K/s) =",heat(klon/2,:)
1239c       print*,"cool (K/s) =",cool(klon/2,:)
1240c       print*,"dtrad1 (K/s) =",dtrad(1,:)
1241c       print*,"dtrad2 (K/s) =",dtrad(klon/2,:)
1242c       print*,"dtrad3 (K/s) =",dtrad(klon,:)
1243           
1244      ENDIF
1245      itaprad = itaprad + 1
1246c====================================================================
1247c
1248c Ajouter la tendance des rayonnements (tous les pas)
1249c
1250      DO k = 1, klev
1251      DO i = 1, klon
1252         t_seri(i,k) = t_seri(i,k) + dtrad(i,k) * dtime
1253      ENDDO
1254      ENDDO
1255 
1256      IF (if_ebil.ge.2) THEN
1257        ztit='after rad'
1258        CALL diagetpq(airephy,ztit,ip_ebil,2,2,dtime
1259     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
1260     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1261        call diagphy(airephy,ztit,ip_ebil
1262     e      , topsw, toplw, solsw, sollw, zero_v
1263     e      , zero_v, zero_v, zero_v, ztsol
1264     e      , d_h_vcol, d_qt, d_ec
1265     s      , fs_bound, fq_bound )
1266      END IF
1267c
1268
1269c====================================================================
1270c+jld ec_conser
1271      DO k = 1, klev
1272      DO i = 1, klon
1273        d_t_ec(i,k)=0.5/cpdet(t_seri(i,k))
1274     $      *(u(i,k)**2+v(i,k)**2-u_seri(i,k)**2-v_seri(i,k)**2)
1275        t_seri(i,k)=t_seri(i,k)+d_t_ec(i,k)
1276        d_t_ec(i,k) = d_t_ec(i,k)/dtime
1277       END DO
1278      END DO
1279         do i=1,klon
1280          flux_ec(i,1) = 0.0
1281          do j=2,klev
1282            flux_ec(i,j) = flux_ec(i,j-1)
1283     . +cpdet(t_seri(i,j-1))/RG*d_t_ec(i,j-1)*delp(i,j-1)
1284          enddo
1285         enddo
1286         
1287c-jld ec_conser
1288c====================================================================
1289
1290      IF (if_ebil.ge.1) THEN
1291        ztit='after physic'
1292        CALL diagetpq(airephy,ztit,ip_ebil,1,1,dtime
1293     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
1294     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1295C     Comme les tendances de la physique sont ajoute dans la dynamique,
1296C     on devrait avoir que la variation d'entalpie par la dynamique
1297C     est egale a la variation de la physique au pas de temps precedent.
1298C     Donc la somme de ces 2 variations devrait etre nulle.
1299        call diagphy(airephy,ztit,ip_ebil
1300     e      , topsw, toplw, solsw, sollw, sens
1301     e      , zero_v, zero_v, zero_v, ztsol
1302     e      , d_h_vcol, d_qt, d_ec
1303     s      , fs_bound, fq_bound )
1304C
1305      d_h_vcol_phy=d_h_vcol
1306C
1307      END IF
1308C
1309c====================================================================
1310c   Calcul  des gravity waves  FLOTT
1311c====================================================================
1312c
1313c      if (ok_orodr.or.ok_gw_nonoro) then
1314cc  CALCUL DE N2
1315c       do i=1,klon
1316c        do k=2,klev
1317c         ztlev(i,k)  = (t_seri(i,k)+t_seri(i,k-1))/2.
1318c         zpklev(i,k) = sqrt(ppk(i,k)*ppk(i,k-1))
1319c       enddo
1320c       enddo
1321c       call t2tpot(klon*klev,ztlev, ztetalev,zpklev)
1322c       call t2tpot(klon*klev,t_seri,ztetalay,ppk)
1323c       do i=1,klon
1324c        do k=2,klev
1325c         zdtetalev(i,k) = ztetalay(i,k)-ztetalay(i,k-1)
1326c         zdzlev(i,k)    = (zphi(i,k)-zphi(i,k-1))/RG
1327c          zn2(i,k) = RG*zdtetalev(i,k)/(ztetalev(i,k)*zdzlev(i,k))
1328c          zn2(i,k) = max(zn2(i,k),1.e-12)  ! securite
1329c       enddo
1330c       enddo
1331c
1332c      endif
1333c     
1334cc ----------------------------ORODRAG
1335c      IF (ok_orodr) THEN
1336cc
1337cc  selection des points pour lesquels le shema est actif:
1338c        igwd=0
1339c        DO i=1,klon
1340c        itest(i)=0
1341cc        IF ((zstd(i).gt.10.0)) THEN
1342c        IF (((zpic(i)-zmea(i)).GT.100.).AND.(zstd(i).GT.10.0)) THEN
1343c          itest(i)=1
1344c          igwd=igwd+1
1345c          idx(igwd)=i
1346c        ENDIF
1347c        ENDDO
1348cc        igwdim=MAX(1,igwd)
1349cc
1350cc A ADAPTER POUR TITAN !!!
1351c        CALL drag_noro(klon,klev,dtime,paprs,pplay,zphi,zn2,
1352c     e                   zmea,zstd, zsig, zgam, zthe,zpic,zval,
1353c     e                   igwd,idx,itest,
1354c     e                   t_seri, u_seri, v_seri,
1355c     s                   zulow, zvlow, zustrdr, zvstrdr,
1356c     s                   d_t_oro, d_u_oro, d_v_oro)
1357c
1358cc  ajout des tendances
1359c           t_seri(:,:) = t_seri(:,:) + d_t_oro(:,:)
1360c           d_t_oro(:,:)= d_t_oro(:,:)/dtime          ! K/s
1361c           u_seri(:,:) = u_seri(:,:) + d_u_oro(:,:)
1362c           d_u_oro(:,:)= d_u_oro(:,:)/dtime          ! (m/s)/s
1363c           v_seri(:,:) = v_seri(:,:) + d_v_oro(:,:)
1364c           d_v_oro(:,:)= d_v_oro(:,:)/dtime          ! (m/s)/s
1365cc   
1366c      ELSE
1367c         d_t_oro = 0.
1368c         d_u_oro = 0.
1369c         d_v_oro = 0.
1370c        zustrdr = 0.
1371c        zvstrdr = 0.
1372cc
1373c      ENDIF ! fin de test sur ok_orodr
1374cc
1375cc ----------------------------OROLIFT
1376c      IF (ok_orolf) THEN
1377cc
1378cc  selection des points pour lesquels le shema est actif:
1379c        igwd=0
1380c        DO i=1,klon
1381c        itest(i)=0
1382c        IF ((zpic(i)-zmea(i)).GT.100.) THEN
1383c          itest(i)=1
1384c          igwd=igwd+1
1385c          idx(igwd)=i
1386c        ENDIF
1387c        ENDDO
1388cc        igwdim=MAX(1,igwd)
1389cc
1390cc A ADAPTER POUR VENUS!!!
1391cc            CALL lift_noro(klon,klev,dtime,paprs,pplay,
1392cc     e                   rlatd,zmea,zstd,zpic,zgam,zthe,zpic,zval,
1393cc     e                   igwd,idx,itest,
1394cc     e                   t_seri, u_seri, v_seri,
1395cc     s                   zulow, zvlow, zustrli, zvstrli,
1396cc     s                   d_t_lif, d_u_lif, d_v_lif               )
1397c
1398cc
1399cc  ajout des tendances
1400c           t_seri(:,:) = t_seri(:,:) + d_t_lif(:,:)
1401c           d_t_lif(:,:)= d_t_lif(:,:)/dtime          ! K/s
1402c           u_seri(:,:) = u_seri(:,:) + d_u_lif(:,:)
1403c           d_u_lif(:,:)= d_u_lif(:,:)/dtime          ! (m/s)/s
1404c           v_seri(:,:) = v_seri(:,:) + d_v_lif(:,:)
1405c           d_v_lif(:,:)= d_v_lif(:,:)/dtime          ! (m/s)/s
1406cc
1407c      ELSE
1408c         d_t_lif = 0.
1409c         d_u_lif = 0.
1410c         d_v_lif = 0.
1411c         zustrli = 0.
1412c         zvstrli = 0.
1413cc
1414c      ENDIF ! fin de test sur ok_orolf
1415c
1416cc ---------------------------- NON-ORO GRAVITY WAVES
1417c       IF(ok_gw_nonoro) then
1418c
1419c      call flott_gwd_ran(klon,klev,dtime,pplay,zn2,
1420c     e               t_seri, u_seri, v_seri,
1421c     o               zustrhi,zvstrhi,
1422c     o               d_t_hin, d_u_hin, d_v_hin)
1423c
1424cc  ajout des tendances
1425c
1426c         t_seri(:,:) = t_seri(:,:) + d_t_hin(:,:)
1427c         d_t_hin(:,:)= d_t_hin(:,:)/dtime          ! K/s
1428c         u_seri(:,:) = u_seri(:,:) + d_u_hin(:,:)
1429c         d_u_hin(:,:)= d_u_hin(:,:)/dtime          ! (m/s)/s
1430c         v_seri(:,:) = v_seri(:,:) + d_v_hin(:,:)
1431c         d_v_hin(:,:)= d_v_hin(:,:)/dtime          ! (m/s)/s
1432c
1433c      ELSE
1434c         d_t_hin = 0.
1435c         d_u_hin = 0.
1436c         d_v_hin = 0.
1437c         zustrhi = 0.
1438c         zvstrhi = 0.
1439c
1440c      ENDIF ! fin de test sur ok_gw_nonoro
1441c
1442c====================================================================
1443c Transport de ballons
1444c====================================================================
1445      if (ballons.eq.1) then
1446         CALL ballon(30,pdtphys,rjourvrai,gmtime,rlatd,rlond,
1447c    C               t,pplay,u,v,pphi)   ! alt above surface (smoothed for GCM)
1448     C               t,pplay,u,v,zphi)   ! alt above planet average radius
1449      endif !ballons
1450
1451c====================================================================
1452c Bilan de mmt angulaire
1453c====================================================================
1454      if (bilansmc.eq.1) then
1455CMODDEB FLOTT
1456C  CALCULER LE BILAN DE MOMENT ANGULAIRE (DIAGNOSTIQUE)
1457C STRESS NECESSAIRES: COUCHE LIMITE ET TOUTE LA PHYSIQUE
1458
1459      DO i = 1, klon
1460        zustrph(i)=0.
1461        zvstrph(i)=0.
1462        zustrcl(i)=0.
1463        zvstrcl(i)=0.
1464      ENDDO
1465      DO k = 1, klev
1466      DO i = 1, klon
1467       zustrph(i)=zustrph(i)+(u_seri(i,k)-u(i,k))/dtime*
1468     c            (paprs(i,k)-paprs(i,k+1))/rg
1469       zvstrph(i)=zvstrph(i)+(v_seri(i,k)-v(i,k))/dtime*
1470     c            (paprs(i,k)-paprs(i,k+1))/rg
1471       zustrcl(i)=zustrcl(i)+d_u_vdf(i,k)*
1472     c            (paprs(i,k)-paprs(i,k+1))/rg
1473       zvstrcl(i)=zvstrcl(i)+d_v_vdf(i,k)*
1474     c            (paprs(i,k)-paprs(i,k+1))/rg
1475      ENDDO
1476      ENDDO
1477
1478      CALL aaam_bud (27,klon,klev,rjourvrai,gmtime,
1479     C               ra,rg,romega,
1480     C               rlatd,rlond,pphis,
1481     C               zustrdr,zustrli,zustrcl,
1482     C               zvstrdr,zvstrli,zvstrcl,
1483     C               paprs,u,v)
1484                     
1485CCMODFIN FLOTT
1486      endif !bilansmc
1487
1488c=======================================================================
1489c   SORTIES
1490c=======================================================================
1491
1492c Convertir les incrementations en tendances
1493c
1494      DO k = 1, klev
1495      DO i = 1, klon
1496         d_u(i,k) = ( u_seri(i,k) - u(i,k) ) / dtime
1497         d_v(i,k) = ( v_seri(i,k) - v(i,k) ) / dtime
1498         d_t(i,k) = ( t_seri(i,k) - t(i,k) ) / dtime
1499      ENDDO
1500      ENDDO
1501c     print*,"vnatphy=",v(705,:)
1502c     print*,"unatphy=",u(705,:)
1503c
1504      DO iq = 1, nqmax
1505      DO  k = 1, klev
1506      DO  i = 1, klon
1507         d_qx(i,k,iq) = ( tr_seri(i,k,iq) - qx(i,k,iq) ) / dtime
1508      ENDDO
1509      ENDDO
1510      ENDDO
1511     
1512c------------------------
1513c Calcul moment cinetique
1514c------------------------
1515c TEST
1516c     mangtot = 0.0
1517c     DO k = 1, klev
1518c     DO i = 1, klon
1519c       mang(i,k) = RA*cos(rlatd(i)*RPI/180.)
1520c    .     *(u_seri(i,k)+RA*cos(rlatd(i)*RPI/180.)*ROMEGA)
1521c    .     *airephy(i)*delp(i,k)/RG
1522c       mangtot=mangtot+mang(i,k)
1523c     ENDDO
1524c     ENDDO
1525c     print*,"Moment cinetique total = ",mangtot
1526c
1527c------------------------
1528c
1529c Sauvegarder les valeurs de t et u a la fin de la physique:
1530c
1531      DO k = 1, klev
1532      DO i = 1, klon
1533         u_ancien(i,k) = u_seri(i,k)
1534         t_ancien(i,k) = t_seri(i,k)
1535      ENDDO
1536      ENDDO
1537c
1538c=============================================================
1539c   Ecriture des sorties
1540c=============================================================
1541     
1542#ifdef CPP_IOIPSL
1543
1544#ifdef histmth
1545#include "write_histmth.h"
1546#endif
1547
1548#ifdef histday
1549#include "write_histday.h"
1550#endif
1551
1552#ifdef histins
1553#include "write_histins.h"
1554#endif
1555             
1556#endif
1557
1558c ---------- Sorties testphys1d -------------
1559     
1560      if (klon.eq.1) then
1561        call writeg1d(klon,klev,t_seri,"Temp","Temperature")
1562        call writeg1d(klon,1,ftsol,"tsurf","Surface Temp")
1563      DO k = 1, klev
1564      DO i = 1, klon
1565c        tmpout(i,k) = real(heat(i,k))
1566         tmpout(i,k) = heat(i,k)
1567      ENDDO
1568      ENDDO
1569        call writeg1d(klon,klev,tmpout,
1570     .                 "heat","Solar heating")
1571      DO k = 1, klev
1572      DO i = 1, klon
1573c        tmpout(i,k) = real(dtrad(i,k))
1574         tmpout(i,k) = dtrad(i,k)
1575      ENDDO
1576      ENDDO
1577        call writeg1d(klon,klev,tmpout,
1578     .                 "dtrad","IR cooling")
1579        call writeg1d(klon,klev,lwnet,"lwnet","Net LW flux")
1580        call writeg1d(klon,klev,swnet,"swnet","Net SW flux")
1581        call writeg1d(klon,klev,fluxt,"flux_vdf","Turbulent flux")
1582        call writeg1d(klon,klev,flux_ajs,"flux_ajs","Dry adjust. flux")
1583        call writeg1d(klon,klev,flux_ec,"flux_ec","Ec flux")
1584c       call writeg1d(klon,1,solsw,"surfsw","Net SW flux at surface")
1585c       call writeg1d(klon,1,sollw,"surflw","Net LW flux at surface")
1586        call writeg1d(klon,1,radsol,"surfnet","Net flux at surface")
1587        call writeg1d(klon,klev,d_t_vdf,"dt_vdf","DT from clmain")
1588        call writeg1d(klon,klev,d_t_ajs,"dt_ajs","DT from ajsec")
1589        call writeg1d(klon,klev,d_t_ec,"dt_ec","DT from Ec")
1590      endif
1591     
1592c====================================================================
1593c Si c'est la fin, il faut conserver l'etat de redemarrage
1594c====================================================================
1595c
1596      IF (lafin) THEN
1597         itau_phy = itau_phy + itap
1598c REMETTRE TOUS LES PARAMETRES POUR OROGW... A FAIRE POUR TITAN
1599         CALL phyredem ("restartphy.nc",dtime,radpas,chimpas,
1600     .      rlatd, rlond, ftsol, ftsoil,
1601     .      falbe,
1602     .      solsw, sollw,dlw,
1603     .      radsol,reservoir,
1604c     .      zmea,zstd,zsig,zgam,zthe,zpic,zval,
1605     .      t_ancien)
1606     
1607c--------------FLOTT
1608CMODEB LOTT
1609C  FERMETURE DU FICHIER FORMATTE CONTENANT LES COMPOSANTES
1610C  DU BILAN DE MOMENT ANGULAIRE.
1611      if (bilansmc.eq.1) then
1612        write(*,*)'Fermeture de aaam_bud.out (FL Vous parle)'
1613        close(27)                                     
1614        close(28)                                     
1615      endif !bilansmc
1616CMODFIN
1617c-------------
1618c--------------SLEBONNOIS
1619C  FERMETURE DES FICHIERS FORMATTES CONTENANT LES POSITIONS ET VITESSES
1620C  DES BALLONS
1621      if (ballons.eq.1) then
1622        write(*,*)'Fermeture des ballons*.out'
1623        close(30)                                     
1624        close(31)                                     
1625        close(32)                                     
1626        close(33)                                     
1627        close(34)                                     
1628      endif !ballons
1629c-------------
1630      ENDIF
1631     
1632
1633      RETURN
1634      END
1635
1636
1637
1638***********************************************************************
1639***********************************************************************
1640***********************************************************************
1641***********************************************************************
1642***********************************************************************
1643***********************************************************************
1644***********************************************************************
1645***********************************************************************
1646
1647      SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit)
1648      IMPLICIT none
1649c
1650c Tranformer une variable de la grille physique a
1651c la grille d'ecriture
1652c
1653      INTEGER nfield,nlon,iim,jjmp1, jjm
1654      REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield)
1655c
1656      INTEGER i, n, ig
1657c
1658      jjm = jjmp1 - 1
1659      DO n = 1, nfield
1660         DO i=1,iim
1661            ecrit(i,n) = fi(1,n)
1662            ecrit(i+jjm*iim,n) = fi(nlon,n)
1663         ENDDO
1664         DO ig = 1, nlon - 2
1665           ecrit(iim+ig,n) = fi(1+ig,n)
1666         ENDDO
1667      ENDDO
1668      RETURN
1669      END
Note: See TracBrowser for help on using the repository browser.