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

Last change on this file since 824 was 815, checked in by slebonnois, 13 years ago

SL: petites modifs Titan et Venus pour tableau controle dans la physique ; pour Titan, petits details lies a raz_date ; modif chemin ioipsl sur gnome ; + elimination d'un warning etrange dans gcm.F

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