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

Last change on this file since 784 was 690, checked in by acolaitis, 13 years ago

Code re-organization in diverse parts of the GCM code. These are NOT cosmetic changes, but are needed for compilation of the Mesoscale model in NESTED configuration

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