xref: /phasta/phSolver/AMG/ramg_smoother.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen!************************************************************
2*59599516SKenneth E. Jansen!      ramg_calcPPEv
3*59599516SKenneth E. Jansen!      calculate u = PPE*v ( parallelly using lesLib calls)
4*59599516SKenneth E. Jansen!************************************************************
5*59599516SKenneth E. Jansen      subroutine ramg_calcPPEv(colm,rowp,lhsK,lhsP,
6*59599516SKenneth E. Jansen     &                   ilwork,BC,iBC,iper,
7*59599516SKenneth E. Jansen     &                         u,v)
8*59599516SKenneth E. Jansen      use ramg_data
9*59599516SKenneth E. Jansen      include "common.h"
10*59599516SKenneth E. Jansen      include "mpif.h"
11*59599516SKenneth E. Jansen
12*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg+1) :: colm
13*59599516SKenneth E. Jansen      integer,intent(in),dimension(nnz_tot) :: rowp
14*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(9,nnz_tot) :: lhsK
15*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(4,nnz_tot) :: lhsP
16*59599516SKenneth E. Jansen      integer,intent(in),dimension(nlwork)  :: ilwork
17*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg)    :: iBC,iper
18*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC
19*59599516SKenneth E. Jansen
20*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg) :: v
21*59599516SKenneth E. Jansen      real(kind=8),intent(inout),dimension(nshg) :: u
22*59599516SKenneth E. Jansen
23*59599516SKenneth E. Jansen      real(kind=8),dimension(nshg,3) :: utmp
24*59599516SKenneth E. Jansen      real(kind=8),dimension(nshg,4) :: utmp4
25*59599516SKenneth E. Jansen      integer                         :: i,j
26*59599516SKenneth E. Jansen
27*59599516SKenneth E. Jansen      call commOut(v,ilwork,1,iper,iBC,BC)
28*59599516SKenneth E. Jansen      call fLesSparseApG(colm,rowp,lhsP,v,utmp,nshg,nnz_tot)
29*59599516SKenneth E. Jansen      call commIn(utmp,ilwork,3,iper,iBC,BC)
30*59599516SKenneth E. Jansen
31*59599516SKenneth E. Jansen      do i=1,3
32*59599516SKenneth E. Jansen         utmp(:,i) = utmp(:,i)*ramg_flowDiag%p(:,i)**2
33*59599516SKenneth E. Jansen      enddo
34*59599516SKenneth E. Jansen
35*59599516SKenneth E. Jansen      do i=1,3
36*59599516SKenneth E. Jansen         utmp4(:,i) = -utmp(:,i)
37*59599516SKenneth E. Jansen      enddo
38*59599516SKenneth E. Jansen
39*59599516SKenneth E. Jansen      utmp4(:,4) = v
40*59599516SKenneth E. Jansen
41*59599516SKenneth E. Jansen      call commOut(utmp4,ilwork,4,iper,iBC,BC)
42*59599516SKenneth E. Jansen      call fLesSparseApNGtC(colm,rowp,lhsP,utmp4,u,nshg,nnz_tot)
43*59599516SKenneth E. Jansen      call commIn(u,ilwork,1,iper,iBC,BC)
44*59599516SKenneth E. Jansen!     There is a slight modify at lesSparse.f:
45*59599516SKenneth E. Jansen!     row(20*nNodes) ==> row(nnz_tot)
46*59599516SKenneth E. Jansen
47*59599516SKenneth E. Jansen      end subroutine ! ramg_calcPPEv
48*59599516SKenneth E. Jansen
49*59599516SKenneth E. Jansen
50*59599516SKenneth E. Jansen!************************************************************
51*59599516SKenneth E. Jansen!      ramg_jacobi
52*59599516SKenneth E. Jansen!      calculate u = D^-1 * [ -( L+U ) v + r ]
53*59599516SKenneth E. Jansen!
54*59599516SKenneth E. Jansen!      and also the residule in and out resout = | r-Au |
55*59599516SKenneth E. Jansen!                                       resin = | r-Av |
56*59599516SKenneth E. Jansen!
57*59599516SKenneth E. Jansen!************************************************************
58*59599516SKenneth E. Jansen      subroutine ramg_jacobi(Acolm,Arowp,Alhs,Anshg,Annz_tot,!A
59*59599516SKenneth E. Jansen     &                      r,u,v)
60*59599516SKenneth E. Jansen
61*59599516SKenneth E. Jansen      include "common.h"
62*59599516SKenneth E. Jansen      integer,intent(in)                   :: Anshg,Annz_tot
63*59599516SKenneth E. Jansen      integer,intent(in),dimension(Anshg+1) :: Acolm
64*59599516SKenneth E. Jansen      integer,intent(in),dimension(Annz_tot) :: Arowp
65*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(Annz_tot) :: Alhs
66*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(Anshg) :: v,r
67*59599516SKenneth E. Jansen      real(kind=8),intent(inout),dimension(Anshg) :: u
68*59599516SKenneth E. Jansen      real(kind=8)                        :: tmp,tmpin,tmpout
69*59599516SKenneth E. Jansen
70*59599516SKenneth E. Jansen      integer                              :: i,j,k,p
71*59599516SKenneth E. Jansen
72*59599516SKenneth E. Jansen      real(kind=8)                     :: damp_jacobi
73*59599516SKenneth E. Jansen
74*59599516SKenneth E. Jansen      ! Useless function, will be removed
75*59599516SKenneth E. Jansen
76*59599516SKenneth E. Jansen      u = 0
77*59599516SKenneth E. Jansen      damp_jacobi = 1.0/ramg_relax
78*59599516SKenneth E. Jansen      do i=1,Anshg
79*59599516SKenneth E. Jansen         do k = Acolm(i)+1,Acolm(i+1)-1
80*59599516SKenneth E. Jansen            j = Arowp(k)
81*59599516SKenneth E. Jansen            tmp = Alhs(k)*v(j)
82*59599516SKenneth E. Jansen            u(i) = u(i) - tmp
83*59599516SKenneth E. Jansen         enddo
84*59599516SKenneth E. Jansen         u(i) = damp_jacobi*u(i) + r(i)
85*59599516SKenneth E. Jansen         u(i) = u(i)/Alhs(Acolm(i))
86*59599516SKenneth E. Jansen      enddo
87*59599516SKenneth E. Jansen
88*59599516SKenneth E. Jansen      end subroutine ! ramg_jacobi
89*59599516SKenneth E. Jansen
90*59599516SKenneth E. Jansen!***************************************************************
91*59599516SKenneth E. Jansen!     ramg_boundary_mls1:
92*59599516SKenneth E. Jansen!     Do polynomial smoothing on boundary nodes only.
93*59599516SKenneth E. Jansen!     remember this can always be simplified to boundary jacobi
94*59599516SKenneth E. Jansen!**************************************************************
95*59599516SKenneth E. Jansen      subroutine ramg_boundary_mls1(acolm,arowp,alhs,
96*59599516SKenneth E. Jansen     &                      anshg,annz_tot,!A
97*59599516SKenneth E. Jansen     &                      r,u,v,clevel,pflag,
98*59599516SKenneth E. Jansen     &                      ilwork,BC,iBC,iper)
99*59599516SKenneth E. Jansen      use ramg_data
100*59599516SKenneth E. Jansen      include "common.h"
101*59599516SKenneth E. Jansen      include "mpif.h"
102*59599516SKenneth E. Jansen      include "auxmpi.h"
103*59599516SKenneth E. Jansen
104*59599516SKenneth E. Jansen      integer,intent(in)                   :: anshg,annz_tot
105*59599516SKenneth E. Jansen      integer,intent(in),dimension(anshg+1) :: acolm
106*59599516SKenneth E. Jansen      integer,intent(in),dimension(annz_tot) :: arowp
107*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(annz_tot) :: alhs
108*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(anshg) :: v,r
109*59599516SKenneth E. Jansen      real(kind=8),intent(inout),dimension(anshg) :: u
110*59599516SKenneth E. Jansen      integer,intent(in)                      :: clevel
111*59599516SKenneth E. Jansen      integer,intent(in),dimension(nlwork)        :: ilwork
112*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg)            :: iBC,iper
113*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC
114*59599516SKenneth E. Jansen      logical,dimension(anshg),intent(inout) :: pflag
115*59599516SKenneth E. Jansen
116*59599516SKenneth E. Jansen      integer i,j,k,p,q
117*59599516SKenneth E. Jansen      real(kind=8),dimension(anshg)    :: tmp,v2,aux
118*59599516SKenneth E. Jansen      real(kind=8) :: cf1
119*59599516SKenneth E. Jansen      real(kind=8) :: cpusec(2)
120*59599516SKenneth E. Jansen
121*59599516SKenneth E. Jansen      if (numpe.eq.1) return
122*59599516SKenneth E. Jansen
123*59599516SKenneth E. Jansen      call cpu_time(cpusec(1))
124*59599516SKenneth E. Jansen
125*59599516SKenneth E. Jansen      tmp = 0
126*59599516SKenneth E. Jansen      v2 = v
127*59599516SKenneth E. Jansen
128*59599516SKenneth E. Jansen      ! tmp = A*x
129*59599516SKenneth E. Jansen      call ramg_calcAv_b(clevel,tmp,v2,ilwork,BC,iBC,iper,1,0)
130*59599516SKenneth E. Jansen!      call ramg_calcAv(acolm,arowp,alhs,anshg,annz_tot,
131*59599516SKenneth E. Jansen!     &                 tmp,v2,1)
132*59599516SKenneth E. Jansen      ! tmp = A*x-b=r
133*59599516SKenneth E. Jansen      tmp = tmp-r
134*59599516SKenneth E. Jansen      ! aux = A*r ! not necessary
135*59599516SKenneth E. Jansen      ! call ramg_calcAv_b(clevel,aux,tmp,ilwork,BC,iBC,iper,1,0)
136*59599516SKenneth E. Jansen
137*59599516SKenneth E. Jansen      cf1 = 1.1*mlsCf(clevel,1)
138*59599516SKenneth E. Jansen
139*59599516SKenneth E. Jansen      !v2 = mlsCf(clevel,1)*tmp!+mlsCf(clevel,2)*aux
140*59599516SKenneth E. Jansen      !v2 = 1.1*v2
141*59599516SKenneth E. Jansen      !v2 = tmp
142*59599516SKenneth E. Jansen
143*59599516SKenneth E. Jansen      do i=1,anshg
144*59599516SKenneth E. Jansen      p = amg_paraext(clevel)%p(i)
145*59599516SKenneth E. Jansen      q = amg_paramap(clevel)%p(i)
146*59599516SKenneth E. Jansen      if (p.ne.(myrank+1)) then
147*59599516SKenneth E. Jansen          pflag(i) = .true.
148*59599516SKenneth E. Jansen      if ((p.gt.0).or.(p.ne.q)) then ! master boundary or extended
149*59599516SKenneth E. Jansen         u(i) = v(i)-cf1*tmp(i)
150*59599516SKenneth E. Jansen      endif
151*59599516SKenneth E. Jansen      endif
152*59599516SKenneth E. Jansen      enddo
153*59599516SKenneth E. Jansen
154*59599516SKenneth E. Jansen      call cpu_time(cpusec(2))
155*59599516SKenneth E. Jansen      ramg_time(20) = ramg_time(20)+cpusec(2)-cpusec(1)
156*59599516SKenneth E. Jansen
157*59599516SKenneth E. Jansen      end subroutine !ramg_boundary_mls1
158*59599516SKenneth E. Jansen
159*59599516SKenneth E. Jansen
160*59599516SKenneth E. Jansen!************************************************************
161*59599516SKenneth E. Jansen!      ramg_gauss
162*59599516SKenneth E. Jansen!
163*59599516SKenneth E. Jansen!      forward and backward, defined in fwdbck
164*59599516SKenneth E. Jansen!      Gauss-Seidel smoothing, u = gaussseidel(M,r)
165*59599516SKenneth E. Jansen!
166*59599516SKenneth E. Jansen!************************************************************
167*59599516SKenneth E. Jansen      subroutine ramg_gauss(acolm,arowp,alhs,anshg,annz_tot,!A
168*59599516SKenneth E. Jansen     &                      r,u,v,fwdbck,clevel,
169*59599516SKenneth E. Jansen     &                      ilwork,BC,iBC,iper)
170*59599516SKenneth E. Jansen      use ramg_data
171*59599516SKenneth E. Jansen      include "common.h"
172*59599516SKenneth E. Jansen      include "mpif.h"
173*59599516SKenneth E. Jansen      include "auxmpi.h"
174*59599516SKenneth E. Jansen
175*59599516SKenneth E. Jansen      integer,intent(in)                   :: anshg,annz_tot
176*59599516SKenneth E. Jansen      integer,intent(in),dimension(anshg+1) :: acolm
177*59599516SKenneth E. Jansen      integer,intent(in),dimension(annz_tot) :: arowp
178*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(annz_tot) :: alhs
179*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(anshg) :: v,r
180*59599516SKenneth E. Jansen      real(kind=8),intent(inout),dimension(anshg) :: u
181*59599516SKenneth E. Jansen      integer,intent(in)                      :: fwdbck!1=fwd,2=bck
182*59599516SKenneth E. Jansen      integer,intent(in)                      :: clevel
183*59599516SKenneth E. Jansen      integer,intent(in),dimension(nlwork)        :: ilwork
184*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg)            :: iBC,iper
185*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC
186*59599516SKenneth E. Jansen
187*59599516SKenneth E. Jansen
188*59599516SKenneth E. Jansen      real(kind=8)                        :: tmp
189*59599516SKenneth E. Jansen      real(kind=8),dimension(anshg) :: u2
190*59599516SKenneth E. Jansen      integer                             :: istr,iend,iint
191*59599516SKenneth E. Jansen      integer                              :: i,j,k,p,ki,kj,kp
192*59599516SKenneth E. Jansen      logical                             :: postsmooth,p2
193*59599516SKenneth E. Jansen      logical,dimension(anshg) :: pflag
194*59599516SKenneth E. Jansen!      integer            :: itag,iacc,iother,numseg,isgbeg,itkbeg
195*59599516SKenneth E. Jansen!      integer            :: numtask
196*59599516SKenneth E. Jansen
197*59599516SKenneth E. Jansen      u = v
198*59599516SKenneth E. Jansen
199*59599516SKenneth E. Jansen      postsmooth = (fwdbck.eq.1)
200*59599516SKenneth E. Jansen      if (postsmooth) then ! post smoothing : 1->nshg
201*59599516SKenneth E. Jansen          istr = 1
202*59599516SKenneth E. Jansen          iend = anshg
203*59599516SKenneth E. Jansen          iint = 1
204*59599516SKenneth E. Jansen      else ! pre smoothing : nshg->1
205*59599516SKenneth E. Jansen          istr = anshg
206*59599516SKenneth E. Jansen          iend = 1
207*59599516SKenneth E. Jansen          iint = -1
208*59599516SKenneth E. Jansen      end if
209*59599516SKenneth E. Jansen
210*59599516SKenneth E. Jansen      pflag = .false.
211*59599516SKenneth E. Jansen      ! boundary
212*59599516SKenneth E. Jansen!      if (postsmooth) then
213*59599516SKenneth E. Jansen          call ramg_boundary_mls1(acolm,arowp,alhs,anshg,annz_tot,
214*59599516SKenneth E. Jansen     &         r,u,v,clevel,pflag,ilwork,BC,iBC,iper)
215*59599516SKenneth E. Jansen!      endif
216*59599516SKenneth E. Jansen!      if (.false.) then
217*59599516SKenneth E. Jansen
218*59599516SKenneth E. Jansen      do i=istr,iend,iint
219*59599516SKenneth E. Jansen         ki = CF_map(clevel)%p(i)
220*59599516SKenneth E. Jansen         p2 = .not.(pflag(ki)) ! nodes that have not been treated
221*59599516SKenneth E. Jansen!         if (numpe.gt.1) then
222*59599516SKenneth E. Jansen!            p2=p2.and.(amg_paraext(clevel)%p(ki).eq.(myrank+1))
223*59599516SKenneth E. Jansen!         endif
224*59599516SKenneth E. Jansen         if (p2) then
225*59599516SKenneth E. Jansen         tmp = 0
226*59599516SKenneth E. Jansen         do k = acolm(ki)+1,acolm(ki+1)-1
227*59599516SKenneth E. Jansen            j = arowp(k)
228*59599516SKenneth E. Jansen            if ((postsmooth).or.(pflag(j))) then
229*59599516SKenneth E. Jansen            ! save some if u(j)=0
230*59599516SKenneth E. Jansen            tmp = tmp + alhs(k)*u(j)
231*59599516SKenneth E. Jansen            endif
232*59599516SKenneth E. Jansen         enddo
233*59599516SKenneth E. Jansen         u(ki) = r(ki) - tmp
234*59599516SKenneth E. Jansen         pflag(ki) = .true.
235*59599516SKenneth E. Jansen         !u(ki) = u(ki)/alhs(acolm(ki))
236*59599516SKenneth E. Jansen         ! The above line is commented out because we know that
237*59599516SKenneth E. Jansen         ! A(i,i)=1.0, this is the result of prescaling of the matrix.
238*59599516SKenneth E. Jansen         endif !p2
239*59599516SKenneth E. Jansen      enddo
240*59599516SKenneth E. Jansen!      if (.not.postsmooth) then
241*59599516SKenneth E. Jansen      u2 = u
242*59599516SKenneth E. Jansen          call ramg_boundary_mls1(acolm,arowp,alhs,anshg,annz_tot,
243*59599516SKenneth E. Jansen     &         r,u,u2,clevel,pflag,ilwork,BC,iBC,iper)
244*59599516SKenneth E. Jansen!      endif
245*59599516SKenneth E. Jansen
246*59599516SKenneth E. Jansen      end subroutine ! ramg_gauss
247*59599516SKenneth E. Jansen
248*59599516SKenneth E. Jansen!**********************************************
249*59599516SKenneth E. Jansen!      ramg dense direct solver routines.
250*59599516SKenneth E. Jansen!
251*59599516SKenneth E. Jansen!      ramg_sparse2dense
252*59599516SKenneth E. Jansen!      Sparse to dense form
253*59599516SKenneth E. Jansen!
254*59599516SKenneth E. Jansen!      ramg_direct_LU ( setup only )
255*59599516SKenneth E. Jansen!      Outside package for direct solve sparse
256*59599516SKenneth E. Jansen!        matrix by LU
257*59599516SKenneth E. Jansen!**********************************************
258*59599516SKenneth E. Jansen
259*59599516SKenneth E. Jansen      subroutine ramg_direct_LU(Acolm,Arowp,Alhs,Anshg,Annz_tot)
260*59599516SKenneth E. Jansen      use ramg_data
261*59599516SKenneth E. Jansen
262*59599516SKenneth E. Jansen      integer,intent(in)                   :: Anshg,Annz_tot
263*59599516SKenneth E. Jansen      integer,intent(in),dimension(Anshg+1) :: Acolm
264*59599516SKenneth E. Jansen      integer,intent(in),dimension(Annz_tot) :: Arowp
265*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(Annz_tot) :: Alhs
266*59599516SKenneth E. Jansen
267*59599516SKenneth E. Jansen!      real(kind=8),dimension(Anshg,Anshg)      :: mtxA
268*59599516SKenneth E. Jansen!      integer,dimension(Anshg)             :: indx
269*59599516SKenneth E. Jansen      real(kind=8)                     :: d
270*59599516SKenneth E. Jansen
271*59599516SKenneth E. Jansen      integer :: i,j
272*59599516SKenneth E. Jansen
273*59599516SKenneth E. Jansen      !write(*,*) "ramg_setup_flag",ramg_setup_flag
274*59599516SKenneth E. Jansen      if (ramg_setup_flag .ne. 0) return
275*59599516SKenneth E. Jansen
276*59599516SKenneth E. Jansen      !write(*,*)'again'
277*59599516SKenneth E. Jansen      if (allocated(cmtxA)) deallocate(cmtxA)
278*59599516SKenneth E. Jansen      if (allocated(cindx)) deallocate(cindx)
279*59599516SKenneth E. Jansen
280*59599516SKenneth E. Jansen      allocate(cmtxA(Anshg,Anshg))
281*59599516SKenneth E. Jansen      allocate(cindx(Anshg))
282*59599516SKenneth E. Jansen
283*59599516SKenneth E. Jansen      if (myrank.eq.master) then
284*59599516SKenneth E. Jansen      write(*,*)"Start to setup LU decomposition at coarsest level"
285*59599516SKenneth E. Jansen      endif
286*59599516SKenneth E. Jansen      call ramg_sparse2dense(Acolm,Arowp,Alhs,Anshg,Annz_tot,cmtxA)
287*59599516SKenneth E. Jansen
288*59599516SKenneth E. Jansen!      Asol = Arhs
289*59599516SKenneth E. Jansen
290*59599516SKenneth E. Jansen      call ludcmp(cmtxA,Anshg,Anshg,cindx,d)
291*59599516SKenneth E. Jansen!      call lubksb(mtxA,Anshg,Anshg,indx,Asol)
292*59599516SKenneth E. Jansen      if (myrank.eq.master) then
293*59599516SKenneth E. Jansen          write(*,*)"End of setup LU decomposition at coarsest level"
294*59599516SKenneth E. Jansen      endif
295*59599516SKenneth E. Jansen
296*59599516SKenneth E. Jansen      end subroutine ! ramg_direct_LU
297*59599516SKenneth E. Jansen
298*59599516SKenneth E. Jansen      subroutine ramg_sparse2dense(Acolm,Arowp,Alhs,Anshg,Annz_tot,
299*59599516SKenneth E. Jansen     &                             mtxA)
300*59599516SKenneth E. Jansen
301*59599516SKenneth E. Jansen      integer,intent(in)                   :: Anshg,Annz_tot
302*59599516SKenneth E. Jansen      integer,intent(in),dimension(Anshg+1) :: Acolm
303*59599516SKenneth E. Jansen      integer,intent(in),dimension(Annz_tot) :: Arowp
304*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(Annz_tot) :: Alhs
305*59599516SKenneth E. Jansen      real(kind=8),intent(inout),dimension(Anshg,Anshg) :: mtxA
306*59599516SKenneth E. Jansen
307*59599516SKenneth E. Jansen      integer                       :: i,j,k,ip,jp,kp
308*59599516SKenneth E. Jansen
309*59599516SKenneth E. Jansen      mtxA = 0
310*59599516SKenneth E. Jansen      do i=1,Anshg
311*59599516SKenneth E. Jansen         do j= Acolm(i),Acolm(i+1)-1
312*59599516SKenneth E. Jansen            k = Arowp(j)
313*59599516SKenneth E. Jansen            mtxA(i,k) = Alhs(j)
314*59599516SKenneth E. Jansen         enddo
315*59599516SKenneth E. Jansen      enddo
316*59599516SKenneth E. Jansen
317*59599516SKenneth E. Jansen      end subroutine !ramg_sparse2dense
318*59599516SKenneth E. Jansen
319*59599516SKenneth E. Jansen
320*59599516SKenneth E. Jansen!**********************************************
321*59599516SKenneth E. Jansen!      ramg_direct_solve
322*59599516SKenneth E. Jansen!**********************************************
323*59599516SKenneth E. Jansen      subroutine ramg_direct_solve(
324*59599516SKenneth E. Jansen     &             Arhs,Asol,
325*59599516SKenneth E. Jansen     &             colm,rowp,lhsK,lhsP,
326*59599516SKenneth E. Jansen     &             ilwork,BC,iBC,iper)
327*59599516SKenneth E. Jansen      use ramg_data
328*59599516SKenneth E. Jansen      include "common.h"
329*59599516SKenneth E. Jansen      include "mpif.h"
330*59599516SKenneth E. Jansen      include "auxmpi.h"
331*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(amg_nshg(ramg_levelx))
332*59599516SKenneth E. Jansen     &                                                     ::Arhs
333*59599516SKenneth E. Jansen      real(kind=8),intent(inout),dimension(amg_nshg(ramg_levelx))
334*59599516SKenneth E. Jansen     &                                                     :: Asol
335*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg+1)       :: colm
336*59599516SKenneth E. Jansen      integer,intent(in),dimension(nnz_tot)      :: rowp
337*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(9,nnz_tot)  :: lhsK
338*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(4,nnz_tot)  :: lhsP
339*59599516SKenneth E. Jansen
340*59599516SKenneth E. Jansen      integer,intent(in),dimension(nlwork)        :: ilwork
341*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg)            :: iBC,iper
342*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC
343*59599516SKenneth E. Jansen
344*59599516SKenneth E. Jansen      real(kind=8)              :: mres_i,mres_n,mres_o
345*59599516SKenneth E. Jansen      real(kind=8),dimension(amg_nshg(ramg_levelx)) :: myV
346*59599516SKenneth E. Jansen      integer                      :: titer,cl
347*59599516SKenneth E. Jansen
348*59599516SKenneth E. Jansen!      Asol=Arhs
349*59599516SKenneth E. Jansen!      return
350*59599516SKenneth E. Jansen      cl = ramg_levelx
351*59599516SKenneth E. Jansen
352*59599516SKenneth E. Jansen      Asol = 0
353*59599516SKenneth E. Jansen      call ramg_L2_norm_p(cl,Arhs,1,mres_i)
354*59599516SKenneth E. Jansen      !mres_i=sqrt(flesDot1(Arhs,amg_nshg(cl),1))
355*59599516SKenneth E. Jansen      mres_o = 1e-7*mres_i
356*59599516SKenneth E. Jansen      titer = 0
357*59599516SKenneth E. Jansen
358*59599516SKenneth E. Jansen      !write(*,*)'mcheck direct solve begin,',myrank
359*59599516SKenneth E. Jansen      if (iamg_c_solver.eq.0) then
360*59599516SKenneth E. Jansen          ! Do two smoothing steps
361*59599516SKenneth E. Jansen          call ramg_smoother(cl,myV,Asol,Arhs,
362*59599516SKenneth E. Jansen     &              colm,rowp,lhsK,lhsP,
363*59599516SKenneth E. Jansen     &              ilwork,BC,iBC,iper,2)
364*59599516SKenneth E. Jansen          call ramg_smoother(cl,Asol,myV,Arhs,
365*59599516SKenneth E. Jansen     &              colm,rowp,lhsK,lhsP,
366*59599516SKenneth E. Jansen     &              ilwork,BC,iBC,iper,1)
367*59599516SKenneth E. Jansen      else if ( iamg_c_solver .eq. 1) then ! direct solve using G-S
368*59599516SKenneth E. Jansen          IF (.FALSE.) THEN ! smoother solve to exact
369*59599516SKenneth E. Jansen          mres_n = mres_i
370*59599516SKenneth E. Jansen          do while ( (mres_n.gt.mres_o).and.(titer.lt.1000) )
371*59599516SKenneth E. Jansen          call ramg_smoother(cl,myV,Asol,Arhs,
372*59599516SKenneth E. Jansen     &              colm,rowp,lhsK,lhsP,
373*59599516SKenneth E. Jansen     &              ilwork,BC,iBC,iper,2)
374*59599516SKenneth E. Jansen          call ramg_smoother(cl,Asol,myV,Arhs,
375*59599516SKenneth E. Jansen     &              colm,rowp,lhsK,lhsP,
376*59599516SKenneth E. Jansen     &              ilwork,BC,iBC,iper,1)
377*59599516SKenneth E. Jansen          call ramg_calcAv_g(cl,myV,Asol,colm,rowp,lhsK,lhsP,
378*59599516SKenneth E. Jansen     &         ilwork,BC,iBC,iper,1)
379*59599516SKenneth E. Jansen          myV = myV-Arhs
380*59599516SKenneth E. Jansen          mres_n=sqrt(flesDot1(myV,amg_nshg(cl),1))
381*59599516SKenneth E. Jansen          titer = titer + 1
382*59599516SKenneth E. Jansen          enddo
383*59599516SKenneth E. Jansen          ELSE ! cg solve to exact
384*59599516SKenneth E. Jansen              call ramg_CG(Asol,Arhs,cl,
385*59599516SKenneth E. Jansen     &        colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper,titer)
386*59599516SKenneth E. Jansen          ENDIF
387*59599516SKenneth E. Jansen      else ! direct solve using dense-matrix LU-decomposition
388*59599516SKenneth E. Jansen          Asol = Arhs
389*59599516SKenneth E. Jansen          call lubksb(cmtxA,amg_nshg(cl),amg_nshg(cl),cindx,Asol)
390*59599516SKenneth E. Jansen      endif
391*59599516SKenneth E. Jansen
392*59599516SKenneth E. Jansen      end subroutine ! ramg_direct_solve
393*59599516SKenneth E. Jansen
394*59599516SKenneth E. Jansen
395*59599516SKenneth E. Jansen      subroutine ramg_smoother(level,xnew,xold,xrhs,
396*59599516SKenneth E. Jansen     &                         colm,rowp,lhsK,lhsP,
397*59599516SKenneth E. Jansen     &                         ilwork,BC,iBC,iper,fwdbck)
398*59599516SKenneth E. Jansen      use ramg_data
399*59599516SKenneth E. Jansen      include "common.h"
400*59599516SKenneth E. Jansen      include "mpif.h"
401*59599516SKenneth E. Jansen      include "auxmpi.h"
402*59599516SKenneth E. Jansen
403*59599516SKenneth E. Jansen      integer level
404*59599516SKenneth E. Jansen      real(kind=8),dimension(amg_nshg(level)),intent(in) ::xrhs,xold
405*59599516SKenneth E. Jansen      real(kind=8),dimension(amg_nshg(level)),intent(inout) :: xnew
406*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg+1)       :: colm
407*59599516SKenneth E. Jansen      integer,intent(in),dimension(nnz_tot)      :: rowp
408*59599516SKenneth E. Jansen
409*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(9,nnz_tot)    :: lhsK
410*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(4,nnz_tot)    :: lhsP
411*59599516SKenneth E. Jansen
412*59599516SKenneth E. Jansen      integer,intent(in),dimension(nlwork)   :: ilwork
413*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg)     :: iBC,iper
414*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC
415*59599516SKenneth E. Jansen      integer fwdbck
416*59599516SKenneth E. Jansen
417*59599516SKenneth E. Jansen      integer i,k
418*59599516SKenneth E. Jansen      logical :: jacobi
419*59599516SKenneth E. Jansen
420*59599516SKenneth E. Jansen      ! NOTICE: fwdbck=2 : pre-smoothing
421*59599516SKenneth E. Jansen      ! NOTICE: fwdbck=1 : post-smoothing
422*59599516SKenneth E. Jansen
423*59599516SKenneth E. Jansen      jacobi = .false. !(ramg_relax.gt.0)
424*59599516SKenneth E. Jansen      ! Disable Jacobi, if needed, using 1-st order MLS instead
425*59599516SKenneth E. Jansen
426*59599516SKenneth E. Jansen      k = 1
427*59599516SKenneth E. Jansen      do i=1,k
428*59599516SKenneth E. Jansen          if ( iamg_smoother .eq. 3 ) then !.or. (level.gt.1) ) then
429*59599516SKenneth E. Jansen              call ramg_mls_apply(xnew,xold,xrhs,level,
430*59599516SKenneth E. Jansen     &                  colm,rowp,lhsK,lhsP,
431*59599516SKenneth E. Jansen     &                  ilwork,BC,iBC,iper,fwdbck)
432*59599516SKenneth E. Jansen          else if ( iamg_smoother .eq. 1 ) then
433*59599516SKenneth E. Jansen          call ramg_gauss(amg_A_colm(level)%p,amg_A_rowp(level)%p,
434*59599516SKenneth E. Jansen     &                    amg_A_lhs(level)%p,amg_nshg(level),
435*59599516SKenneth E. Jansen     &                    amg_nnz(level),xrhs,xnew,xold,
436*59599516SKenneth E. Jansen     &                    fwdbck,level,
437*59599516SKenneth E. Jansen     &                    ilwork,BC,iBC,iper)
438*59599516SKenneth E. Jansen          else if ( ( iamg_smoother .eq. 2) ) then !.and. (level.eq.1) ) then
439*59599516SKenneth E. Jansen              call ramg_cheby_apply(xnew,xold,xrhs,level,
440*59599516SKenneth E. Jansen     &             colm,rowp,lhsK,lhsP,
441*59599516SKenneth E. Jansen     &             ilwork,BC,iBC,iper,fwdbck)
442*59599516SKenneth E. Jansen          else if (jacobi) then
443*59599516SKenneth E. Jansen             call ramg_jacobi(amg_A_colm(level)%p,
444*59599516SKenneth E. Jansen     &                 amg_A_rowp(level)%p,amg_A_lhs(level)%p,
445*59599516SKenneth E. Jansen     &                 amg_nshg(level),amg_nnz(level),
446*59599516SKenneth E. Jansen     &                 xrhs,xnew,xold)
447*59599516SKenneth E. Jansen          endif
448*59599516SKenneth E. Jansen      enddo
449*59599516SKenneth E. Jansen
450*59599516SKenneth E. Jansen      end subroutine !ramg_smoother
451