source: trunk/LMDZ.VENUS/libf/phyvenus/physiq.F @ 1460

Last change on this file since 1460 was 1460, checked in by slebonnois, 9 years ago

SL: bug dans la prise en compte de d_tr_vdf dans la physique (Venus)

File size: 53.3 KB
RevLine 
[3]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,
[97]7     .            paprs,pplay,ppk,pphi,pphis,presnivs,
[3]8     .            u,v,t,qx,
[1301]9     .            flxmw,
[3]10     .            d_u, d_v, d_t, d_qx, d_ps)
11
12c======================================================================
13c
14c Modifications pour la physique de Venus
15c     S. Lebonnois (LMD/CNRS) Septembre 2005
16c
17c ---------------------------------------------------------------------
18c Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
19c
20c Objet: Moniteur general de la physique du modele
21cAA      Modifications quant aux traceurs :
22cAA                  -  uniformisation des parametrisations ds phytrac
23cAA                  -  stockage des moyennes des champs necessaires
24cAA                     en mode traceur off-line
25c    modif   ( P. Le Van ,  12/10/98 )
26c
27c  Arguments:
28c
29c nlon----input-I-nombre de points horizontaux
30c nlev----input-I-nombre de couches verticales
31c nqmax---input-I-nombre de traceurs
32c debut---input-L-variable logique indiquant le premier passage
33c lafin---input-L-variable logique indiquant le dernier passage
34c rjour---input-R-numero du jour de l'experience
[953]35c gmtime--input-R-fraction de la journee (0 a 1)
[3]36c pdtphys-input-R-pas d'integration pour la physique (seconde)
37c paprs---input-R-pression pour chaque inter-couche (en Pa)
38c pplay---input-R-pression pour le mileu de chaque couche (en Pa)
39c ppk  ---input-R-fonction d'Exner au milieu de couche
40c pphi----input-R-geopotentiel de chaque couche (g z) (reference sol)
41c pphis---input-R-geopotentiel du sol
42c presnivs-input_R_pressions approximat. des milieux couches ( en PA)
43c u-------input-R-vitesse dans la direction X (de O a E) en m/s
44c v-------input-R-vitesse Y (de S a N) en m/s
45c t-------input-R-temperature (K)
46c qx------input-R-mass mixing ratio traceurs (kg/kg)
47c d_t_dyn-input-R-tendance dynamique pour "t" (K/s)
[1301]48c flxmw---input-R-flux de masse vertical en kg/s
[3]49c
50c d_u-----output-R-tendance physique de "u" (m/s/s)
51c d_v-----output-R-tendance physique de "v" (m/s/s)
52c d_t-----output-R-tendance physique de "t" (K/s)
53c d_qx----output-R-tendance physique de "qx" (kg/kg/s)
54c d_ps----output-R-tendance physique de la pression au sol
55c======================================================================
[101]56      USE ioipsl
[849]57!      USE histcom ! not needed; histcom is included in ioipsl
[101]58      USE infotrac
59      USE control_mod
60      use dimphy
61      USE comgeomphy
[892]62      USE mod_phys_lmdz_para, only : is_parallel,jj_nb
63      USE phys_state_var_mod ! Variables sauvegardees de la physique
[973]64      USE write_field_phy
[892]65      USE iophy
[1442]66      USE cpdet_mod, only: cpdet, t2tpot
[1310]67      USE chemparam_mod
68      USE conc
69      USE compo_hedin83_mod2
70      use moyzon_mod, only: tmoy
[1305]71      use ieee_arithmetic
[1443]72      use temps_mod, only: itau_phy,day_ref,annee_ref
73      use logic_mod, only: iflag_trac
[101]74      IMPLICIT none
75c======================================================================
76c   CLEFS CPP POUR LES IO
77c   =====================
78c#define histhf
79#define histday
80#define histmth
81#define histins
82c======================================================================
[3]83#include "dimensions.h"
84      integer jjmp1
85      parameter (jjmp1=jjm+1-1/jjm)
86#include "dimsoil.h"
87#include "clesphys.h"
88#include "iniprint.h"
89#include "timerad.h"
[815]90#include "tabcontrol.h"
[1310]91#include "nirdata.h"
92#include "hedin.h"
[3]93c======================================================================
94      LOGICAL ok_journe ! sortir le fichier journalier
95      save ok_journe
96c      PARAMETER (ok_journe=.true.)
97c
98      LOGICAL ok_mensuel ! sortir le fichier mensuel
99      save ok_mensuel
100c      PARAMETER (ok_mensuel=.true.)
101c
102      LOGICAL ok_instan ! sortir le fichier instantane
103      save ok_instan
104c      PARAMETER (ok_instan=.true.)
105c
106c======================================================================
107c
108c Variables argument:
109c
110      INTEGER nlon
111      INTEGER nlev
112      INTEGER nqmax
113      REAL rjourvrai
114      REAL gmtime
115      REAL pdtphys
116      LOGICAL debut, lafin
117      REAL paprs(klon,klev+1)
118      REAL pplay(klon,klev)
119      REAL pphi(klon,klev)
120      REAL pphis(klon)
121      REAL presnivs(klev)
122
123! ADAPTATION GCM POUR CP(T)
124      REAL ppk(klon,klev)
125
126      REAL u(klon,klev)
127      REAL v(klon,klev)
128      REAL t(klon,klev)
129      REAL qx(klon,klev,nqmax)
130
131      REAL d_u_dyn(klon,klev)
132      REAL d_t_dyn(klon,klev)
133
[1301]134      REAL flxmw(klon,klev)
[3]135
136      REAL d_u(klon,klev)
137      REAL d_v(klon,klev)
138      REAL d_t(klon,klev)
139      REAL d_qx(klon,klev,nqmax)
140      REAL d_ps(klon)
141
142      logical ok_hf
143      real ecrit_hf
144      integer nid_hf
145      save ok_hf, ecrit_hf, nid_hf
146
147#ifdef histhf
148      data ok_hf,ecrit_hf/.true.,0.25/
149#else
150      data ok_hf/.false./
151#endif
152
153c Variables propres a la physique
154c
[101]155      INTEGER,save :: itap        ! compteur pour la physique
[3]156      REAL delp(klon,klev)        ! epaisseur d'une couche
[1301]157      REAL omega(klon,klev)       ! vitesse verticale en Pa/s
158
[3]159     
160      INTEGER igwd,idx(klon),itest(klon)
161c
162c  Diagnostiques 2D de drag_noro, lift_noro et gw_nonoro
163
164      REAL zulow(klon),zvlow(klon)
165      REAL zustrdr(klon), zvstrdr(klon)
166      REAL zustrli(klon), zvstrli(klon)
167      REAL zustrhi(klon), zvstrhi(klon)
168
169c Pour calcul GW drag oro et nonoro: CALCUL de N2:
170      real zdzlev(klon,klev)
171      real ztlev(klon,klev),zpklev(klon,klev)
172      real ztetalay(klon,klev),ztetalev(klon,klev)
173      real zdtetalev(klon,klev)
174      real zn2(klon,klev) ! BV^2 at plev
175
176c Pour les bilans de moment angulaire,
177      integer bilansmc
178c Pour le transport de ballons
179      integer ballons
180c j'ai aussi besoin
181c du stress de couche limite a la surface:
182
183      REAL zustrcl(klon),zvstrcl(klon)
184
185c et du stress total c de la physique:
186
187      REAL zustrph(klon),zvstrph(klon)
188
189c Variables locales:
190c
191      REAL cdragh(klon) ! drag coefficient pour T and Q
192      REAL cdragm(klon) ! drag coefficient pour vent
193c
194cAA  Pour  TRACEURS
195cAA
[105]196      REAL,save,allocatable :: source(:,:)
[3]197      REAL ycoefh(klon,klev)    ! coef d'echange pour phytrac
198      REAL yu1(klon)            ! vents dans la premiere couche U
199      REAL yv1(klon)            ! vents dans la premiere couche V
200
201      REAL sens(klon), dsens(klon) ! chaleur sensible et sa derivee
202      REAL ve(klon) ! integr. verticale du transport meri. de l'energie
203      REAL vq(klon) ! integr. verticale du transport meri. de l'eau
204      REAL ue(klon) ! integr. verticale du transport zonal de l'energie
205      REAL uq(klon) ! integr. verticale du transport zonal de l'eau
206c
[1305]207      REAL Fsedim(klon,klev+1)  ! Flux de sedimentation (kg.m-2)
[3]208
209c======================================================================
210c
211c Declaration des procedures appelees
212c
213      EXTERNAL ajsec     ! ajustement sec
214      EXTERNAL clmain    ! couche limite
215      EXTERNAL hgardfou  ! verifier les temperatures
216c     EXTERNAL orbite    ! calculer l'orbite
217      EXTERNAL phyetat0  ! lire l'etat initial de la physique
218      EXTERNAL phyredem  ! ecrire l'etat de redemarrage de la physique
219      EXTERNAL radlwsw   ! rayonnements solaire et infrarouge
220      EXTERNAL suphec    ! initialiser certaines constantes
221      EXTERNAL transp    ! transport total de l'eau et de l'energie
222      EXTERNAL abort_gcm
223      EXTERNAL printflag
224      EXTERNAL zenang
225      EXTERNAL diagetpq
226      EXTERNAL conf_phys
227      EXTERNAL diagphy
228      EXTERNAL mucorr
229      EXTERNAL phytrac
[1310]230      EXTERNAL nirco2abs
231      EXTERNAL nir_leedat
232      EXTERNAL nltecool
233      EXTERNAL nlte_tcool
234      EXTERNAL nlte_setup
235      EXTERNAL blendrad
236      EXTERNAL nlthermeq
237      EXTERNAL euvheat
238      EXTERNAL param_read
239      EXTERNAL param_read_e107
240      EXTERNAL conduction
241      EXTERNAL molvis
[1442]242      EXTERNAL moldiff_red
243
[3]244c
245c Variables locales
246c
247CXXX PB
248      REAL fluxt(klon,klev)   ! flux turbulent de chaleur
249      REAL fluxu(klon,klev)   ! flux turbulent de vitesse u
250      REAL fluxv(klon,klev)   ! flux turbulent de vitesse v
251c
252      REAL flux_dyn(klon,klev)  ! flux de chaleur produit par la dynamique
253      REAL flux_ajs(klon,klev)  ! flux de chaleur ajustement sec
254      REAL flux_ec(klon,klev)   ! flux de chaleur Ec
255c
256      REAL tmpout(klon,klev)  ! K s-1
257
258      INTEGER itaprad
259      SAVE itaprad
260c
261      REAL dist, rmu0(klon), fract(klon)
262      REAL zdtime, zlongi
263c
[1442]264      INTEGER i, k, iq, ig, j, ll, ilon, ilat, ilev
[3]265c
266      REAL zphi(klon,klev)
[1301]267      REAL zzlev(klon,klev+1),zzlay(klon,klev),z1,z2
[1310]268      real tsurf(klon)
[1301]269
[1310]270c va avec nlte_tcool
271      INTEGER ierr_nlte
272      REAL    varerr
273
[3]274c Variables du changement
275c
276c ajs: ajustement sec
277c vdf: couche limite (Vertical DiFfusion)
278      REAL d_t_ajs(klon,klev), d_tr_ajs(klon,klev,nqmax)
279      REAL d_u_ajs(klon,klev), d_v_ajs(klon,klev)
280c
281      REAL d_ts(klon)
282c
283      REAL d_u_vdf(klon,klev), d_v_vdf(klon,klev)
284      REAL d_t_vdf(klon,klev), d_tr_vdf(klon,klev,nqmax)
285c
286CMOD LOTT: Tendances Orography Sous-maille
287      REAL d_u_oro(klon,klev), d_v_oro(klon,klev)
288      REAL d_t_oro(klon,klev)
289      REAL d_u_lif(klon,klev), d_v_lif(klon,klev)
290      REAL d_t_lif(klon,klev)
291C          Tendances Ondes de G non oro (runs strato).
292      REAL d_u_hin(klon,klev), d_v_hin(klon,klev)
293      REAL d_t_hin(klon,klev)
294
[1310]295c Tendencies due to radiative scheme   [K/s]
296c     d_t_rad,dtsw,dtlw,d_t_nirco2,d_t_nlte,d_t_euv
297c are not computed at each physical timestep
298c therefore, they are defined and saved in phys_state_var_mod
299
300c Tendencies due to molecular viscosity and conduction
301      real d_t_conduc(klon,klev)     ! [K/s]
302      real d_u_molvis(klon,klev)     ! (m/s) /s
303      real d_v_molvis(klon,klev)     ! (m/s) /s
304
[1442]305c Tendencies due to molecular diffusion
306      real d_q_moldif(klon,klev,nqmax)
307
[3]308c
309c Variables liees a l'ecriture de la bande histoire physique
310c
311      INTEGER ecrit_mth
312      SAVE ecrit_mth   ! frequence d'ecriture (fichier mensuel)
313c
314      INTEGER ecrit_day
315      SAVE ecrit_day   ! frequence d'ecriture (fichier journalier)
316c
317      INTEGER ecrit_ins
318      SAVE ecrit_ins   ! frequence d'ecriture (fichier instantane)
319c
320      integer itau_w   ! pas de temps ecriture = itap + itau_phy
321
322c Variables locales pour effectuer les appels en serie
323c
324      REAL t_seri(klon,klev)
325      REAL u_seri(klon,klev), v_seri(klon,klev)
326c
[1305]327      REAL :: tr_seri(klon,klev,nqmax)
328      REAL :: d_tr(klon,klev,nqmax)
329
[1442]330c Champ de modification de la temperature par rapport a VIRAII
331      REAL delta_temp(klon,klev)
332c      SAVE delta_temp
333      REAL mat_dtemp(33,50)
334      SAVE mat_dtemp
335
[1305]336c Variables tendance sedimentation
337
338      REAL :: d_tr_sed(klon,klev,2)
339      REAL :: d_tr_ssed(klon)
[3]340c
341c pour ioipsl
342      INTEGER nid_day, nid_mth, nid_ins
343      SAVE nid_day, nid_mth, nid_ins
344      INTEGER nhori, nvert, idayref
[97]345      REAL zsto, zout, zsto1, zsto2, zero
[3]346      parameter (zero=0.0e0)
347      real zjulian
348      save zjulian
349
350      CHARACTER*2  str2
351      character*20 modname
352      character*80 abort_message
353      logical ok_sync
354
355      character*30 nom_fichier
356      character*10 varname
357      character*40 vartitle
358      character*20 varunits
359C     Variables liees au bilan d'energie et d'enthalpi
360      REAL ztsol(klon)
361      REAL      h_vcol_tot, h_dair_tot, h_qw_tot, h_ql_tot
362     $        , h_qs_tot, qw_tot, ql_tot, qs_tot , ec_tot
363      SAVE      h_vcol_tot, h_dair_tot, h_qw_tot, h_ql_tot
364     $        , h_qs_tot, qw_tot, ql_tot, qs_tot , ec_tot
365      REAL      d_h_vcol, d_h_dair, d_qt, d_qw, d_ql, d_qs, d_ec
366      REAL      d_h_vcol_phy
367      REAL      fs_bound, fq_bound
368      SAVE      d_h_vcol_phy
369      REAL      zero_v(klon),zero_v2(klon,klev)
370      CHARACTER*15 ztit
371      INTEGER   ip_ebil  ! PRINT level for energy conserv. diag.
372      SAVE      ip_ebil
373      DATA      ip_ebil/2/
374      INTEGER   if_ebil ! level for energy conserv. dignostics
375      SAVE      if_ebil
376c+jld ec_conser
377      REAL d_t_ec(klon,klev)    ! tendance du a la conversion Ec -> E thermique
378c-jld ec_conser
379
380c TEST VENUS...
381      REAL mang(klon,klev)    ! moment cinetique
382      REAL mangtot            ! moment cinetique total
383
384c Declaration des constantes et des fonctions thermodynamiques
385c
386#include "YOMCST.h"
387
388c======================================================================
389c INITIALISATIONS
390c================
391
392      modname = 'physiq'
393      ok_sync=.TRUE.
394
[119]395      bilansmc = 0
[3]396      ballons  = 0
[892]397! NE FONCTIONNENT PAS ENCORE EN PARALLELE !!!
398      if (is_parallel) then
399        bilansmc = 0
400        ballons  = 0
401      endif
[3]402
403      IF (if_ebil.ge.1) THEN
404        DO i=1,klon
405          zero_v(i)=0.
406        END DO
407        DO i=1,klon
408         DO j=1,klev
409          zero_v2(i,j)=0.
410         END DO
411        END DO
412      END IF
413     
414c PREMIER APPEL SEULEMENT
415c========================
416      IF (debut) THEN
[105]417         allocate(source(klon,nqmax))
[3]418
419         CALL suphec ! initialiser constantes et parametres phys.
420
421         IF (if_ebil.ge.1) d_h_vcol_phy=0.
422c
423c appel a la lecture du physiq.def
424c
425         call conf_phys(ok_journe, ok_mensuel,
426     .                  ok_instan,
427     .                  if_ebil)
428
[892]429         call phys_state_var_init
[3]430c
[1310]431c Initialising Hedin model for upper atm
432c   (to be revised when coupled to chemistry) :
433         call conc_init
434c
[3]435c Initialiser les compteurs:
436c
437         itap    = 0
438         itaprad = 0
439c         
440c Lecture startphy.nc :
441c
[892]442         CALL phyetat0 ("startphy.nc")
[3]443
[815]444c dtime est defini dans tabcontrol.h et lu dans startphy
[150]445c pdtphys est calcule a partir des nouvelles conditions:
446c Reinitialisation du pas de temps physique quand changement
447         IF (ABS(dtime-pdtphys).GT.0.001) THEN
448            WRITE(lunout,*) 'Pas physique a change',dtime,
449     .                        pdtphys
450c           abort_message='Pas physique n est pas correct '
451c           call abort_gcm(modname,abort_message,1)
[152]452c----------------
453c pour initialiser convenablement le time_counter, il faut tenir compte
454c du changement de dtime en changeant itau_phy (point de depart)
455            itau_phy = NINT(itau_phy*dtime/pdtphys)
456c----------------
[150]457            dtime=pdtphys
458         ENDIF
459
[3]460         radpas = NINT( RDAY/pdtphys/nbapp_rad)
461
[815]462         CALL printflag( ok_journe,ok_instan )
[3]463c
464c---------
465c FLOTT
466       IF (ok_orodr) THEN
467         DO i=1,klon
468         rugoro(i) = MAX(1.0e-05, zstd(i)*zsig(i)/2.0)
469         ENDDO
470         CALL SUGWD(klon,klev,paprs,pplay)
471         DO i=1,klon
472         zuthe(i)=0.
473         zvthe(i)=0.
474         if(zstd(i).gt.10.)then
475           zuthe(i)=(1.-zgam(i))*cos(zthe(i))
476           zvthe(i)=(1.-zgam(i))*sin(zthe(i))
477         endif
478         ENDDO
479       ENDIF
480
481      if (bilansmc.eq.1) then
482C  OUVERTURE D'UN FICHIER FORMATTE POUR STOCKER LES COMPOSANTES
483C  DU BILAN DE MOMENT ANGULAIRE.
484      open(27,file='aaam_bud.out',form='formatted')
485      open(28,file='fields_2d.out',form='formatted')
486      write(*,*)'Ouverture de aaam_bud.out (FL Vous parle)'
487      write(*,*)'Ouverture de fields_2d.out (FL Vous parle)'
488      endif !bilansmc
489
490c--------------SLEBONNOIS
491C  OUVERTURE DES FICHIERS FORMATTES CONTENANT LES POSITIONS ET VITESSES
492C  DES BALLONS
493      if (ballons.eq.1) then
494      open(30,file='ballons-lat.out',form='formatted')
495      open(31,file='ballons-lon.out',form='formatted')
496      open(32,file='ballons-u.out',form='formatted')
497      open(33,file='ballons-v.out',form='formatted')
498      open(34,file='ballons-alt.out',form='formatted')
499      write(*,*)'Ouverture des ballons*.out'
500      endif !ballons
501c-------------
502
503c---------
504C TRACEURS
505C source dans couche limite
[1442]506         source(:,:) = 0.0 ! pas de source, pour l'instant
[3]507c---------
[1310]508
509c---------
510c INITIALIZE THERMOSPHERIC PARAMETERS
511c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
512
513         if (callthermos) then
514            if(solvarmod.eq.0) call param_read
515            if(solvarmod.eq.1) call param_read_e107
516         endif
517
518c Initialisation (recomputed in concentration2)
519       do ig=1,klon
520         do j=1,klev
521            rnew(ig,j)=R
522            cpnew(ig,j)=cpdet(tmoy(j))
[1442]523            mmean(ig,j)=RMD
[1310]524            akknew(ig,j)=1.e-4
525          enddo
526c        stop
527
528        enddo 
529     
530      IF(callthermos.or.callnlte.or.callnirco2) THEN 
531         call compo_hedin83_init2
532      ENDIF
[1442]533      if (callnlte.and.nltemodel.eq.2) call nlte_setup
534      if (callnirco2.and.nircorr.eq.1) call nir_leedat         
[1310]535c---------
536     
[3]537c
538c Verifications:
539c
540         IF (nlon .NE. klon) THEN
541            WRITE(lunout,*)'nlon et klon ne sont pas coherents', nlon,
542     .                      klon
543            abort_message='nlon et klon ne sont pas coherents'
544            call abort_gcm(modname,abort_message,1)
545         ENDIF
546         IF (nlev .NE. klev) THEN
547            WRITE(lunout,*)'nlev et klev ne sont pas coherents', nlev,
548     .                       klev
549            abort_message='nlev et klev ne sont pas coherents'
550            call abort_gcm(modname,abort_message,1)
551         ENDIF
552c
[808]553         IF (dtime*REAL(radpas).GT.(RDAY*0.25).AND.cycle_diurne)
[3]554     $    THEN
555           WRITE(lunout,*)'Nbre d appels au rayonnement insuffisant'
556           WRITE(lunout,*)"Au minimum 4 appels par jour si cycle diurne"
557           abort_message='Nbre d appels au rayonnement insuffisant'
558           call abort_gcm(modname,abort_message,1)
559         ENDIF
560c
561         WRITE(lunout,*)"Clef pour la convection seche, iflag_ajs=",
562     .                   iflag_ajs
563c
[1120]564         ecrit_mth = NINT(RDAY/dtime*ecriphy)  ! tous les ecritphy jours
[3]565         IF (ok_mensuel) THEN
566         WRITE(lunout,*)'La frequence de sortie mensuelle est de ',
567     .                   ecrit_mth
568         ENDIF
[97]569
[3]570         ecrit_day = NINT(RDAY/dtime *1.0)  ! tous les jours
571         IF (ok_journe) THEN
572         WRITE(lunout,*)'La frequence de sortie journaliere est de ',
573     .                   ecrit_day
574         ENDIF
[97]575
[1120]576         ecrit_ins = NINT(RDAY/dtime*ecriphy)  ! Fraction de jour reglable
[3]577         IF (ok_instan) THEN
578         WRITE(lunout,*)'La frequence de sortie instant. est de ',
579     .                   ecrit_ins
580         ENDIF
581
582c Initialisation des sorties
583c===========================
584
585#ifdef CPP_IOIPSL
586
587#ifdef histhf
588#include "ini_histhf.h"
589#endif
590
591#ifdef histday
592#include "ini_histday.h"
593#endif
594
595#ifdef histmth
596#include "ini_histmth.h"
597#endif
598
599#ifdef histins
600#include "ini_histins.h"
601#endif
602
603#endif
604
605c
606c Initialiser les valeurs de u pour calculs tendances
607c (pour T, c'est fait dans phyetat0)
608c
609      DO k = 1, klev
610      DO i = 1, klon
611         u_ancien(i,k) = u(i,k)
612      ENDDO
613      ENDDO
614
[1305]615c---------
616c       Ecriture fichier initialisation
617c       PRINT*,'Ecriture Initial_State.csv'
618c       OPEN(88,file='Trac_Point.csv',
619c     & form='formatted')
620c---------
621     
622c---------
623c       Initialisation des parametres des nuages
624c===============================================
625     
626      if ((nlon .EQ. 1) .AND. ok_cloud) then
627        PRINT*,'Open profile_cloud_parameters.csv'
628        OPEN(66,file='profile_cloud_parameters.csv',
629     &   form='formatted')
630      endif
631
632      if ((nlon .EQ. 1) .AND. ok_sedim) then
633        PRINT*,'Open profile_cloud_sedim.csv'
634        OPEN(77,file='profile_cloud_sedim.csv',
635     &   form='formatted')
636      endif
637           
[1442]638      if ((nlon .GT. 1) .AND. (ok_chem.OR.ok_cloud)) then
[1305]639c !!! DONC 3D !!!
640        CALL chemparam_ini()
641      endif
[3]642         
[1305]643      if ((nlon .GT. 1) .AND. ok_cloud) then
644c !!! DONC 3D !!!
645        CALL cloud_ini(nlon,nlev)
646      endif
[1442]647
648c======================================================================
649c      Lecture du fichier DeltaT
650c======================================================================
651
652c     ATTENTION tout ce qui suit est pour un 48*32*50
653
654      if (ok_deltatemp) then
655
656      print*,'lecture de VenusDeltaT.txt '
657      open(99, form = 'formatted', status = 'old', file =
658     &     'VenusDeltaT.dat')
659      print*,'Ouverture de VenusDeltaT.txt '
660
661      DO ilev = 1, klev
662         read(99,'(33(1x,e13.6))') (mat_dtemp(ilat,ilev),ilat=1,33)
663         print*,'lecture de VenusDeltaT.txt ligne:',ilev
664      ENDDO
665     
666      close(99)
667      print*,'FIN lecture de VenusDeltaT.txt ok.'
668
669      DO k = 1, klev
670      DO i = 1, klon     
671         ilat=(rlatd(i)/5.625) + 17.
672         delta_temp(i,k)=mat_dtemp(INT(ilat),k)
673      ENDDO
674      ENDDO
675     
676      endif
[1305]677       
[3]678      ENDIF ! debut
679c======================================================================
[1442]680c======================================================================
[3]681
682c Mettre a zero des variables de sortie (pour securite)
683c
684      DO i = 1, klon
685         d_ps(i) = 0.0
686      ENDDO
687      DO k = 1, klev
688      DO i = 1, klon
689         d_t(i,k) = 0.0
690         d_u(i,k) = 0.0
691         d_v(i,k) = 0.0
692      ENDDO
693      ENDDO
694      DO iq = 1, nqmax
695      DO k = 1, klev
696      DO i = 1, klon
697         d_qx(i,k,iq) = 0.0
698      ENDDO
699      ENDDO
700      ENDDO
701c
702c Ne pas affecter les valeurs entrees de u, v, h, et q
703c
704      DO k = 1, klev
705      DO i = 1, klon
706         t_seri(i,k)  = t(i,k)
707         u_seri(i,k)  = u(i,k)
708         v_seri(i,k)  = v(i,k)
709      ENDDO
710      ENDDO
711      DO iq = 1, nqmax
712      DO  k = 1, klev
713      DO  i = 1, klon
714         tr_seri(i,k,iq) = qx(i,k,iq)
715      ENDDO
716      ENDDO
717      ENDDO
718C
719      DO i = 1, klon
720          ztsol(i) = ftsol(i)
721      ENDDO
722C
723      IF (if_ebil.ge.1) THEN
724        ztit='after dynamic'
725        CALL diagetpq(airephy,ztit,ip_ebil,1,1,dtime
726     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
727     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
728C     Comme les tendances de la physique sont ajoute dans la dynamique,
729C     on devrait avoir que la variation d'entalpie par la dynamique
730C     est egale a la variation de la physique au pas de temps precedent.
731C     Donc la somme de ces 2 variations devrait etre nulle.
732        call diagphy(airephy,ztit,ip_ebil
733     e      , zero_v, zero_v, zero_v, zero_v, zero_v
734     e      , zero_v, zero_v, zero_v, ztsol
735     e      , d_h_vcol+d_h_vcol_phy, d_qt, 0.
736     s      , fs_bound, fq_bound )
737      END IF
738
739c====================================================================
740c Diagnostiquer la tendance dynamique
741c
742      IF (ancien_ok) THEN
743         DO k = 1, klev
744         DO i = 1, klon
745            d_u_dyn(i,k) = (u_seri(i,k)-u_ancien(i,k))/dtime
746            d_t_dyn(i,k) = (t_seri(i,k)-t_ancien(i,k))/dtime
747         ENDDO
748         ENDDO
749
750! ADAPTATION GCM POUR CP(T)
751         do i=1,klon
752          flux_dyn(i,1) = 0.0
753          do j=2,klev
754            flux_dyn(i,j) = flux_dyn(i,j-1)
755     . +cpdet(t_seri(i,j-1))/RG*d_t_dyn(i,j-1)*(paprs(i,j-1)-paprs(i,j))
756          enddo
757         enddo
758         
759      ELSE
760         DO k = 1, klev
761         DO i = 1, klon
762            d_u_dyn(i,k) = 0.0
763            d_t_dyn(i,k) = 0.0
764         ENDDO
765         ENDDO
766         ancien_ok = .TRUE.
767      ENDIF
768c====================================================================
[1301]769
770c Calcule de vitesse verticale a partir de flux de masse verticale
771      DO k = 1, klev
772       DO i = 1, klon
773        omega(i,k) = RG*flxmw(i,k) / airephy(i)
774       END DO
775      END DO
776
[3]777c
778c Ajouter le geopotentiel du sol:
779c
780      DO k = 1, klev
781      DO i = 1, klon
782         zphi(i,k) = pphi(i,k) + pphis(i)
783      ENDDO
784      ENDDO
[1301]785
786c   calcul du geopotentiel aux niveaux intercouches
787c   ponderation des altitudes au niveau des couches en dp/p
788
789      DO k=1,klev
790         DO i=1,klon
[1310]791            zzlay(i,k)=zphi(i,k)/RG        ! [m]
[1301]792         ENDDO
793      ENDDO
794      DO i=1,klon
[1310]795         zzlev(i,1)=pphis(i)/RG            ! [m]
[1301]796      ENDDO
797      DO k=2,klev
798         DO i=1,klon
799            z1=(pplay(i,k-1)+paprs(i,k))/(pplay(i,k-1)-paprs(i,k))
800            z2=(paprs(i,k)  +pplay(i,k))/(paprs(i,k)  -pplay(i,k))
801            zzlev(i,k)=(z1*zzlay(i,k-1)+z2*zzlay(i,k))/(z1+z2)
802         ENDDO
803      ENDDO
804      DO i=1,klon
805         zzlev(i,klev+1)=zzlay(i,klev)+(zzlay(i,klev)-zzlev(i,klev))
806      ENDDO
807
[3]808c====================================================================
809c
810c Verifier les temperatures
811c
812      CALL hgardfou(t_seri,ftsol,'debutphy')
813c====================================================================
814c
815c Incrementer le compteur de la physique
816c
817      itap   = itap + 1
818
819c====================================================================
820c Orbite et eclairement
821c====================================================================
822
823c Pour VENUS, on fixe l'obliquite a 0 et l'eccentricite a 0.
824c donc pas de variations de Ls, ni de dist.
825c La seule chose qui compte, c'est la rotation de la planete devant
826c le Soleil...
827     
828      zlongi = 0.0
829      dist   = 0.72  ! en UA
830
831c Si on veut remettre l'obliquite a 3 degres et/ou l'eccentricite
832c a sa valeur, et prendre en compte leur evolution,
833c il faudra refaire un orbite.F...
834c     CALL orbite(zlongi,dist)
835
836      IF (cycle_diurne) THEN
[808]837        zdtime=dtime*REAL(radpas) ! pas de temps du rayonnement (s)
[3]838        CALL zenang(zlongi,gmtime,zdtime,rlatd,rlond,rmu0,fract)
839      ELSE
840        call mucorr(klon,zlongi,rlatd,rmu0,fract)
841      ENDIF
842     
843c====================================================================
844c   Calcul  des tendances traceurs
845c====================================================================
846
847      if (iflag_trac.eq.1) then
[1160]848
849       if (tr_scheme.eq.1) then
850! Case 1: pseudo-chemistry with relaxation toward fixed profile
[1305]851
[1160]852         call phytrac_relax (debut,lafin,nqmax,
853     I                   nlon,nlev,dtime,pplay,
[3]854     O                   tr_seri)
[1160]855
856       elseif (tr_scheme.eq.2) then
857! Case 2: surface emission
858! For the moment, inspired from Mars version
859! However, the variable 'source' could be used in physiq
860! so the call to phytrac_emiss could be to initialise it.
[1305]861
[1160]862         call phytrac_emiss ( (rjourvrai+gmtime)*RDAY,
863     I                   debut,lafin,nqmax,
864     I                   nlon,nlev,dtime,paprs,
865     I                   rlatd,rlond,
866     O                   tr_seri)
[3]867
[1305]868       elseif (tr_scheme.eq.3) then  ! identical to ok_chem.or.ok_cloud
869! Case 3: Full chemistry and/or clouds
870
[1442]871         if (ok_deltatemp) then
872!           PRINT*,'Def de delta_temp'
873           DO k = 1, klev
874           DO i = 1, klon     
875             ilat=(rlatd(i)/5.625) + 17.
876!             PRINT*,INT(ilat),rlatd(i),mat_dtemp(INT(ilat),k)
877             delta_temp(i,k)=mat_dtemp(INT(ilat),k)
878           ENDDO
879           ENDDO
880     
881         endif
882
[1452]883!         Pas appel concentrations car BUG ... comprend pas
884!         valeur de mmean trop hautes 5,25E+29
885!         call concentrations2(pplay,t_seri,d_t,tr_seri, nqmax)
886
[1442]887         if (ok_deltatemp) then
888! Utilisation du champ de temperature modifie         
889           call phytrac_chimie(                             
[1305]890     I             debut,
891     I             gmtime,
892     I             nqmax,
[1442]893     I             klon,
[1305]894     I             rlatd,
895     I             rlond,
896     I             nlev,
897     I             dtime,
[1442]898     I             t_seri+delta_temp,
899     I             pplay,
900     O             tr_seri)
901         else
902         
903           call phytrac_chimie(                             
904     I             debut,
905     I             gmtime,
906     I             nqmax,
907     I             klon,
908     I             rlatd,
909     I             rlond,
910     I             nlev,
911     I             dtime,
912     I             t_seri,
913     I             pplay,
914     O             tr_seri)
915         endif
[1305]916
917c        CALL WriteField_phy('Pression',pplay,nlev)
918c        CALL WriteField_phy('PressionBnd',paprs,nlev+1)
919c        CALL WriteField_phy('Temp',t_seri,nlev)
920c        IF (ok_cloud) THEN
921c          CALL WriteField_phy('NBRTOT',NBRTOT,nlev)
922c        ENDIF
923c        CALL WriteField_phy('SAl',tr_seri(:,:,i_h2so4liq),nlev)
924c        CALL WriteField_phy('SAg',tr_seri(:,:,i_h2so4),nlev)
925
926         if (ok_sedim) then
927
[1442]928         if (ok_deltatemp) then
929! Utilisation du champ de temperature modifie 
[1305]930           CALL new_cloud_sedim(
[1442]931     I                 klon,
932     I               nlev,
933     I               dtime,
934     I                 pplay,
935     I               paprs,
936     I               t_seri+delta_temp,
937     I                 tr_seri,
938     O               d_tr_sed,
939     O               d_tr_ssed,
940     I               nqmax,
941     O                 Fsedim)
942         else
943         
944           CALL new_cloud_sedim(
945     I                 klon,
946     I               nlev,
947     I               dtime,
948     I                 pplay,
949     I               paprs,
950     I               t_seri,
951     I                 tr_seri,
952     O               d_tr_sed,
953     O               d_tr_ssed,
954     I               nqmax,
955     O                 Fsedim)
[1305]956
[1442]957         endif
[1305]958
[1442]959        DO k = 1, klev
960         DO i = 1, klon
961
[1305]962c--------------------
963c   Ce test est necessaire pour eviter Xliq=NaN   
964        IF (ieee_is_nan(d_tr_sed(i,k,1)).OR.
965     &  ieee_is_nan(d_tr_sed(i,k,2))) THEN
966        PRINT*,'sedim NaN PROBLEM'
967        PRINT*,'d_tr_sed Nan?',d_tr_sed(i,k,:),'Temp',t_seri(i,k)
968        PRINT*,'lat-lon',i,'level',k,'dtime',dtime
[1442]969        PRINT*,'F_sed',Fsedim(i,k)
[1305]970        PRINT*,'==============================================='
971                d_tr_sed(i,k,:)=0.
972        ENDIF
973c--------------------
974
975        tr_seri(i,k,i_h2so4liq) = tr_seri(i,k,i_h2so4liq)+
976     &                            d_tr_sed(i,k,1)
977        tr_seri(i,k,i_h2oliq)   = tr_seri(i,k,i_h2oliq)+
978     &                            d_tr_sed(i,k,2)
979        d_tr_sed(i,k,:) = d_tr_sed(i,k,:) / dtime
980        Fsedim(i,k)     = Fsedim(i,k) / dtime
981     
982          ENDDO
[1442]983         ENDDO
[1305]984     
985        Fsedim(:,klev+1) = 0.
986
987         endif ! ok_sedim
988
989       endif   ! tr_scheme
990      endif    ! iflag_trac
991
[3]992c
993c====================================================================
994c Appeler la diffusion verticale (programme de couche limite)
995c====================================================================
996
997c-------------------------------
998c VENUS TEST: on ne tient pas compte des calculs de clmain mais on force
999c l'equilibre radiatif du sol
1000      if (1.eq.0) then
1001              if (debut) then
1002                print*,"ATTENTION, CLMAIN SHUNTEE..."
1003              endif
1004
1005      DO i = 1, klon
1006         sens(i) = 0.0e0 ! flux de chaleur sensible au sol
1007         fder(i) = 0.0e0
1008         dlw(i)  = 0.0e0
1009      ENDDO
1010
1011c Incrementer la temperature du sol
1012c
1013      DO i = 1, klon
1014         d_ts(i)  = dtime * radsol(i)/22000. !valeur calculee par GCM pour I=200
1015         ftsol(i) = ftsol(i) + d_ts(i)
1016         do j=1,nsoilmx
1017           ftsoil(i,j)=ftsol(i)
1018         enddo
1019      ENDDO
1020
1021c-------------------------------
1022      else
1023c-------------------------------
1024
1025      fder = dlw
1026
1027! ADAPTATION GCM POUR CP(T)
[808]1028
[3]1029      CALL clmain(dtime,itap,
1030     e            t_seri,u_seri,v_seri,
1031     e            rmu0,
1032     e            ftsol,
1033     $            ftsoil,
1034     $            paprs,pplay,ppk,radsol,falbe,
1035     e            solsw, sollw, sollwdown, fder,
1036     e            rlond, rlatd, cuphy, cvphy,   
1037     e            debut, lafin,
1038     s            d_t_vdf,d_u_vdf,d_v_vdf,d_ts,
1039     s            fluxt,fluxu,fluxv,cdragh,cdragm,
1040     s            dsens,
1041     s            ycoefh,yu1,yv1)
1042
1043CXXX Incrementation des flux
1044      DO i = 1, klon
1045         sens(i) = - fluxt(i,1) ! flux de chaleur sensible au sol
1046         fder(i) = dlw(i) + dsens(i)
1047      ENDDO
1048CXXX
1049
1050      DO k = 1, klev
1051      DO i = 1, klon
1052         t_seri(i,k) = t_seri(i,k) + d_t_vdf(i,k)
1053         d_t_vdf(i,k)= d_t_vdf(i,k)/dtime          ! K/s
1054         u_seri(i,k) = u_seri(i,k) + d_u_vdf(i,k)
1055         d_u_vdf(i,k)= d_u_vdf(i,k)/dtime          ! (m/s)/s
1056         v_seri(i,k) = v_seri(i,k) + d_v_vdf(i,k)
1057         d_v_vdf(i,k)= d_v_vdf(i,k)/dtime          ! (m/s)/s
1058      ENDDO
1059      ENDDO
1060
1061C TRACEURS
1062
1063      if (iflag_trac.eq.1) then
[101]1064         DO k = 1, klev
1065         DO i = 1, klon
1066            delp(i,k) = paprs(i,k)-paprs(i,k+1)
1067         ENDDO
1068         ENDDO
[1442]1069   
[3]1070         DO iq=1, nqmax
[1442]1071     
[3]1072             CALL cltrac(dtime,ycoefh,t_seri,
[1442]1073     s               tr_seri(:,:,iq),source(:,iq),
[3]1074     e               paprs, pplay,delp,
[1442]1075     s               d_tr_vdf(:,:,iq))
1076     
[3]1077             tr_seri(:,:,iq) = tr_seri(:,:,iq) + d_tr_vdf(:,:,iq)
1078             d_tr_vdf(:,:,iq)= d_tr_vdf(:,:,iq)/dtime          ! /s
[1442]1079
1080         ENDDO !nqmax
[3]1081
[1442]1082       endif
1083
[3]1084      IF (if_ebil.ge.2) THEN
1085        ztit='after clmain'
1086        CALL diagetpq(airephy,ztit,ip_ebil,2,1,dtime
1087     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
1088     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1089         call diagphy(airephy,ztit,ip_ebil
1090     e      , zero_v, zero_v, zero_v, zero_v, sens
1091     e      , zero_v, zero_v, zero_v, ztsol
1092     e      , d_h_vcol, d_qt, d_ec
1093     s      , fs_bound, fq_bound )
1094      END IF
1095C
1096c
1097c Incrementer la temperature du sol
1098c
1099      DO i = 1, klon
1100         ftsol(i) = ftsol(i) + d_ts(i)
1101      ENDDO
1102
1103c Calculer la derive du flux infrarouge
1104c
1105      DO i = 1, klon
1106            dlw(i) = - 4.0*RSIGMA*ftsol(i)**3
1107      ENDDO
1108
1109c-------------------------------
1110      endif  ! fin du VENUS TEST
1111
[973]1112      ! tests: output tendencies
1113!      call writefield_phy('physiq_d_t_vdf',d_t_vdf,klev)
1114!      call writefield_phy('physiq_d_u_vdf',d_u_vdf,klev)
1115!      call writefield_phy('physiq_d_v_vdf',d_v_vdf,klev)
1116!      call writefield_phy('physiq_d_ts',d_ts,1)
1117
[3]1118c
1119c Appeler l'ajustement sec
1120c
1121c===================================================================
1122c Convection seche
1123c===================================================================
1124c
1125      d_t_ajs(:,:)=0.
1126      d_u_ajs(:,:)=0.
1127      d_v_ajs(:,:)=0.
1128      d_tr_ajs(:,:,:)=0.
1129c
1130      IF(prt_level>9)WRITE(lunout,*)
1131     .    'AVANT LA CONVECTION SECHE , iflag_ajs='
1132     s   ,iflag_ajs
1133
1134      if(iflag_ajs.eq.0) then
1135c  Rien
1136c  ====
1137         IF(prt_level>9)WRITE(lunout,*)'pas de convection'
1138
1139      else if(iflag_ajs.eq.1) then
1140
1141c  Ajustement sec
1142c  ==============
1143         IF(prt_level>9)WRITE(lunout,*)'ajsec'
1144
1145! ADAPTATION GCM POUR CP(T)
1146         CALL ajsec(paprs, pplay, ppk, t_seri, u_seri, v_seri, nqmax,
1147     .              tr_seri, d_t_ajs, d_u_ajs, d_v_ajs, d_tr_ajs)
1148
1149! ADAPTATION GCM POUR CP(T)
1150         do i=1,klon
1151          flux_ajs(i,1) = 0.0
1152          do j=2,klev
1153            flux_ajs(i,j) = flux_ajs(i,j-1)
1154     .        + cpdet(t_seri(i,j-1))/RG*d_t_ajs(i,j-1)/dtime
1155     .                                 *(paprs(i,j-1)-paprs(i,j))
1156          enddo
1157         enddo
1158         
1159         t_seri(:,:) = t_seri(:,:) + d_t_ajs(:,:)
1160         d_t_ajs(:,:)= d_t_ajs(:,:)/dtime          ! K/s
1161         u_seri(:,:) = u_seri(:,:) + d_u_ajs(:,:)
1162         d_u_ajs(:,:)= d_u_ajs(:,:)/dtime          ! (m/s)/s
1163         v_seri(:,:) = v_seri(:,:) + d_v_ajs(:,:)
1164         d_v_ajs(:,:)= d_v_ajs(:,:)/dtime          ! (m/s)/s
[1301]1165
1166         if (iflag_trac.eq.1) then
[3]1167           tr_seri(:,:,:) = tr_seri(:,:,:) + d_tr_ajs(:,:,:)
1168           d_tr_ajs(:,:,:)= d_tr_ajs(:,:,:)/dtime  ! /s
[1301]1169         endif
[3]1170      endif
[973]1171
1172      ! tests: output tendencies
1173!      call writefield_phy('physiq_d_t_ajs',d_t_ajs,klev)
1174!      call writefield_phy('physiq_d_u_ajs',d_u_ajs,klev)
1175!      call writefield_phy('physiq_d_v_ajs',d_v_ajs,klev)
[3]1176c
1177      IF (if_ebil.ge.2) THEN
1178        ztit='after dry_adjust'
1179        CALL diagetpq(airephy,ztit,ip_ebil,2,2,dtime
1180     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
1181     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1182        call diagphy(airephy,ztit,ip_ebil
1183     e      , zero_v, zero_v, zero_v, zero_v, sens
1184     e      , zero_v, zero_v, zero_v, ztsol
1185     e      , d_h_vcol, d_qt, d_ec
1186     s      , fs_bound, fq_bound )
1187      END IF
1188
[1442]1189
[3]1190c====================================================================
1191c RAYONNEMENT
1192c====================================================================
1193
[1442]1194c------------------------------------
[1310]1195c    . Compute radiative tendencies :
1196c------------------------------------
1197c====================================================================
[3]1198      IF (MOD(itaprad,radpas).EQ.0) THEN
[1310]1199c====================================================================
[3]1200
[808]1201      dtimerad = dtime*REAL(radpas)  ! pas de temps du rayonnement (s)
[3]1202c     PRINT*,'dtimerad,dtime,radpas',dtimerad,dtime,radpas
1203           
[1310]1204
[1442]1205c------------------------------------
1206c    . Compute mean mass, cp and R :
1207c------------------------------------
1208
[1310]1209      if(callthermos) then
[1452]1210         call concentrations2(pplay,t_seri,d_t,tr_seri, nqmax)
[1442]1211
[1310]1212      endif
1213
1214
[1442]1215cc!!! ADD key callhedin
[1310]1216
[1442]1217      IF(callnlte.or.callthermos) THEN                                 
1218         call compo_hedin83_mod(pplay,rmu0,   
1219     &           co2vmr_gcm,covmr_gcm,ovmr_gcm,n2vmr_gcm,nvmr_gcm)
1220
1221         IF(ok_chem) then
1222 
1223CC  !! GG : Using only mayor species tracers abundances to compute NLTE heating/cooling
1224
1225CC               Conversion [mmr] ---> [vmr]
1226       
1227                 co2vmr_gcm(:,:) = tr_seri(1:nlon,1:nlev,i_co2)*
1228     &                             mmean(1:nlon,1:nlev)/M_tr(i_co2)
1229                 covmr_gcm(:,:)  = tr_seri(1:nlon,1:nlev,i_co)*
1230     &                              mmean(1:nlon,1:nlev)/M_tr(i_co)
1231                 ovmr_gcm(:,:)   = tr_seri(1:nlon,1:nlev,i_o)*
1232     &                             mmean(1:nlon,1:nlev)/M_tr(i_o)
1233                 n2vmr_gcm(:,:)   = tr_seri(1:nlon,1:nlev,i_n2)*
1234     &                             mmean(1:nlon,1:nlev)/M_tr(i_n2)
1235
1236         ENDIF
1237
1238       ENDIF   
1239
1240c
1241c   NLTE cooling from CO2 emission
1242c   ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1243
[1310]1244        IF(callnlte) THEN                                 
1245            if(nltemodel.eq.0.or.nltemodel.eq.1) then
[1442]1246                CALL nltecool(klon, klev, nqmax, pplay*9.869e-6, t_seri,
1247     $                    tr_seri, d_t_nlte)
[1310]1248            else if(nltemodel.eq.2) then                               
[1442]1249               CALL nlte_tcool(klon,klev,pplay*9.869e-6,             
[1310]1250     $               t_seri,zzlay,co2vmr_gcm, n2vmr_gcm, covmr_gcm,
1251     $               ovmr_gcm,d_t_nlte,ierr_nlte,varerr )
1252                  if(ierr_nlte.gt.0) then
1253                     write(*,*)
1254     $                'WARNING: nlte_tcool output with error message',
1255     $                'ierr_nlte=',ierr_nlte,'varerr=',varerr
1256                     write(*,*)'I will continue anyway'
1257                  endif
1258
1259             endif
1260             
1261        ELSE
1262 
1263          d_t_nlte(:,:)=0.
1264
1265        ENDIF   
1266
1267c      Find number of layers for LTE radiation calculations
1268
1269      IF(callnlte .or. callnirco2)
1270     $        CALL nlthermeq(klon, klev, paprs, pplay)
1271
[1442]1272c
1273c       LTE radiative transfert / solar / IR matrix
1274c       ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[3]1275      CALL radlwsw
[1301]1276     e            (dist, rmu0, fract, zzlev,
[1310]1277     e             paprs, pplay,ftsol, t_seri)
[808]1278
[1301]1279
[1442]1280c       CO2 near infrared absorption
1281c      ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[1310]1282
1283        d_t_nirco2(:,:)=0.
1284        if (callnirco2) then
[1442]1285           call nirco2abs (klon, klev, pplay, dist, nqmax, tr_seri,
1286     .                 rmu0, fract, d_t_nirco2)
[1310]1287        endif
1288
1289
1290c          Net atmospheric radiative heating rate (K.s-1)
1291c          ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1292
1293        IF(callnlte.or.callnirco2) THEN
1294           CALL blendrad(klon, klev, pplay,heat,
1295     &          cool, d_t_nirco2,d_t_nlte, dtsw, dtlw)
1296        ELSE
1297           dtsw(:,:)=heat(:,:)
1298           dtlw(:,:)=-1*cool(:,:)
1299        ENDIF
1300
1301         DO k=1,klev
1302            DO i=1,klon
1303               d_t_rad(i,k) = dtsw(i,k) + dtlw(i,k)   ! K/s
1304            ENDDO
1305         ENDDO
1306
1307
1308cc---------------------------------------------
1309
1310c          EUV heating rate (K.s-1)
1311c          ~~~~~~~~~~~~~~~~~~~~~~~~
1312
1313        d_t_euv(:,:)=0.
1314
1315        IF (callthermos) THEN
1316
[1442]1317c           call euvheat(klon, klev,t_seri,paprs,pplay,zzlay,
1318c     $          rmu0,pdtphys,gmtime,rjourvrai, co2vmr_gcm, n2vmr_gcm,
1319c     $          covmr_gcm, ovmr_gcm,d_t_euv )
1320           call euvheat(klon, klev, nqmax, t_seri,paprs,pplay,zzlay,
1321     $         rmu0,pdtphys,gmtime,rjourvrai,
1322     $         tr_seri, d_tr, d_t_euv )
1323               
[1310]1324           DO k=1,klev
1325              DO ig=1,klon
1326                 d_t_rad(ig,k)=d_t_rad(ig,k)+d_t_euv(ig,k)
1327               
1328              ENDDO
1329           ENDDO
1330
1331        ENDIF  ! callthermos
1332
[3]1333c====================================================================
[1310]1334        itaprad = 0
1335      ENDIF    ! radpas
1336c====================================================================
[3]1337c
1338c Ajouter la tendance des rayonnements (tous les pas)
1339c
1340      DO k = 1, klev
1341      DO i = 1, klon
[1310]1342         t_seri(i,k) = t_seri(i,k) + d_t_rad(i,k) * dtime
[3]1343      ENDDO
1344      ENDDO
[973]1345
[1310]1346! CONDUCTION  and  MOLECULAR VISCOSITY
1347c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1348
1349        d_t_conduc(:,:)=0.
1350        d_u_molvis(:,:)=0.
1351        d_v_molvis(:,:)=0.
1352
1353        IF (callthermos) THEN
1354
1355           tsurf(:)=t_seri(:,1)
1356           call conduction(klon, klev,pdtphys,
1357     $            pplay,paprs,t_seri,
1358     $            tsurf,zzlev,zzlay,d_t_conduc)
1359
1360            call molvis(klon, klev,pdtphys,
1361     $            pplay,paprs,t_seri,
1362     $            u,tsurf,zzlev,zzlay,d_u_molvis)
1363
1364            call molvis(klon, klev, pdtphys,
1365     $            pplay,paprs,t_seri,
[1442]1366     $            v,tsurf,zzlev,zzlay,d_v_molvis)
[1310]1367
1368            DO k=1,klev
1369               DO ig=1,klon
1370                  t_seri(ig,k)= t_seri(ig,k)+ d_t_conduc(ig,k)*dtime ! [K]
1371                  u_seri(ig,k)= u_seri(ig,k)+ d_u_molvis(ig,k)*dtime ! m/s
1372                  v_seri(ig,k)= v_seri(ig,k)+ d_v_molvis(ig,k)*dtime ! m/s
1373               ENDDO
1374            ENDDO
[1442]1375        ENDIF
[1310]1376
1377
[1442]1378!  --  MOLECULAR DIFFUSION ---
1379
1380          d_q_moldif(:,:,:)=0
1381
1382         IF (callthermos .and. ok_chem) THEN
1383
1384             call moldiff_red(klon, klev, nqmax,
1385     &                   pplay,paprs,t_seri, tr_seri, pdtphys,
1386     &                   zzlay,d_t_euv,d_t_conduc,d_q_moldif)
1387
1388
1389! --- update tendencies tracers ---
1390
1391          DO iq = 1, nqmax
1392           DO k=1,klev
1393              DO ig=1,klon
1394                tr_seri(ig,k,iq)= tr_seri(ig,k,iq)+
1395     &                           d_q_moldif(ig,k,iq)*dtime ! [Kg/kg]?
1396              ENDDO
1397            ENDDO
1398           ENDDO
1399           
1400
1401         ENDIF  ! callthermos & ok_chem
1402
[1310]1403c====================================================================
[973]1404      ! tests: output tendencies
1405!      call writefield_phy('physiq_dtrad',dtrad,klev)
[3]1406 
1407      IF (if_ebil.ge.2) THEN
1408        ztit='after rad'
1409        CALL diagetpq(airephy,ztit,ip_ebil,2,2,dtime
1410     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
1411     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1412        call diagphy(airephy,ztit,ip_ebil
1413     e      , topsw, toplw, solsw, sollw, zero_v
1414     e      , zero_v, zero_v, zero_v, ztsol
1415     e      , d_h_vcol, d_qt, d_ec
1416     s      , fs_bound, fq_bound )
1417      END IF
1418c
1419
1420c====================================================================
1421c   Calcul  des gravity waves  FLOTT
1422c====================================================================
1423c
1424      if (ok_orodr.or.ok_gw_nonoro) then
1425c  CALCUL DE N2
1426       do i=1,klon
1427        do k=2,klev
1428          ztlev(i,k)  = (t_seri(i,k)+t_seri(i,k-1))/2.
1429          zpklev(i,k) = sqrt(ppk(i,k)*ppk(i,k-1))
1430        enddo
1431       enddo
1432       call t2tpot(klon*klev,ztlev, ztetalev,zpklev)
1433       call t2tpot(klon*klev,t_seri,ztetalay,ppk)
1434       do i=1,klon
1435        do k=2,klev
1436          zdtetalev(i,k) = ztetalay(i,k)-ztetalay(i,k-1)
1437          zdzlev(i,k)    = (zphi(i,k)-zphi(i,k-1))/RG
1438          zn2(i,k) = RG*zdtetalev(i,k)/(ztetalev(i,k)*zdzlev(i,k))
1439          zn2(i,k) = max(zn2(i,k),1.e-12)  ! securite
1440        enddo
[808]1441        zn2(i,1) = 1.e-12  ! securite
[3]1442       enddo
1443
1444      endif
1445     
1446c ----------------------------ORODRAG
1447      IF (ok_orodr) THEN
1448c
1449c  selection des points pour lesquels le shema est actif:
1450        igwd=0
1451        DO i=1,klon
1452        itest(i)=0
1453c        IF ((zstd(i).gt.10.0)) THEN
1454        IF (((zpic(i)-zmea(i)).GT.100.).AND.(zstd(i).GT.10.0)) THEN
1455          itest(i)=1
1456          igwd=igwd+1
1457          idx(igwd)=i
1458        ENDIF
1459        ENDDO
1460c        igwdim=MAX(1,igwd)
1461c
1462c A ADAPTER POUR VENUS!!!
[1306]1463        CALL drag_noro(klon,klev,dtime,paprs,pplay,pphi,zn2,
[3]1464     e                   zmea,zstd, zsig, zgam, zthe,zpic,zval,
1465     e                   igwd,idx,itest,
1466     e                   t_seri, u_seri, v_seri,
1467     s                   zulow, zvlow, zustrdr, zvstrdr,
1468     s                   d_t_oro, d_u_oro, d_v_oro)
1469
[808]1470c       print*,"d_u_oro=",d_u_oro(klon/2,:)
[3]1471c  ajout des tendances
1472           t_seri(:,:) = t_seri(:,:) + d_t_oro(:,:)
1473           d_t_oro(:,:)= d_t_oro(:,:)/dtime          ! K/s
1474           u_seri(:,:) = u_seri(:,:) + d_u_oro(:,:)
1475           d_u_oro(:,:)= d_u_oro(:,:)/dtime          ! (m/s)/s
1476           v_seri(:,:) = v_seri(:,:) + d_v_oro(:,:)
1477           d_v_oro(:,:)= d_v_oro(:,:)/dtime          ! (m/s)/s
1478c   
1479      ELSE
1480         d_t_oro = 0.
1481         d_u_oro = 0.
1482         d_v_oro = 0.
1483         zustrdr = 0.
1484         zvstrdr = 0.
1485c
1486      ENDIF ! fin de test sur ok_orodr
1487c
[973]1488      ! tests: output tendencies
1489!      call writefield_phy('physiq_d_t_oro',d_t_oro,klev)
1490!      call writefield_phy('physiq_d_u_oro',d_u_oro,klev)
1491!      call writefield_phy('physiq_d_v_oro',d_v_oro,klev)
1492
[3]1493c ----------------------------OROLIFT
1494      IF (ok_orolf) THEN
[808]1495       print*,"ok_orolf NOT IMPLEMENTED !"
1496       stop
[3]1497c
1498c  selection des points pour lesquels le shema est actif:
1499        igwd=0
1500        DO i=1,klon
1501        itest(i)=0
1502        IF ((zpic(i)-zmea(i)).GT.100.) THEN
1503          itest(i)=1
1504          igwd=igwd+1
1505          idx(igwd)=i
1506        ENDIF
1507        ENDDO
1508c        igwdim=MAX(1,igwd)
1509c
1510c A ADAPTER POUR VENUS!!!
1511c            CALL lift_noro(klon,klev,dtime,paprs,pplay,
1512c     e                   rlatd,zmea,zstd,zpic,zgam,zthe,zpic,zval,
1513c     e                   igwd,idx,itest,
1514c     e                   t_seri, u_seri, v_seri,
1515c     s                   zulow, zvlow, zustrli, zvstrli,
1516c     s                   d_t_lif, d_u_lif, d_v_lif               )
1517
1518c
1519c  ajout des tendances
1520           t_seri(:,:) = t_seri(:,:) + d_t_lif(:,:)
1521           d_t_lif(:,:)= d_t_lif(:,:)/dtime          ! K/s
1522           u_seri(:,:) = u_seri(:,:) + d_u_lif(:,:)
1523           d_u_lif(:,:)= d_u_lif(:,:)/dtime          ! (m/s)/s
1524           v_seri(:,:) = v_seri(:,:) + d_v_lif(:,:)
1525           d_v_lif(:,:)= d_v_lif(:,:)/dtime          ! (m/s)/s
1526c
1527      ELSE
1528         d_t_lif = 0.
1529         d_u_lif = 0.
1530         d_v_lif = 0.
1531         zustrli = 0.
1532         zvstrli = 0.
1533c
1534      ENDIF ! fin de test sur ok_orolf
1535
1536c ---------------------------- NON-ORO GRAVITY WAVES
1537       IF(ok_gw_nonoro) then
1538
1539      call flott_gwd_ran(klon,klev,dtime,pplay,zn2,
1540     e               t_seri, u_seri, v_seri,
1541     o               zustrhi,zvstrhi,
1542     o               d_t_hin, d_u_hin, d_v_hin)
1543
1544c  ajout des tendances
1545
1546         t_seri(:,:) = t_seri(:,:) + d_t_hin(:,:)
1547         d_t_hin(:,:)= d_t_hin(:,:)/dtime          ! K/s
1548         u_seri(:,:) = u_seri(:,:) + d_u_hin(:,:)
1549         d_u_hin(:,:)= d_u_hin(:,:)/dtime          ! (m/s)/s
1550         v_seri(:,:) = v_seri(:,:) + d_v_hin(:,:)
1551         d_v_hin(:,:)= d_v_hin(:,:)/dtime          ! (m/s)/s
1552
1553      ELSE
1554         d_t_hin = 0.
1555         d_u_hin = 0.
1556         d_v_hin = 0.
1557         zustrhi = 0.
1558         zvstrhi = 0.
1559
1560      ENDIF ! fin de test sur ok_gw_nonoro
1561
[973]1562      ! tests: output tendencies
1563!      call writefield_phy('physiq_d_t_hin',d_t_hin,klev)
1564!      call writefield_phy('physiq_d_u_hin',d_u_hin,klev)
1565!      call writefield_phy('physiq_d_v_hin',d_v_hin,klev)
1566
[3]1567c====================================================================
1568c Transport de ballons
1569c====================================================================
1570      if (ballons.eq.1) then
[953]1571         CALL ballon(30,pdtphys,rjourvrai,gmtime*RDAY,rlatd,rlond,
[3]1572c    C               t,pplay,u,v,pphi)   ! alt above surface (smoothed for GCM)
1573     C               t,pplay,u,v,zphi)   ! alt above planet average radius
1574      endif !ballons
1575
1576c====================================================================
1577c Bilan de mmt angulaire
1578c====================================================================
1579      if (bilansmc.eq.1) then
1580CMODDEB FLOTT
1581C  CALCULER LE BILAN DE MOMENT ANGULAIRE (DIAGNOSTIQUE)
1582C STRESS NECESSAIRES: COUCHE LIMITE ET TOUTE LA PHYSIQUE
1583
1584      DO i = 1, klon
1585        zustrph(i)=0.
1586        zvstrph(i)=0.
1587        zustrcl(i)=0.
1588        zvstrcl(i)=0.
1589      ENDDO
1590      DO k = 1, klev
1591      DO i = 1, klon
1592       zustrph(i)=zustrph(i)+(u_seri(i,k)-u(i,k))/dtime*
1593     c            (paprs(i,k)-paprs(i,k+1))/rg
1594       zvstrph(i)=zvstrph(i)+(v_seri(i,k)-v(i,k))/dtime*
1595     c            (paprs(i,k)-paprs(i,k+1))/rg
1596       zustrcl(i)=zustrcl(i)+d_u_vdf(i,k)*
1597     c            (paprs(i,k)-paprs(i,k+1))/rg
1598       zvstrcl(i)=zvstrcl(i)+d_v_vdf(i,k)*
1599     c            (paprs(i,k)-paprs(i,k+1))/rg
1600      ENDDO
1601      ENDDO
1602
[953]1603      CALL aaam_bud (27,klon,klev,rjourvrai,gmtime*RDAY,
[3]1604     C               ra,rg,romega,
1605     C               rlatd,rlond,pphis,
1606     C               zustrdr,zustrli,zustrcl,
1607     C               zvstrdr,zvstrli,zvstrcl,
1608     C               paprs,u,v)
1609                     
1610CCMODFIN FLOTT
1611      endif !bilansmc
1612
1613c====================================================================
1614c====================================================================
1615c Calculer le transport de l'eau et de l'energie (diagnostique)
1616c
1617c  A REVOIR POUR VENUS...
1618c
1619c     CALL transp (paprs,ftsol,
1620c    e                   t_seri, q_seri, u_seri, v_seri, zphi,
1621c    s                   ve, vq, ue, uq)
1622c
1623c====================================================================
1624c+jld ec_conser
1625      DO k = 1, klev
1626      DO i = 1, klon
1627        d_t_ec(i,k)=0.5/cpdet(t_seri(i,k))
1628     $      *(u(i,k)**2+v(i,k)**2-u_seri(i,k)**2-v_seri(i,k)**2)
1629        t_seri(i,k)=t_seri(i,k)+d_t_ec(i,k)
1630        d_t_ec(i,k) = d_t_ec(i,k)/dtime
1631       END DO
1632      END DO
1633         do i=1,klon
1634          flux_ec(i,1) = 0.0
1635          do j=2,klev
1636            flux_ec(i,j) = flux_ec(i,j-1)
1637     . +cpdet(t_seri(i,j-1))/RG*d_t_ec(i,j-1)*(paprs(i,j-1)-paprs(i,j))
1638          enddo
1639         enddo
1640         
1641c-jld ec_conser
1642c====================================================================
1643      IF (if_ebil.ge.1) THEN
1644        ztit='after physic'
1645        CALL diagetpq(airephy,ztit,ip_ebil,1,1,dtime
1646     e      , t_seri,zero_v2,zero_v2,zero_v2,u_seri,v_seri,paprs,pplay
1647     s      , d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1648C     Comme les tendances de la physique sont ajoute dans la dynamique,
1649C     on devrait avoir que la variation d'entalpie par la dynamique
1650C     est egale a la variation de la physique au pas de temps precedent.
1651C     Donc la somme de ces 2 variations devrait etre nulle.
1652        call diagphy(airephy,ztit,ip_ebil
1653     e      , topsw, toplw, solsw, sollw, sens
1654     e      , zero_v, zero_v, zero_v, ztsol
1655     e      , d_h_vcol, d_qt, d_ec
1656     s      , fs_bound, fq_bound )
1657C
1658      d_h_vcol_phy=d_h_vcol
1659C
1660      END IF
1661C
1662c=======================================================================
1663c   SORTIES
1664c=======================================================================
1665
1666c Convertir les incrementations en tendances
1667c
1668      DO k = 1, klev
1669      DO i = 1, klon
1670         d_u(i,k) = ( u_seri(i,k) - u(i,k) ) / dtime
1671         d_v(i,k) = ( v_seri(i,k) - v(i,k) ) / dtime
1672         d_t(i,k) = ( t_seri(i,k) - t(i,k) ) / dtime
1673      ENDDO
1674      ENDDO
1675c
1676      DO iq = 1, nqmax
1677      DO  k = 1, klev
1678      DO  i = 1, klon
1679         d_qx(i,k,iq) = ( tr_seri(i,k,iq) - qx(i,k,iq) ) / dtime
1680      ENDDO
1681      ENDDO
1682      ENDDO
1683     
1684c------------------------
1685c Calcul moment cinetique
1686c------------------------
1687c TEST VENUS...
1688c     mangtot = 0.0
1689c     DO k = 1, klev
1690c     DO i = 1, klon
1691c       mang(i,k) = RA*cos(rlatd(i)*RPI/180.)
1692c    .     *(u_seri(i,k)+RA*cos(rlatd(i)*RPI/180.)*ROMEGA)
1693c    .     *airephy(i)*(paprs(i,k)-paprs(i,k+1))/RG
1694c       mangtot=mangtot+mang(i,k)
1695c     ENDDO
1696c     ENDDO
1697c     print*,"Moment cinetique total = ",mangtot
1698c
1699c------------------------
1700c
1701c Sauvegarder les valeurs de t et u a la fin de la physique:
1702c
1703      DO k = 1, klev
1704      DO i = 1, klon
1705         u_ancien(i,k) = u_seri(i,k)
1706         t_ancien(i,k) = t_seri(i,k)
1707      ENDDO
1708      ENDDO
1709c
1710c=============================================================
1711c   Ecriture des sorties
1712c=============================================================
1713     
1714#ifdef CPP_IOIPSL
1715
1716#ifdef histhf
1717#include "write_histhf.h"
1718#endif
1719
1720#ifdef histday
1721#include "write_histday.h"
1722#endif
1723
1724#ifdef histmth
1725#include "write_histmth.h"
1726#endif
1727
1728#ifdef histins
1729#include "write_histins.h"
1730#endif
1731
1732#endif
1733
1734c====================================================================
1735c Si c'est la fin, il faut conserver l'etat de redemarrage
1736c====================================================================
1737c
1738      IF (lafin) THEN
1739         itau_phy = itau_phy + itap
[892]1740         CALL phyredem ("restartphy.nc")
[3]1741     
1742c--------------FLOTT
1743CMODEB LOTT
1744C  FERMETURE DU FICHIER FORMATTE CONTENANT LES COMPOSANTES
1745C  DU BILAN DE MOMENT ANGULAIRE.
1746      if (bilansmc.eq.1) then
1747        write(*,*)'Fermeture de aaam_bud.out (FL Vous parle)'
1748        close(27)                                     
1749        close(28)                                     
1750      endif !bilansmc
1751CMODFIN
1752c-------------
1753c--------------SLEBONNOIS
1754C  FERMETURE DES FICHIERS FORMATTES CONTENANT LES POSITIONS ET VITESSES
1755C  DES BALLONS
1756      if (ballons.eq.1) then
1757        write(*,*)'Fermeture des ballons*.out'
1758        close(30)                                     
1759        close(31)                                     
1760        close(32)                                     
1761        close(33)                                     
1762        close(34)                                     
1763      endif !ballons
1764c-------------
1765      ENDIF
1766     
1767      RETURN
1768      END
1769
1770
1771
1772***********************************************************************
1773***********************************************************************
1774***********************************************************************
1775***********************************************************************
1776***********************************************************************
1777***********************************************************************
1778***********************************************************************
1779***********************************************************************
1780
1781      SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit)
1782      IMPLICIT none
1783c
1784c Tranformer une variable de la grille physique a
1785c la grille d'ecriture
1786c
1787      INTEGER nfield,nlon,iim,jjmp1, jjm
1788      REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield)
1789c
1790      INTEGER i, n, ig
1791c
1792      jjm = jjmp1 - 1
1793      DO n = 1, nfield
1794         DO i=1,iim
1795            ecrit(i,n) = fi(1,n)
1796            ecrit(i+jjm*iim,n) = fi(nlon,n)
1797         ENDDO
1798         DO ig = 1, nlon - 2
1799           ecrit(iim+ig,n) = fi(1+ig,n)
1800         ENDDO
1801      ENDDO
1802      RETURN
1803      END
Note: See TracBrowser for help on using the repository browser.