source: trunk/LMDZ.MARS/libf/aeronomars/moldiff.F @ 317

Last change on this file since 317 was 38, checked in by emillour, 14 years ago

Ajout du modè Martien (mon LMDZ.MARS.BETA, du 28/01/2011) dans le rértoire mars, pour pouvoir suivre plus facilement les modifs.
EM

File size: 17.3 KB
Line 
1      subroutine moldiff(pplay,pplev,pt,pdt,pq,pdq,ptimestep,
2     &                   zzlay,pdteuv,pdtconduc,pdqdiff)
3
4
5      implicit none
6
7#include "dimensions.h"
8#include "dimphys.h"
9#include "comcstfi.h"
10#include "callkeys.h"
11#include "comdiurn.h"
12#include "chimiedata.h"
13#include "tracer.h"
14#include "conc.h"
15
16
17c
18c Input/Output
19c
20      real ptimestep
21      real pplay(ngridmx,nlayermx)
22      real zzlay(ngridmx,nlayermx)
23      real pplev(ngridmx,nlayermx+1)
24      real pq(ngridmx,nlayermx,nqmx)
25      real pdq(ngridmx,nlayermx,nqmx)
26      real pt(ngridmx,nlayermx)
27      real pdt(ngridmx,nlayermx)
28      real pdteuv(ngridmx,nlayermx)
29      real pdtconduc(ngridmx,nlayermx)
30      real pdqdiff(ngridmx,nlayermx,nqmx)
31c
32c Local
33c
34       real hco2(ncomptot),ho
35
36      integer ig,nz,l,n,nn,iq
37      real del1,del2, tmean ,dalfinvdz, d
38      real hh,dcoef,dcoef1,ptfac, ntot, dens, dens2, dens3
39      real hp(nlayermx)
40      real tt(nlayermx)
41      real qq(nlayermx,ncomptot)
42      real dmmeandz(nlayermx)
43      real qnew(nlayermx,ncomptot)
44      real zlocal(nlayermx)
45      real alf(ncomptot-1,ncomptot-1)
46      real alfinv(nlayermx,ncomptot-1,ncomptot-1)
47      real indx(ncomptot-1)
48      real b(nlayermx,ncomptot-1)
49      real y(ncomptot-1,ncomptot-1)
50      real aa(nlayermx,ncomptot-1,ncomptot-1)
51      real bb(nlayermx,ncomptot-1,ncomptot-1)
52      real cc(nlayermx,ncomptot-1,ncomptot-1)
53      real atri(nlayermx-2)
54      real btri(nlayermx-2)
55      real ctri(nlayermx-2)
56      real rtri(nlayermx-2)
57      real qtri(nlayermx-2)
58      real alfdiag(ncomptot-1)
59      real wi(ncomptot), flux(ncomptot), pote
60
61cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
62c     tracer numbering in the molecular diffusion
63cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
64c  Atomic oxygen must always be the LAST species of the list as
65c it is the dominant species at high altitudes. 
66      integer,parameter :: i_co   = 1
67      integer,parameter :: i_n2   = 2
68      integer,parameter :: i_o2   = 3
69      integer,parameter :: i_co2  = 4
70      integer,parameter :: i_h2   = 5
71      integer,parameter :: i_h    = 6
72      integer,parameter :: i_oh   = 7
73      integer,parameter :: i_ho2  = 8
74      integer,parameter :: i_h2o  = 9
75      integer,parameter :: i_h2o2 = 10
76      integer,parameter :: i_o1d  = 11
77      integer,parameter :: i_o3   = 12
78      integer,parameter :: i_ar   = 13
79      integer,parameter :: i_o    = 14
80
81! Tracer indexes in the GCM:
82      integer,save :: g_co2=0
83      integer,save :: g_co=0
84      integer,save :: g_o=0
85      integer,save :: g_o1d=0
86      integer,save :: g_o2=0
87      integer,save :: g_o3=0
88      integer,save :: g_h=0
89      integer,save :: g_h2=0
90      integer,save :: g_oh=0
91      integer,save :: g_ho2=0
92      integer,save :: g_h2o2=0
93      integer,save :: g_n2=0
94      integer,save :: g_ar=0
95      integer,save :: g_h2o=0
96
97      integer,save :: gcmind(ncomptot) ! array of GCM indexes
98      integer ierr
99
100      logical,save :: firstcall=.true.
101      real abfac(ncomptot)
102      real,save :: dij(ncomptot,ncomptot)
103
104! Initializations at first call
105      if (firstcall) then
106        call moldiffcoeff(dij)
107        print*,'MOLDIFF  EXO'
108       
109        ! identify the indexes of the tracers we'll need
110        g_co2=igcm_co2
111        if (g_co2.eq.0) then
112          write(*,*) "moldiff: Error; no CO2 tracer !!!"
113          stop
114        endif
115        g_co=igcm_co
116        if (g_co.eq.0) then
117          write(*,*) "moldiff: Error; no CO tracer !!!"
118          stop
119        endif
120        g_o=igcm_o
121        if (g_o.eq.0) then
122          write(*,*) "moldiff: Error; no O tracer !!!"
123          stop
124        endif
125        g_o1d=igcm_o1d
126        if (g_o1d.eq.0) then
127          write(*,*) "moldiff: Error; no O1D tracer !!!"
128          stop
129        endif
130        g_o2=igcm_o2
131        if (g_o2.eq.0) then
132          write(*,*) "moldiff: Error; no O2 tracer !!!"
133          stop
134        endif
135        g_o3=igcm_o3
136        if (g_o3.eq.0) then
137          write(*,*) "moldiff: Error; no O3 tracer !!!"
138          stop
139        endif
140        g_h=igcm_h
141        if (g_h.eq.0) then
142          write(*,*) "moldiff: Error; no H tracer !!!"
143          stop
144        endif
145        g_h2=igcm_h2
146        if (g_h2.eq.0) then
147          write(*,*) "moldiff: Error; no H2 tracer !!!"
148          stop
149        endif
150        g_oh=igcm_oh
151        if (g_oh.eq.0) then
152          write(*,*) "moldiff: Error; no OH tracer !!!"
153          stop
154        endif
155        g_ho2=igcm_ho2
156        if (g_ho2.eq.0) then
157          write(*,*) "moldiff: Error; no HO2 tracer !!!"
158          stop
159        endif
160        g_h2o2=igcm_h2o2
161        if (g_h2o2.eq.0) then
162          write(*,*) "moldiff: Error; no H2O2 tracer !!!"
163          stop
164        endif
165        g_n2=igcm_n2
166        if (g_n2.eq.0) then
167          write(*,*) "moldiff: Error; no N2 tracer !!!"
168          stop
169        endif
170        g_ar=igcm_ar
171        if (g_ar.eq.0) then
172          write(*,*) "moldiff: Error; no AR tracer !!!"
173          stop
174        endif
175        g_h2o=igcm_h2o_vap
176        if (g_h2o.eq.0) then
177          write(*,*) "moldiff: Error; no water vapor tracer !!!"
178          stop
179        endif
180
181cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
182c   fill array to relate local indexes to gcm indexes
183cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
184
185        gcmind(i_co)    = g_co
186        gcmind(i_n2)    = g_n2
187        gcmind(i_o2)    = g_o2
188        gcmind(i_co2)   = g_co2
189        gcmind(i_h2)    = g_h2
190        gcmind(i_h)     = g_h
191        gcmind(i_oh)    = g_oh
192        gcmind(i_ho2)   = g_ho2
193        gcmind(i_h2o)   = g_h2o
194        gcmind(i_h2o2)  = g_h2o2
195        gcmind(i_o1d)   = g_o1d
196        gcmind(i_o3)    = g_o3
197        gcmind(i_o)     = g_o
198        gcmind(i_ar)    = g_ar
199
200        firstcall= .false.
201      endif ! of if (firstcall)
202
203
204
205c
206cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
207
208      nz=nlayermx
209
210      do ig=1,ngridmx
211
212        do l=2,nz-1
213          tt(l)=pt(ig,l)+pdt(ig,l)*ptimestep
214     &                  +pdteuv(ig,l)*ptimestep
215     &                  +pdtconduc(ig,l)*ptimestep
216
217          do nn=1,ncomptot
218            qq(l,nn)=pq(ig,l,gcmind(nn))+pdq(ig,l,gcmind(nn))*ptimestep
219            qq(l,nn)=max(qq(l,nn),1.e-30)
220          enddo
221          hp(l)=-log(pplay(ig,l+1)/pplay(ig,l-1))
222          dmmeandz(l)=(mmean(ig,l+1)-mmean(ig,l-1))/hp(l)
223        enddo
224
225        tt(1)=pt(ig,1)  +pdt(ig,1)*ptimestep
226     &                  +pdteuv(ig,1)*ptimestep
227     &                  +pdtconduc(ig,1)*ptimestep
228        tt(nz)=pt(ig,nz)+pdt(ig,nz)*ptimestep
229     &                  +pdteuv(ig,nz)*ptimestep
230     &                  +pdtconduc(ig,nz)*ptimestep
231
232        do nn=1,ncomptot
233          qq(1,nn)=pq(ig,1,gcmind(nn))+pdq(ig,1,gcmind(nn))*ptimestep
234          qq(nz,nn)=pq(ig,nz,gcmind(nn))+pdq(ig,nz,gcmind(nn))*ptimestep
235          qq(1,nn)=max(qq(1,nn),1.e-30)
236          qq(nz,nn)=max(qq(nz,nn),1.e-30)
237        enddo
238        hp(1)=-log(pplay(ig,2)/pplay(ig,1))
239        dmmeandz(1)=(-3.*mmean(ig,1)+4.*mmean(ig,2)
240     &               -mmean(ig,3))/(2.*hp(1))
241        hp(nz)=-log(pplay(ig,nz)/pplay(ig,nz-1))
242        dmmeandz(nz)=(3.*mmean(ig,nz)-4.*mmean(ig,nz-1)
243     &                +mmean(ig,nz-2))/(2.*hp(nz))
244c
245c Setting-up matrix of alfa coefficients from Dickinson and Ridley 1972
246c
247      do l=1,nz
248       if(abs(dmmeandz(l)) .lt. 1.e-5) dmmeandz(l)=0.0
249        hh=rnew(ig,l)*tt(l)/g
250        ptfac=(1.e5/pplay(ig,l))*(tt(l)/273)**1.75
251        ntot=pplay(ig,l)/(1.381e-23*tt(l))                      ! in #/m3
252
253        do nn=1,ncomptot-1            ! rows
254          alfdiag(nn)=0.
255          dcoef1=dij(nn,i_o)*ptfac
256          do n=1,ncomptot-1           ! columns
257            y(nn,n)=0.
258            dcoef=dij(nn,n)*ptfac
259            alf(nn,n)=qq(l,nn)/ntot/1.66e-27
260     &         *(1./(mmol(gcmind(n))*dcoef)-1./(mmol(g_o)*dcoef1))
261            alfdiag(nn)=alfdiag(nn)
262     &       +(1./(mmol(gcmind(n))*dcoef)-1./(mmol(g_o)*dcoef1))
263     &        *qq(l,n)
264          enddo
265          dcoef=dij(nn,nn)*ptfac
266          alfdiag(nn)=alfdiag(nn)
267     &       -(1./(mmol(gcmind(nn))*dcoef)-1./(mmol(g_o)*dcoef1))
268     &          *qq(l,nn)
269          alf(nn,nn)=-(alfdiag(nn)
270     &                 +1./(mmol(g_o)*dcoef1))/ntot/1.66e-27
271          y(nn,nn)=1.
272          b(l,nn)=-(dmmeandz(l)/mmean(ig,l)
273     &              +mmol(gcmind(nn))/mmean(ig,l)-1.)
274        enddo
275
276c
277c Inverting the alfa matrix
278c
279        call ludcmp(alf,ncomptot-1,ncomptot-1,indx,d,ierr)
280
281c       TEMPORAIRE *****************************
282        if (ierr.ne.0) then
283            write(*,*) 'In moldiff: Problem in LUDCMP with matrix alf'
284            write(*,*) 'Singular matrix ?'
285c           write(*,*) 'Matrix alf = ', alf
286            write(*,*) 'ig, l=',ig, l
287            write(*,*) 'No molecular diffusion this time !'
288            call zerophys(ngridmx*nlayermx*nqmx,pdqdiff)
289            return
290c           stop
291        end if
292c       *******************************************
293        do n=1,ncomptot-1
294          call lubksb(alf,ncomptot-1,ncomptot-1,indx,y(1,n))
295          do nn=1,ncomptot-1
296            alfinv(l,nn,n)=y(nn,n)/hh
297          enddo
298        enddo
299      enddo
300
301c
302c Calculating coefficients of the system
303c
304
305c      zlocal(1)=-log(pplay(ig,1)/pplev(ig,1))* Rnew(ig,1)*tt(1)/g
306      zlocal(1)=zzlay(ig,1)
307
308      do l=2,nz-1
309        del1=hp(l)*pplay(ig,l)/(g*ptimestep)
310        del2=(hp(l)/2)**2*pplay(ig,l)/(g*ptimestep)
311        do nn=1,ncomptot-1
312          do n=1,ncomptot-1
313            dalfinvdz=(alfinv(l+1,nn,n)-alfinv(l-1,nn,n))/hp(l)
314            aa(l,nn,n)=-dalfinvdz/del1+alfinv(l,nn,n)/del2
315     &                +alfinv(l-1,nn,n)*b(l-1,n)/del1   
316            bb(l,nn,n)=-2.*alfinv(l,nn,n)/del2
317            cc(l,nn,n)=dalfinvdz/del1+alfinv(l,nn,n)/del2
318     &                -alfinv(l+1,nn,n)*b(l+1,n)/del1   
319          enddo
320        enddo
321
322c        tmean=tt(l)
323c        if(tt(l).ne.tt(l-1))
324c     &        tmean=(tt(l)-tt(l-1))/log(tt(l)/tt(l-1))
325c        zlocal(l)= zlocal(l-1)
326c     &         -log(pplay(ig,l)/pplay(ig,l-1))*rnew(ig,l)*tmean/g
327      zlocal(l)=zzlay(ig,l)
328      enddo
329
330c      zlocal(nz)= zlocal(nz-1)
331c     &         -log(pplay(ig,nz)/pplay(ig,nz-1))*rnew(ig,nz)*tmean/g
332      zlocal(nz)=zzlay(ig,nz)
333       
334ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
335c Escape velocity from Jeans equation for the flux of H and H2
336c (Hunten 1973, eq. 5)
337     
338      do n=1,ncomptot
339        wi(n)=1.
340        flux(n)=0.
341        abfac(n)=1.
342      enddo
343
344       dens=pplay(ig,nz)/(rnew(ig,nz)*tt(nz))
345c
346c For H:
347c
348       pote=(3398000.+zlocal(nz))/
349     &         (1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h)*g))
350       wi(i_h)=sqrt(2.*1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h)))
351     &             /(2.*sqrt(3.1415))*(1.+pote)*exp(-pote)
352       flux(i_h)=qq(nz,i_h)*dens/(1.6605e-27*mmol(g_h))*wi(i_h)
353       flux(i_h)=flux(i_h)*1.6606e-27
354       abfac(i_h)=0.
355c
356c For H2:
357c
358       pote=(3398000.+zlocal(nz))/
359     &         (1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h2)*g))
360       wi(i_h2)=sqrt(2.*1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h2)))
361     &              /(2.*sqrt(3.1415))*(1.+pote)*exp(-pote)
362       flux(i_h2)=qq(nz,i_h2)*dens/(1.6605e-27*mmol(g_h2))*wi(i_h2)
363       flux(i_h2)=flux(i_h2)*1.6606e-27
364       abfac(i_h2)=0.
365
366c ********* TEMPORAIRE : no escape for h and h2
367c     do n=1,ncomptot
368c       wi(n)=1.
369c       flux(n)=0.
370c       abfac(n)=1.
371c     enddo
372c ********************************************
373
374
375ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
376
377c
378c Setting coefficients for tridiagonal matrix and solving the system
379c
380
381      do nn=1,ncomptot-1
382        do l=2,nz-1
383          atri(l-1)=aa(l,nn,nn)
384          btri(l-1)=bb(l,nn,nn)+1.
385          ctri(l-1)=cc(l,nn,nn)
386          rtri(l-1)=qq(l,nn)
387          do n=1,ncomptot-1
388            rtri(l-1)=rtri(l-1)-(aa(l,nn,n)*qq(l-1,n)
389     &                           +bb(l,nn,n)*qq(l,n)
390     &                           +cc(l,nn,n)*qq(l+1,n))
391          enddo
392          rtri(l-1)=rtri(l-1)+(aa(l,nn,nn)*qq(l-1,nn)
393     &                         +bb(l,nn,nn)*qq(l,nn)
394     &                         +cc(l,nn,nn)*qq(l+1,nn))
395        enddo
396
397c
398c Boundary conditions:
399c       Escape flux for H and H2 at top
400c       Diffusive equilibrium for the other species at top
401c       Perfect mixing for all at bottom
402c
403
404        rtri(nz-2)=rtri(nz-2)
405     &             -ctri(nz-2)*flux(nn)*mmol(gcmind(nn))/(dens*wi(nn))
406
407        atri(nz-2)=atri(nz-2)
408     &             -abfac(nn)*ctri(nz-2)/(3.-2.*hp(nz)*b(nz,nn))
409        btri(nz-2)=btri(nz-2)
410     &             +abfac(nn)*4.*ctri(nz-2)/(3.-2.*hp(nz)*b(nz,nn))
411
412c        rtri(1)=rtri(1)-atri(1)*qq(1,nn)
413        btri(1)=btri(1)+atri(1)
414
415        call tridag(atri,btri,ctri,rtri,qtri,nz-2)
416
417        do l=2,nz-1
418c          qnew(l,nn)=qtri(l-1)
419          qnew(l,nn)=max(qtri(l-1),1.e-30)
420        enddo
421
422        qnew(nz,nn)=flux(nn)*mmol(gcmind(nn))/(dens*wi(nn))
423     &               +abfac(nn)*(4.*qnew(nz-1,nn)-qnew(nz-2,nn))
424     &                /(3.-2.*hp(nz)*b(nz,nn))
425c        qnew(1,nn)=qq(1,nn)
426        qnew(1,nn)=qnew(2,nn)
427         
428        qnew(nz,nn)=max(qnew(nz,nn),1.e-30)
429        qnew(1,nn)=max(qnew(1,nn),1.e-30)
430
431      enddo     ! loop on species
432
433      DO l=1,nz
434        if(zlocal(l).gt.65000.)then
435        pdqdiff(ig,l,g_o)=0.
436        do n=1,ncomptot-1
437          pdqdiff(ig,l,gcmind(n))=(qnew(l,n)-qq(l,n))
438          pdqdiff(ig,l,g_o)=pdqdiff(ig,l,g_o)-(qnew(l,n)-qq(l,n))
439          pdqdiff(ig,l,gcmind(n))=pdqdiff(ig,l,gcmind(n))/ptimestep
440        enddo
441                pdqdiff(ig,l,g_o)=pdqdiff(ig,l,g_o)/ptimestep
442        endif
443      ENDDO
444
445c      do l=2,nz
446c        do n=1,ncomptot-1
447c          hco2(n)=qnew(l,n)*pplay(ig,l)/(rnew(ig,l)*tt(l)) /
448c     &      (qnew(l-1,n)*pplay(ig,l-1)/(rnew(ig,l-1)*tt(l-1)))
449c          hco2(n)=-(zlocal(l)-zlocal(l-1))/log(hco2(n))/1000.
450c        enddo
451c        write(225,*),l,pt(1,l),(hco2(n),n=1,6)
452c        write(226,*),l,pt(1,l),(hco2(n),n=7,12)
453c      enddo
454
455      enddo             ! ig loop
456
457      return
458      end
459
460c    ********************************************************************
461c    ********************************************************************
462c    ********************************************************************
463 
464      subroutine tridag(a,b,c,r,u,n)
465      parameter (nmax=100)
466c      dimension gam(nmax),a(n),b(n),c(n),r(n),u(n)
467      real gam(nmax),a(n),b(n),c(n),r(n),u(n)
468      if(b(1).eq.0.)then
469        stop 'tridag: error: b(1)=0 !!! '
470      endif
471      bet=b(1)
472      u(1)=r(1)/bet
473      do 11 j=2,n
474        gam(j)=c(j-1)/bet
475        bet=b(j)-a(j)*gam(j)
476        if(bet.eq.0.) then
477          stop 'tridag: error: bet=0 !!! '
478        endif
479        u(j)=(r(j)-a(j)*u(j-1))/bet
48011    continue
481      do 12 j=n-1,1,-1
482        u(j)=u(j)-gam(j+1)*u(j+1)
48312    continue
484      return
485      end
486
487c    ********************************************************************
488c    ********************************************************************
489c    ********************************************************************
490
491      SUBROUTINE LUBKSB(A,N,NP,INDX,B)
492
493      implicit none
494
495      integer i,j,n,np,ii,ll
496      real sum
497      real a(np,np),indx(np),b(np)
498
499c      DIMENSION A(NP,NP),INDX(N),B(N)
500      II=0
501      DO 12 I=1,N
502        LL=INDX(I)
503        SUM=B(LL)
504        B(LL)=B(I)
505        IF (II.NE.0)THEN
506          DO 11 J=II,I-1
507            SUM=SUM-A(I,J)*B(J)
50811        CONTINUE
509        ELSE IF (SUM.NE.0.) THEN
510          II=I
511        ENDIF
512        B(I)=SUM
51312    CONTINUE
514      DO 14 I=N,1,-1
515        SUM=B(I)
516        IF(I.LT.N)THEN
517          DO 13 J=I+1,N
518            SUM=SUM-A(I,J)*B(J)
51913        CONTINUE
520        ENDIF
521        B(I)=SUM/A(I,I)
52214    CONTINUE
523      RETURN
524      END
525
526c    ********************************************************************
527c    ********************************************************************
528c    ********************************************************************
529
530      SUBROUTINE LUDCMP(A,N,NP,INDX,D,ierr)
531
532      implicit none
533
534      integer n,np,nmax,i,j,k,imax
535      real d,tiny,aamax
536      real a(np,np),indx(np)
537      integer ierr  ! error =0 if OK, =1 if problem
538
539      PARAMETER (NMAX=100,TINY=1.0E-20)
540c      DIMENSION A(NP,NP),INDX(N),VV(NMAX)
541      real sum,vv(nmax),dum
542
543      D=1.
544      DO 12 I=1,N
545        AAMAX=0.
546        DO 11 J=1,N
547          IF (ABS(A(I,J)).GT.AAMAX) AAMAX=ABS(A(I,J))
54811      CONTINUE
549        IF (AAMAX.EQ.0.) then
550            write(*,*) 'In moldiff: Problem in LUDCMP with matrix A'
551            write(*,*) 'Singular matrix ?'
552c           write(*,*) 'Matrix A = ', A
553c           TO DEBUG :
554            ierr =1
555            return
556c           stop
557        END IF
558
559        VV(I)=1./AAMAX
56012    CONTINUE
561      DO 19 J=1,N
562        IF (J.GT.1) THEN
563          DO 14 I=1,J-1
564            SUM=A(I,J)
565            IF (I.GT.1)THEN
566              DO 13 K=1,I-1
567                SUM=SUM-A(I,K)*A(K,J)
56813            CONTINUE
569              A(I,J)=SUM
570            ENDIF
57114        CONTINUE
572        ENDIF
573        AAMAX=0.
574        DO 16 I=J,N
575          SUM=A(I,J)
576          IF (J.GT.1)THEN
577            DO 15 K=1,J-1
578              SUM=SUM-A(I,K)*A(K,J)
57915          CONTINUE
580            A(I,J)=SUM
581          ENDIF
582          DUM=VV(I)*ABS(SUM)
583          IF (DUM.GE.AAMAX) THEN
584            IMAX=I
585            AAMAX=DUM
586          ENDIF
58716      CONTINUE
588        IF (J.NE.IMAX)THEN
589          DO 17 K=1,N
590            DUM=A(IMAX,K)
591            A(IMAX,K)=A(J,K)
592            A(J,K)=DUM
59317        CONTINUE
594          D=-D
595          VV(IMAX)=VV(J)
596        ENDIF
597        INDX(J)=IMAX
598        IF(J.NE.N)THEN
599          IF(A(J,J).EQ.0.)A(J,J)=TINY
600          DUM=1./A(J,J)
601          DO 18 I=J+1,N
602            A(I,J)=A(I,J)*DUM
60318        CONTINUE
604        ENDIF
60519    CONTINUE
606      IF(A(N,N).EQ.0.)A(N,N)=TINY
607      ierr =0
608      RETURN
609      END
610
Note: See TracBrowser for help on using the repository browser.