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

Last change on this file since 778 was 495, checked in by slebonnois, 14 years ago

Mise a jour physique Titan, ajout des forces de marees (dans la dynamique, sous flag titan). SL.

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