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

Last change on this file since 3289 was 3158, checked in by csegonne, 13 months ago

MARS PCM
Cleaning of conduction.F, euvheat.F90, moldiff.F and molvis.F, some commented lines referring to a local calculation of layers/levels altitudes have been removed.

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