source: LMDZ6/trunk/libf/phylmd/calltherm.F90 @ 4627

Last change on this file since 4627 was 4590, checked in by fhourdin, 17 months ago

Passage des thermiques a la nouvelle norme.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 18.1 KB
Line 
1!
2! $Id: calltherm.F90 4590 2023-06-29 01:03:15Z acozic $
3!
4      subroutine calltherm(dtime  &
5     &      ,pplay,paprs,pphi,weak_inversion  &
6     &      ,u_seri,v_seri,t_seri,q_seri,zqsat,debut  &
7     &      ,d_u_ajs,d_v_ajs,d_t_ajs,d_q_ajs  &
8     &      ,fm_therm,entr_therm,detr_therm,zqasc,clwcon0,lmax,ratqscth,  &
9     &       ratqsdiff,zqsatth,ale_bl,alp_bl,lalim_conv,wght_th, &
10     &       zmax0,f0,zw2,fraca,ztv,zpspsk,ztla,zthl &
11!!! nrlmd le 10/04/2012
12     &      ,pbl_tke,pctsrf,omega,airephy &
13     &      ,zlcl_th,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &
14     &      ,n2,s2,ale_bl_stat &
15     &      ,therm_tke_max,env_tke_max &
16     &      ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &
17     &      ,alp_bl_conv,alp_bl_stat &
18!!! fin nrlmd le 10/04/2012
19     &      ,zqla,ztva &
20#ifdef ISO         
21     &      ,xt_seri,d_xt_ajs &
22#ifdef DIAGISO         
23     &      ,q_the,xt_the &
24#endif
25#endif         
26     &   )
27
28      USE dimphy
29      USE indice_sol_mod
30      USE print_control_mod, ONLY: prt_level,lunout
31      USE lmdz_thermcell_alp, ONLY: thermcell_alp
32      USE lmdz_thermcell_main, ONLY: thermcell_main
33      USE lmdz_thermcell_old, ONLY: thermcell, thermcell_2002, thermcell_eau, calcul_sec, thermcell_sec
34#ifdef ISO
35      use infotrac_phy, ONLY: ntiso
36#ifdef ISOVERIF
37      USE isotopes_mod, ONLY: iso_eau,iso_HDO
38      USE isotopes_verif_mod, ONLY: iso_verif_aberrant_enc_vect2D, &
39        iso_verif_egalite_vect2D
40#endif   
41#endif
42
43      implicit none
44      include "clesphys.h"
45      include "thermcell_old.h"
46
47
48!IM 140508
49      INTEGER, SAVE ::  itap
50!$OMP THREADPRIVATE(itap)
51      REAL dtime
52      LOGICAL debut
53      LOGICAL logexpr0, logexpr2(klon,klev), logexpr1(klon)
54      REAL fact(klon)
55      INTEGER nbptspb
56
57      REAL u_seri(klon,klev),v_seri(klon,klev)
58      REAL t_seri(klon,klev),q_seri(klon,klev),qmemoire(klon,klev)
59      REAL weak_inversion(klon)
60      REAL paprs(klon,klev+1)
61      REAL pplay(klon,klev)
62      REAL pphi(klon,klev)
63      real zlev(klon,klev+1)
64!test: on sort lentr et a* pour alimenter KE
65      REAL wght_th(klon,klev)
66      INTEGER lalim_conv(klon)
67      REAL zw2(klon,klev+1),fraca(klon,klev+1)
68
69!FH Update Thermiques
70      REAL d_t_ajs(klon,klev), d_q_ajs(klon,klev)
71      REAL d_u_ajs(klon,klev),d_v_ajs(klon,klev)
72      real fm_therm(klon,klev+1)
73      real entr_therm(klon,klev),detr_therm(klon,klev)
74
75!********************************************************
76!     declarations
77      LOGICAL flag_bidouille_stratocu
78      real fmc_therm(klon,klev+1),zqasc(klon,klev)
79      real zqla(klon,klev)
80      real zqta(klon,klev)
81      real ztv(klon,klev),ztva(klon,klev)
82      real zpspsk(klon,klev)
83      real ztla(klon,klev)
84      real zthl(klon,klev)
85      real wmax_sec(klon)
86      real zmax_sec(klon)
87      real f_sec(klon)
88      real detrc_therm(klon,klev)
89! FH WARNING : il semble que ces save ne servent a rien
90!     save fmc_therm, detrc_therm
91      real clwcon0(klon,klev)
92      real zqsat(klon,klev)
93      real zw_sec(klon,klev+1)
94      integer lmix_sec(klon)
95      integer lmax(klon)
96      real ratqscth(klon,klev)
97      real ratqsdiff(klon,klev)
98      real zqsatth(klon,klev) 
99!nouvelles variables pour la convection
100      real ale_bl(klon)
101      real alp_bl(klon)
102      real ale(klon)
103      real alp(klon)
104!RC
105      !on garde le zmax du pas de temps precedent
106      real zmax0(klon), f0(klon)
107
108!!! nrlmd le 10/04/2012
109      real pbl_tke(klon,klev+1,nbsrf)
110      real pctsrf(klon,nbsrf)
111      real omega(klon,klev)
112      real airephy(klon)
113      real zlcl_th(klon),fraca0(klon),w0(klon),w_conv(klon)
114      real therm_tke_max0(klon),env_tke_max0(klon)
115      real n2(klon),s2(klon)
116      real ale_bl_stat(klon)
117      real therm_tke_max(klon,klev),env_tke_max(klon,klev)
118      real alp_bl_det(klon),alp_bl_fluct_m(klon),alp_bl_fluct_tke(klon),alp_bl_conv(klon),alp_bl_stat(klon)
119!!! fin nrlmd le 10/04/2012
120
121!********************************************************
122
123      real, dimension(klon) :: pcon
124      real, dimension(klon,klev) :: rhobarz,wth3
125      integer,dimension(klon) :: lalim
126      real, dimension(klon,klev+1) :: fm
127      real, dimension(klon,klev) :: alim_star
128      real, dimension(klon) :: zmax
129
130
131
132
133! variables locales
134      REAL d_t_the(klon,klev), d_q_the(klon,klev)
135      REAL d_u_the(klon,klev),d_v_the(klon,klev)
136!
137      real zfm_therm(klon,klev+1),zdt
138      real zentr_therm(klon,klev),zdetr_therm(klon,klev)
139! FH A VERIFIER : SAVE INUTILES
140!      save zentr_therm,zfm_therm
141
142      character (len=20) :: modname='calltherm'
143      character (len=80) :: abort_message
144
145      integer i,k,isplit
146      logical, save :: first=.true.
147      logical :: new_thermcell
148
149#ifdef ISO
150      REAL xt_seri(ntiso,klon,klev),xtmemoire(ntiso,klon,klev)
151      REAL d_xt_ajs(ntiso,klon,klev)
152      real d_xt_the(ntiso,klon,klev)
153#ifdef DIAGISO
154      real q_the(klon,klev)
155      real xt_the(ntiso,klon,klev)
156#endif
157      real qprec(klon,klev)
158      integer ixt
159#endif
160
161
162!$OMP THREADPRIVATE(first)
163!********************************************************
164      if (first) then
165        itap=0
166        first=.false.
167      endif
168
169! Incrementer le compteur de la physique
170     itap   = itap + 1
171
172!  Modele du thermique
173!  ===================
174!         print*,'thermiques: WARNING on passe t au lieu de t_seri'
175
176
177! On prend comme valeur initiale des thermiques la valeur du pas
178! de temps precedent
179         zfm_therm(:,:)=fm_therm(:,:)
180         zdetr_therm(:,:)=detr_therm(:,:)
181         zentr_therm(:,:)=entr_therm(:,:)
182
183! On reinitialise les flux de masse a zero pour le cumul en
184! cas de splitting
185         fm_therm(:,:)=0.
186         entr_therm(:,:)=0.
187         detr_therm(:,:)=0.
188
189         ale_bl(:)=0.
190         alp_bl(:)=0.
191         if (prt_level.ge.10) then
192          print*,'thermV4 nsplit: ',nsplit_thermals,' weak_inversion'
193         endif
194
195!   tests sur les valeurs negatives de l'eau
196         logexpr0=prt_level.ge.10
197         nbptspb=0
198         do k=1,klev
199            do i=1,klon
200! Attention teste abderr 19-03-09
201!               logexpr2(i,k)=.not.q_seri(i,k).ge.0.
202                logexpr2(i,k)=.not.q_seri(i,k).ge.1.e-15
203               if (logexpr2(i,k)) then
204#ifdef ISO
205                qprec(i,k)=q_seri(i,k)
206#endif
207                q_seri(i,k)=1.e-15
208                nbptspb=nbptspb+1
209#ifdef ISO
210                do ixt=1,ntiso
211                  xt_seri(ixt,i,k)=1.e-15*(xt_seri(ixt,i,k)/qprec(i,k))
212                  ! xt_seri(ixt,i,k)=1.e-15*(Rdefault(index_iso(ixt)))
213                enddo
214#endif
215               endif
216!               if (logexpr0) &
217!    &             print*,'WARN eau<0 avant therm i=',i,'  k=',k  &
218!    &         ,' dq,q',d_q_the(i,k),q_seri(i,k)
219            enddo
220         enddo
221         if(nbptspb.GT.0) print*,'Number of points with q_seri(i,k)<=0 ',nbptspb   
222
223
224         new_thermcell=iflag_thermals>=15.and.iflag_thermals<=18
225#ifdef ISO
226      if (.not.new_thermcell) then
227           CALL abort_gcm('calltherm 234','isos pas prevus ici',1)
228      endif
229#ifdef ISOVERIF
230      if (iso_eau.gt.0) then
231       call iso_verif_egalite_vect2D( &
232     &           xt_seri,q_seri, &
233     &           'calltherm 174',ntiso,klon,klev)
234      endif !if (iso_eau.gt.0) then
235#endif   
236#endif
237         zdt=dtime/REAL(nsplit_thermals)
238
239
240         do isplit=1,nsplit_thermals
241
242          if (iflag_thermals>=1000) then
243            CALL thermcell_2002(klon,klev,zdt,iflag_thermals   &
244     &      ,pplay,paprs,pphi  &
245     &      ,u_seri,v_seri,t_seri,q_seri  &
246     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
247     &      ,zfm_therm,zentr_therm,fraca,zw2  &
248     &      ,r_aspect_thermals,30.,w2di_thermals  &
249     &      ,tau_thermals)
250          else if (iflag_thermals.eq.2) then
251            CALL thermcell_sec(klon,klev,zdt  &
252     &      ,pplay,paprs,pphi,zlev  &
253     &      ,u_seri,v_seri,t_seri,q_seri  &
254     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
255     &      ,zfm_therm,zentr_therm  &
256     &      ,r_aspect_thermals,30.,w2di_thermals  &
257     &      ,tau_thermals)
258          else if (iflag_thermals.eq.3) then
259            CALL thermcell(klon,klev,zdt  &
260     &      ,pplay,paprs,pphi  &
261     &      ,u_seri,v_seri,t_seri,q_seri  &
262     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
263     &      ,zfm_therm,zentr_therm  &
264     &      ,r_aspect_thermals,l_mix_thermals,w2di_thermals  &
265     &      ,tau_thermals)
266          else if (iflag_thermals.eq.10) then
267            CALL thermcell_eau(klon,klev,zdt  &
268     &      ,pplay,paprs,pphi  &
269     &      ,u_seri,v_seri,t_seri,q_seri  &
270     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
271     &      ,zfm_therm,zentr_therm  &
272     &      ,r_aspect_thermals,l_mix_thermals,w2di_thermals  &
273     &      ,tau_thermals)
274          else if (iflag_thermals.eq.11) then
275              abort_message = 'cas non prevu dans calltherm'
276              CALL abort_physic (modname,abort_message,1)
277          else if (iflag_thermals.eq.12) then
278            CALL calcul_sec(klon,klev,zdt  &
279     &      ,pplay,paprs,pphi,zlev  &
280     &      ,u_seri,v_seri,t_seri,q_seri  &
281     &      ,zmax_sec,wmax_sec,zw_sec,lmix_sec  &
282     &      ,r_aspect_thermals,l_mix_thermals,w2di_thermals  &
283     &      ,tau_thermals)
284          else if (iflag_thermals==13.or.iflag_thermals==14) then
285              abort_message = 'thermcellV0_main enleve svn>2084'
286              CALL abort_physic (modname,abort_message,1)
287          else if (new_thermcell) then
288            CALL thermcell_main(itap,klon,klev,zdt  &
289     &      ,pplay,paprs,pphi,debut  &
290     &      ,u_seri,v_seri,t_seri,q_seri  &
291     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
292     &      ,zfm_therm,zentr_therm,zdetr_therm,zqasc,zqla,lmax  &
293     &      ,ratqscth,ratqsdiff,zqsatth  &
294     &      ,zmax0,f0,zw2,fraca,ztv,zpspsk &
295     &      ,ztla,zthl,ztva &
296     &      ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax &
297#ifdef ISO         
298     &      ,xt_seri,d_xt_the &
299#endif         
300     &   )
301
302            CALL thermcell_alp(klon,klev,zdt  &                      ! in
303     &        ,pplay,paprs  &                                        ! in
304     &        ,zfm_therm,zentr_therm,lmax  &                         ! in
305     &        ,pbl_tke,pctsrf,omega,airephy &                        ! in
306     &        ,zw2,fraca &                                           ! in
307     &        ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax &  ! in
308     &        ,ale,alp,lalim_conv,wght_th &                          ! out
309     &        ,zlcl_th,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &! out
310     &        ,n2,s2,ale_bl_stat &                                   ! out
311     &        ,therm_tke_max,env_tke_max &                           ! out
312     &        ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &          ! out
313     &        ,alp_bl_conv,alp_bl_stat &                             ! out
314     &        )
315
316           if (prt_level.gt.10) write(lunout,*)'Apres thermcell_main OK'
317         else
318           abort_message = 'Cas des thermiques non prevu'
319           CALL abort_physic (modname,abort_message,1)
320         endif
321
322! Attention : les noms sont contre intuitif.
323! flag_bidouille_stratocu est .true. si on ne fait pas de bidouille.
324! Il aurait mieux valu avoir un nobidouille_stratocu
325! Et pour simplifier :
326! nobidouille_stratocu=.not.(iflag_thermals==13.or.iflag_thermals=15)
327! Ce serait bien de changer, mai en prenant le temps de vérifier que ca
328! fait bien ce qu'on croit.
329
330       flag_bidouille_stratocu=iflag_thermals<=12.or.iflag_thermals==14.or.iflag_thermals==16.or.iflag_thermals==18
331
332! Calcul a posteriori du niveau max des thermiques pour les schémas qui
333! ne la sortent pas.
334      if (iflag_thermals<=12.or.iflag_thermals>=1000) then
335         lmax(:)=1
336         do k=1,klev-1
337            zdetr_therm(:,k)=zentr_therm(:,k)+zfm_therm(:,k)-zfm_therm(:,k+1)
338         enddo
339         do k=1,klev-1
340            do i=1,klon
341               if (zfm_therm(i,k+1)>0.) lmax(i)=k
342            enddo
343         enddo
344      endif
345
346      fact(:)=0.
347      DO i=1,klon
348       logexpr1(i)=flag_bidouille_stratocu.or.weak_inversion(i).gt.0.5
349       IF(logexpr1(i)) fact(i)=1./REAL(nsplit_thermals)
350      ENDDO
351
352     DO k=1,klev
353!  transformation de la derivee en tendance
354            d_t_the(:,k)=d_t_the(:,k)*dtime*fact(:)
355            d_u_the(:,k)=d_u_the(:,k)*dtime*fact(:)
356            d_v_the(:,k)=d_v_the(:,k)*dtime*fact(:)
357            d_q_the(:,k)=d_q_the(:,k)*dtime*fact(:)
358            fm_therm(:,k)=fm_therm(:,k)  &
359     &      +zfm_therm(:,k)*fact(:)
360            entr_therm(:,k)=entr_therm(:,k)  &
361     &       +zentr_therm(:,k)*fact(:)
362            detr_therm(:,k)=detr_therm(:,k)  &
363     &       +zdetr_therm(:,k)*fact(:)
364#ifdef ISO
365            do ixt=1,ntiso
366              d_xt_the(ixt,:,k)=d_xt_the(ixt,:,k)*dtime*fact(:)
367            enddo
368#endif
369      ENDDO
370       fm_therm(:,klev+1)=0.
371
372
373
374!  accumulation de la tendance
375            d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_the(:,:)
376            d_u_ajs(:,:)=d_u_ajs(:,:)+d_u_the(:,:)
377            d_v_ajs(:,:)=d_v_ajs(:,:)+d_v_the(:,:)
378            d_q_ajs(:,:)=d_q_ajs(:,:)+d_q_the(:,:)
379#ifdef ISO
380            d_xt_ajs(:,:,:)=d_xt_ajs(:,:,:)+d_xt_the(:,:,:)
381#endif
382
383!  incrementation des variables meteo
384            t_seri(:,:) = t_seri(:,:) + d_t_the(:,:)
385            u_seri(:,:) = u_seri(:,:) + d_u_the(:,:)
386            v_seri(:,:) = v_seri(:,:) + d_v_the(:,:)
387            qmemoire(:,:)=q_seri(:,:)
388            q_seri(:,:) = q_seri(:,:) + d_q_the(:,:)
389#ifdef ISO
390            xtmemoire(:,:,:)=xt_seri(:,:,:)
391            xt_seri(:,:,:) = xt_seri(:,:,:) + d_xt_the(:,:,:)
392#ifdef ISOVERIF
393!      write(*,*) 'calltherm 350 tmp: ajout d_xt_the'
394      if (iso_HDO.gt.0) then
395!      i=479
396!      k=4
397!      write(*,*) 'xt_seri(iso_hdo,i,k),q_seri(i,k)=', &
398!     &   xt_seri(iso_hdo,i,k),q_seri(i,k)
399!      write(*,*) 'd_xt_the(iso_hdo,i,k),d_q_the(i,k)=', &
400!     &   d_xt_the(iso_hdo,i,k),d_q_the(i,k)
401      call iso_verif_aberrant_enc_vect2D( &
402     &        xt_seri,q_seri, &
403     &        'calltherm 353, apres ajout d_xt_the',ntiso,klon,klev)
404      endif     
405#endif
406#endif
407           if (prt_level.gt.10) write(lunout,*)'Apres apres thermcell_main OK'
408
409       DO i=1,klon
410            fm_therm(i,klev+1)=0.
411            ale_bl(i)=ale_bl(i)+ale(i)/REAL(nsplit_thermals)
412!            write(22,*)'ALE CALLTHERM',ale_bl(i),ale(i)
413            alp_bl(i)=alp_bl(i)+alp(i)/REAL(nsplit_thermals)
414!            write(23,*)'ALP CALLTHERM',alp_bl(i),alp(i)
415        if(prt_level.GE.10) print*,'calltherm i alp_bl alp ale_bl ale',i,alp_bl(i),alp(i),ale_bl(i),ale(i)
416       ENDDO
417
418!IM 060508 marche pas comme cela !!!        enddo ! isplit
419
420!   tests sur les valeurs negatives de l'eau
421         nbptspb=0
422            DO k = 1, klev
423            DO i = 1, klon
424               logexpr2(i,k)=.not.q_seri(i,k).ge.0.
425               if (logexpr2(i,k)) then
426                q_seri(i,k)=1.e-15
427                nbptspb=nbptspb+1
428#ifdef ISO
429                do ixt=1,ntiso
430                  xt_seri(ixt,i,k)=1.e-15*(xtmemoire(ixt,i,k)/qmemoire(i,k))
431                enddo
432#endif
433!                if (prt_level.ge.10) then
434!                  print*,'WARN eau<0 apres therm i=',i,'  k=',k  &
435!    &         ,' dq,q',d_q_the(i,k),q_seri(i,k),  &
436!    &         'fm=',zfm_therm(i,k),'entr=',entr_therm(i,k)
437                 endif
438            ENDDO
439            ENDDO
440#ifdef ISO
441#ifdef ISOVERIF
442      if (iso_HDO.gt.0) then
443      call iso_verif_aberrant_enc_vect2D( &
444     &        xt_seri,q_seri, &
445     &        'calltherm 393, apres bidouille q<0',ntiso,klon,klev)
446      endif     
447#endif
448#endif
449
450        IF(nbptspb.GT.0) print*,'Number of points with q_seri(i,k)<=0 ',nbptspb   
451! tests sur les valeurs de la temperature
452        nbptspb=0
453            DO k = 1, klev
454            DO i = 1, klon
455               logexpr2(i,k)=t_seri(i,k).lt.50..or.t_seri(i,k).gt.370.
456               if (logexpr2(i,k)) nbptspb=nbptspb+1
457!              if ((t_seri(i,k).lt.50.) .or.  &
458!    &              (t_seri(i,k).gt.370.)) then
459!                 print*,'WARN temp apres therm i=',i,'  k=',k  &
460!    &         ,' t_seri',t_seri(i,k)
461!              CALL abort
462!              endif
463            ENDDO
464            ENDDO
465        IF(nbptspb.GT.0) print*,'Number of points with q_seri(i,k)<=0 ',nbptspb
466         enddo ! isplit
467
468!
469!***************************************************************
470!     calcul du flux ascencant conservatif
471!            print*,'<<<<calcul flux ascendant conservatif'
472
473      fmc_therm=0.
474               do k=1,klev
475            do i=1,klon
476                  if (entr_therm(i,k).gt.0.) then
477                     fmc_therm(i,k+1)=fmc_therm(i,k)+entr_therm(i,k)
478                  else
479                     fmc_therm(i,k+1)=fmc_therm(i,k)
480                  endif
481                  detrc_therm(i,k)=(fmc_therm(i,k+1)-fm_therm(i,k+1))  &
482     &                 -(fmc_therm(i,k)-fm_therm(i,k))
483               enddo
484            enddo
485     
486     
487!****************************************************************
488!     calcul de l'humidite dans l'ascendance
489!      print*,'<<<<calcul de lhumidite dans thermique'
490!CR:on ne le calcule que pour le cas sec
491      if (iflag_thermals.le.11) then     
492      do i=1,klon
493         zqasc(i,1)=q_seri(i,1)
494         do k=2,klev
495            if (fmc_therm(i,k+1).gt.1.e-6) then
496               zqasc(i,k)=(fmc_therm(i,k)*zqasc(i,k-1)  &
497     &              +entr_therm(i,k)*q_seri(i,k))/fmc_therm(i,k+1)
498!CR:test on asseche le thermique
499!               zqasc(i,k)=zqasc(i,k)/2.
500!            else
501!               zqasc(i,k)=q_seri(i,k)
502            endif
503         enddo
504       enddo
505     
506
507!     calcul de l'eau condensee dans l'ascendance
508!             print*,'<<<<calcul de leau condensee dans thermique'
509             do i=1,klon
510                do k=1,klev
511                   clwcon0(i,k)=zqasc(i,k)-zqsat(i,k)
512                   if (clwcon0(i,k).lt.0. .or.   &
513     &             (fm_therm(i,k+1)+detrc_therm(i,k)).lt.1.e-6) then
514                      clwcon0(i,k)=0.
515                   endif
516                enddo
517             enddo
518       else
519              do i=1,klon
520                do k=1,klev
521                   clwcon0(i,k)=zqla(i,k) 
522                   if (clwcon0(i,k).lt.0. .or.   &
523     &             (fm_therm(i,k+1)+detrc_therm(i,k)).lt.1.e-6) then
524                   clwcon0(i,k)=0.
525                   endif
526                enddo
527             enddo
528       endif
529!*******************************************************************   
530
531
532!jyg  Protection contre les temperatures nulles
533          do i=1,klon
534             do k=1,klev
535                if (ztla(i,k) .lt. 1.e-10) fraca(i,k) =0.
536             enddo
537          enddo
538
539
540      return
541
542      end
Note: See TracBrowser for help on using the repository browser.