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

Last change on this file since 4657 was 4657, checked in by fhourdin, 8 months ago

Poursuite nettoyage replauisation

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