source: src/solvers/pcg.F90 @ 13bb187

Revision 13bb187, 32.8 KB checked in by Oleg Batrashev <ogbash@…>, 2 years ago (diff)

Remove coarse preliminary factorisation, move coarse solve into levels>1.

  • Property mode set to 100644
Line 
1! DOUG - Domain decomposition On Unstructured Grids
2! Copyright (C) 1998-2006 Faculty of Computer Science, University of Tartu and
3! Department of Mathematics, University of Bath
4!
5! This library is free software; you can redistribute it and/or
6! modify it under the terms of the GNU Lesser General Public
7! License as published by the Free Software Foundation; either
8! version 2.1 of the License, or (at your option) any later version.
9!
10! This library is distributed in the hope that it will be useful,
11! but WITHOUT ANY WARRANTY; without even the implied warranty of
12! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
13! Lesser General Public License for more details.
14!
15! You should have received a copy of the GNU Lesser General Public
16! License along with this library; if not, write to the Free Software
17! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
18! or contact the authors (University of Tartu, Faculty of Computer Science, Chair
19! of Distributed Systems, Liivi 2, 50409 Tartu, Estonia, http://dougdevel.org,
20! mailto:info(at)dougdevel.org)
21
22!!----------------------------------
23!! Preconditioned Conjugate Gradient
24!!----------------------------------
25module pcg_mod
26
27  use ConvInf_mod
28  use SpMtx_mods
29  use Mesh_class
30  use globals
31
32  implicit none
33
34#include<doug_config.h>
35
36#ifdef D_COMPLEX
37#define float complex
38#else
39#define float real
40#endif
41
42  integer,save :: coarseid=0
43  private :: &
44       preconditioner, &
45       msolve
46
47  ! profiling
48  real(kind=rk), private, save :: time_preconditioner = 0
49
50contains
51
52  !--------------------------------------------
53  ! Preconditioned conjugate gradient method
54  !--------------------------------------------
55  subroutine pcg (A, b, x, Msh, tol_, maxit_, &
56       x0_, solinf, resvects_, CoarseMtx_)
57    implicit none
58
59    type(SpMtx),intent(in out) :: A ! System matrix (sparse)
60    float(kind=rk),dimension(:),pointer :: b ! RHS
61    float(kind=rk),dimension(:),pointer :: x ! Solution
62    ! Mesh - aux data for Ax operation
63    type(Mesh),                       intent(in) :: Msh
64
65    ! Optional arguments
66    ! Tolerance of the method
67    real(kind=rk),          intent(in), optional :: tol_
68    ! Max number of iterations
69    integer,                intent(in), optional :: maxit_
70    ! Initial guess
71    float(kind=rk), dimension(:), intent(in), optional :: x0_
72    ! Solution statistics
73    type(ConvInf),            intent(in out), optional :: solinf
74    ! Fill in the 'resvect' or not
75    logical,                      intent(in), optional :: resvects_
76    type(SpMtx),optional                         :: CoarseMtx_ ! Coarse matrix
77
78    real(kind=rk) :: tol   ! Tolerance
79    integer       :: maxit ! Max number of iterations
80    logical       :: resvects=.false.
81
82    !real(kind=rk), dimension(size(b)) :: r, p, q, m, z, b_k
83    real(kind=rk),dimension(:),pointer :: r, p, q, m, z, b_k
84    real(kind=rk) :: rho_curr, rho_prev, init_norm, res_norm
85    real(kind=rk) :: alpha, beta, tmp, ratio_norm
86    integer :: iter,nfreds
87    integer :: ierr
88    ! + test
89    float(kind=rk), dimension(:), pointer :: arr_copy, gv_tmp
90
91!!$    real(kind=rk) :: r_max
92
93    write(stream,'(/a)') 'Preconditioned conjugate gradient:'
94
95    if (size(b) /= size(x)) &
96         call DOUG_abort('[pcg] : SEVERE : size(b) /= size(x)',-1)
97
98!!$    if (present(x0_))
99!!$         ! use method with initial guess
100
101
102    nfreds=size(b)
103    allocate(r(nfreds))
104    allocate(p(nfreds))
105    allocate(q(nfreds))
106    allocate(m(nfreds))
107    allocate(z(nfreds))
108    allocate(b_k(nfreds))
109
110
111    tol = sctls%solve_tolerance
112    if (present(tol_)) &
113         tol = tol_
114
115    maxit = sctls%solve_maxiters
116    if (maxit==-1) maxit=100000
117    if (present(maxit_)) &
118         maxit = maxit_
119
120    if (present(resvects_).and.(.not.present(solinf))) &
121         call DOUG_abort('[pcg] : SEVERE : "resvects_" must be given along'//&
122         ' with "solinf".',-1)
123
124    ! Allocate convergence info structure
125    if (present(solinf).and.(.not.present(resvects_))) then
126       call ConvInf_Init(solinf)
127    else if (present(solinf)) then
128       if (present(resvects_).and.(resvects_.eqv.(.true.))) then
129         call ConvInf_Init(solinf, maxit)
130         resvects = .true.
131       end if
132    end if
133
134    ! Initialise auxiliary data structures
135    ! to assist with pmvm
136    call pmvmCommStructs_init(A, Msh)
137
138    ! + test
139    if (ismaster()) &
140         allocate(gv_tmp(Msh%ngf))
141    allocate(arr_copy(Msh%nlf))
142
143    ! Build preconditioner
144    !  call preconditioner(A, m, CoarseMtx_)
145!!$    ! + test
146!!$    arr_copy = m
147!!$    call Vect_newToOldPerm(arr_copy)
148!!$    gv_tmp = 0.0_rk
149!!$    call Vect_Gather(arr_copy, gv_tmp, Msh)
150!!$    if (ismaster()) &
151!!$         call Vect_Print(gv_tmp, 'm ::')
152
153!!$    ! + test
154!!$    !b = 1.0_rk
155!!$    arr_copy = b
156!!$    call Vect_newToOldPerm(arr_copy)
157!!$    call Vect_Print(arr_copy, 'b local (original ordering) ::')
158!!$    gv_tmp = 0.0_rk
159!!$    call Vect_Gather(arr_copy, gv_tmp, Msh)
160!!$    if (ismaster()) &
161!!$         call Vect_Print(gv_tmp, 'b global (original ordering) ::')
162
163!!$    ! + test
164!!$    arr_copy = x
165!!$    call Vect_newToOldPerm(arr_copy)
166!!$    call Vect_Print(arr_copy, 'x local (original ordering) ::')
167!!$    gv_tmp = 0.0_rk
168!!$    call Vect_Gather(arr_copy, gv_tmp, Msh)
169!!$    if (ismaster()) &
170!!$         call Vect_Print(gv_tmp, 'x global (original ordering) ::')
171
172    !r = b - SpMtx_pmvm(A, x, Msh)
173    call SpMtx_pmvm(r,A,x,Msh)
174    r = b - r
175    !call Vect_Print(r,'initial residual')
176
177    ! + test
178!!$    arr_copy = r
179!!$    call Vect_newToOldPerm(arr_copy)
180!!$    gv_tmp = 0.0_rk
181!!$    call Vect_Gather(arr_copy, gv_tmp, Msh)
182!!$    if (ismaster()) &
183!!$         call Vect_Print(gv_tmp, 'r ::')
184
185    init_norm = Vect_dot_product(r,r)
186    if (init_norm == 0.0) &
187         init_norm = 1.0
188
189    write(stream,'(a,i6,a,e8.2)') 'maxit = ',maxit,', tol = ',tol
190    write(stream,'(a,e22.15)') 'init_norm = ', dsqrt(init_norm)
191
192    iter = 1
193    ratio_norm = 1.0_rk
194    do while((ratio_norm > tol*tol).and.(iter <= maxit))
195       call msolve(m,r,z)
196       !if (present(CoarseMtx_)) then
197       !  call preconditioner(z,A,r,CoarseMtx_)
198       !  !z=r
199       !else
200       !  call preconditioner(z,A,r)
201       !endif
202
203       !call Vect_Print(z,'preconditioner')
204
205       ! compute current rho
206       rho_curr = Vect_dot_product(r,z)
207       if (iter == 1) then
208          p = z
209       else
210          beta = rho_curr / rho_prev
211          p = z + beta * p
212       end if
213       !q = SpMtx_pmvm(A, p, Msh)
214       call SpMtx_pmvm(q,A,p,Msh)
215       ! compute alpha
216       alpha = rho_curr / Vect_dot_product(p,q)
217       x = x + alpha * p
218       r = r - alpha * q
219       rho_prev = rho_curr
220       ! check
221       res_norm = Vect_dot_product(r,r)
222       ratio_norm = dsqrt(res_norm) / dsqrt(init_norm)
223!!$       ratio_norm = res_norm / init_norm
224       if (ismaster()) &
225            write(stream, '(i5,a,e22.15)') iter,': res_norm=',dsqrt(res_norm)
226            !write(stream, '(i5,a,e22.15,e22.15)') iter,': dsqrt(ratio_norm)=',&
227            !dsqrt(ratio_norm), dsqrt(res_norm)
228       iter = iter + 1
229    end do
230
231    ! Deallocate auxiliary data structures
232    ! helped to assist with pmvm
233    call pmvmCommStructs_destroy()
234
235    deallocate(b_k)
236    deallocate(z)
237    deallocate(m)
238    deallocate(q)
239    deallocate(p)
240    deallocate(r)
241    if (associated(gv_tmp)) deallocate(gv_tmp)
242    if (associated(arr_copy)) deallocate(arr_copy)
243  end subroutine pcg
244
245  subroutine preconditioner_1level(sol,A,rhs,M,res,A_interf_,refactor_)
246    use subsolvers
247    implicit none
248    real(kind=rk),dimension(:),pointer :: sol !< solution
249    type(SpMtx)                        :: A   !< sparse system matrix
250    real(kind=rk),dimension(:),pointer :: rhs !< right hand side
251    type(Mesh),intent(in)              :: M   !< Mesh
252    real(kind=rk),dimension(:),pointer :: res !< residual vector, allocated
253                                              !! here for multiplicative Schwarz
254    type(SpMtx),optional               :: A_interf_  !< matr@interf.
255    logical,intent(inout),optional :: refactor_
256
257    type(SpMtx)                        :: A_tmp
258
259    if (refactor_.and.present(A_interf_)) then!{
260       if (sctls%verbose>9) then
261          !call SpMtx_printMat(A)
262          call SpMtx_printRaw(A)
263          !call SpMtx_printMat(A_interf_)
264          call SpMtx_printRaw(A_interf_)
265          !stop
266       endif
267       if (sctls%input_type==DCTL_INPUT_TYPE_ASSEMBLED) then
268          ! we need the fully sorted entries but A has block struct...
269          !   keep the original values...
270          if (associated(A%M_bound)) then
271             A_tmp=SpMtx_newInit(nnz=A%nnz,nblocks=A%nblocks,&
272                  nrows=A%nrows,ncols=A%ncols,&
273                  indi=A%indi,indj=A%indj,val=A%val,&
274                  arrange_type=A%arrange_type,M_bound=A%M_bound)
275          else
276             A_tmp=SpMtx_newInit(nnz=A%nnz,nblocks=A%nblocks,&
277                  nrows=A%nrows,ncols=A%ncols,&
278                  indi=A%indi,indj=A%indj,val=A%val,&
279                  arrange_type=A%arrange_type)
280          endif
281          A%arrange_type=D_SpMtx_ARRNG_NO
282          !A%nrows=max(A%nrows,A_interf_%nrows)
283          !A%ncols=max(A%nrows,A_interf_%ncols)
284          !if (associated(A%M_bound)) then
285          !  deallocate(A%M_bound)
286          !endif
287          call SpMtx_arrange(A,D_SpMtx_ARRNG_ROWS,sort=.true.)
288          call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
289               A_interf_=A_interf_, &
290               refactor=refactor_) !fine solves
291          if (sctls%method>1) then ! For multiplicative Schwarz method...:
292             refactor_=.false.
293             call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
294                  A_interf_=A_interf_, &
295                  refactor=refactor_) !fine solves
296          endif
297          ! put the original structure and orders back:
298          A%indi=A_tmp%indi
299          A%indj=A_tmp%indj
300          A%val=A_tmp%val
301          if (associated(A_tmp%M_bound)) A%M_bound=A_tmp%M_bound
302          A%nrows=A_tmp%nrows
303          A%ncols=A_tmp%ncols
304          A%arrange_type=A_tmp%arrange_type
305          call SpMtx_Destroy(A_tmp)
306       else
307          call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
308               A_interf_=A_interf_, &
309               refactor=refactor_) !fine solves
310          if (sctls%method>1) then ! For multiplicative Schwarz method...:
311             refactor_=.false.
312             call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
313                  A_interf_=A_interf_, &
314                  refactor=refactor_) !fine solves
315          endif
316       endif
317    elseif (refactor_.and.sctls%input_type==DCTL_INPUT_TYPE_ASSEMBLED) then!}{
318       if (associated(A%M_bound)) then
319          A_tmp=SpMtx_newInit(nnz=size(A%indi),nblocks=A%nblocks,&
320               nrows=A%nrows,ncols=A%ncols,&
321               indi=A%indi,indj=A%indj,val=A%val,&
322               arrange_type=A%arrange_type,M_bound=A%M_bound)
323       else
324          A_tmp=SpMtx_newInit(nnz=size(A%indi),nblocks=A%nblocks,&
325               nrows=A%nrows,ncols=A%ncols,&
326               indi=A%indi,indj=A%indj,val=A%val,&
327               arrange_type=A%arrange_type)
328       endif
329       A_tmp%nnz=A%nnz
330       A%arrange_type=D_SpMtx_ARRNG_NO
331       call SpMtx_arrange(A,D_SpMtx_ARRNG_ROWS,sort=.true.)
332       call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
333            refactor=refactor_) !fine solves
334       if (sctls%method>1) then ! For multiplicative Schwarz method...:
335          refactor_=.false.
336          call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
337               refactor=refactor_) !fine solves
338       endif
339       ! put the original structure and orders back:
340       A%indi=A_tmp%indi
341       A%indj=A_tmp%indj
342       A%val=A_tmp%val
343       if (associated(A_tmp%M_bound)) A%M_bound=A_tmp%M_bound
344       A%nrows=A_tmp%nrows
345       A%ncols=A_tmp%ncols
346       A%arrange_type=A_tmp%arrange_type
347       call SpMtx_Destroy(A_tmp)
348    else!}{
349       call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
350            A_interf_=A_interf_, &
351            refactor=refactor_) !fine solves
352       if (sctls%method>1) then ! For multiplicative Schwarz method...:
353          refactor_=.false.
354          call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
355               A_interf_=A_interf_, &
356               refactor=refactor_) !fine solves
357       endif
358    endif!}
359  end subroutine preconditioner_1level
360
361  subroutine prec2Level_exchangeMatrix(CoarseMtx_)
362    use CoarseAllgathers
363    type(SpMtx) :: CoarseMtx_ !< Coarse matrix
364    integer :: ol
365    logical :: add
366
367    ol=max(sctls%overlap,sctls%smoothers)
368    if (ol==0) then
369      add=.false.
370    else
371      add=.true.
372    endif
373
374    call AllSendCoarseMtx(cdat%LAC,CoarseMtx_,cdat%lg_cfmap,&
375         cdat%ngfc,cdat%nprocs,cdat%send)
376    call AllRecvCoarseMtx(CoarseMtx_,cdat%send,add=add) ! Recieve it
377
378  end subroutine prec2Level_exchangeMatrix
379
380  !-------------------------------
381  !> Make preconditioner
382  !-------------------------------
383  subroutine preconditioner(sol,A,rhs,M, &
384               A_interf_,CoarseMtx_,Restrict,refactor_,bugtrack_)
385    use subsolvers
386    use CoarseMtx_mod
387    use CoarseAllgathers
388    use Vect_mod
389    implicit none
390    real(kind=rk),dimension(:),pointer :: sol !< solution
391    type(SpMtx)                        :: A   !< sparse system matrix
392    real(kind=rk),dimension(:),pointer :: rhs !< right hand side
393    type(Mesh),intent(in)              :: M   !< Mesh
394    real(kind=rk),dimension(:),pointer :: res !< residual vector, allocated
395                                              !! here for multiplicative Schwarz
396    type(SpMtx),optional               :: A_interf_  !< matr@interf.
397    type(SpMtx),optional               :: CoarseMtx_ !< Coarse matrix
398    type(SpMtx),optional               :: Restrict   !< Restriction matrix
399    logical,intent(inout),optional :: refactor_
400    logical,optional                   :: bugtrack_
401    ! ----- local: ------
402    real(kind=rk),dimension(:),pointer,save :: csol,crhs,tmpsol,tmpsol2,clrhs
403    integer :: i,ol
404    logical :: add,bugtrack
405    real(kind=rk) :: t1
406    logical :: isFirstIter
407
408    t1 = MPI_WTime()
409
410    if (present(bugtrack_)) then
411      bugtrack=bugtrack_
412    else
413      bugtrack=.false.
414    endif
415    ! ----------------------------
416    isFirstIter = present(refactor_).and.refactor_
417    if (sctls%method==0) then
418      sol=rhs
419      return
420    endif
421    sol=0.0_rk
422    if (sctls%method>1) then ! For multiplicative Schwarz method...:
423      allocate(res(size(rhs)))
424    endif
425    if (present(CoarseMtx_)) then !{
426      if (.not.present(Restrict)) call DOUG_abort("Restriction matrix needs to be passed along with the coarse matrix!")
427
428      if (isFirstIter.and.cdat%active) then
429        ! First iteration - send matrix
430        call prec2Level_exchangeMatrix(CoarseMtx_)
431      end if
432
433      ! after coarse matrix is received allocate coarse vectors
434      if (isFirstIter.and.sctls%levels>1) then
435        if (associated(crhs)) then
436          if (size(crhs)/=CoarseMtx_%ncols) then
437            deallocate(csol)
438            deallocate(crhs)
439            allocate(crhs(CoarseMtx_%ncols))
440            allocate(csol(CoarseMtx_%nrows))
441          endif
442        else
443          allocate(crhs(CoarseMtx_%ncols))
444          allocate(csol(CoarseMtx_%nrows))
445        endif
446        allocate(clrhs(Restrict%nrows)) ! allocate memory for vector
447      end if
448      if (.not.associated(tmpsol)) then
449        !allocate(tmpsol(A%nrows))
450        allocate(tmpsol(size(rhs)))
451      endif
452
453      if (cdat%active) then
454        ! Send coarse vector
455        call SpMtx_Ax(clrhs,Restrict,rhs,dozero=.true.) ! restrict <RA>
456        if (cdat_vec%active) then
457          call AllSendCoarseVector(clrhs,cdat_vec%nprocs,cdat_vec%cdisps,&
458               cdat_vec%send,useprev=.not.isFirstIter)
459        else
460          call AllSendCoarseVector(clrhs,cdat%nprocs,cdat%cdisps,&
461               cdat%send,useprev=.not.isFirstIter)
462        endif
463      end if ! cdat%active
464
465      ! first level prec
466      if (sctls%input_type==DCTL_INPUT_TYPE_ELEMENTAL.or.numprocs>1) then
467        call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
468                        A_interf_=A_interf_, &
469                        refactor=refactor_) !fine solves
470      else
471        call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
472                        A_interf_=A_interf_,AC=CoarseMtx_, &
473                        refactor=refactor_,Restrict=Restrict) !fine solves
474      endif
475
476      if (sctls%levels>1) then
477        ol=max(sctls%overlap,sctls%smoothers)
478        if (isFirstIter.and..not.cdat%active) then ! 1 processor case
479          if (sctls%smoothers==-1) then
480            allocate(tmpsol2(A%nrows))
481            tmpsol2=0.0_rk
482            call exact_sparse_multismoother(tmpsol2,A,rhs)
483            call SpMtx_Ax(crhs,Restrict,tmpsol2,dozero=.true.) ! restriction
484          else
485            if (sctls%method>1.and.sctls%method/=5) then ! multiplicative Schwarz
486              call SpMtx_Ax(crhs,Restrict,res,dozero=.true.) ! restriction
487            else
488              call SpMtx_Ax(crhs,Restrict,rhs,dozero=.true.) ! restriction
489            endif
490          endif
491        endif
492      end if
493
494      !csol=0.0_rk
495      if (cdat%active) then
496        ! Recieve the vector for solve
497        if (cdat_vec%active) then
498          call AllRecvCoarseVector(crhs,cdat_vec%nprocs,&
499               cdat_vec%cdisps,cdat_vec%glg_cfmap,cdat_vec%send)
500        else
501          call AllRecvCoarseVector(crhs,cdat%nprocs,&
502               cdat%cdisps,cdat%glg_cfmap,cdat%send)
503        endif
504        !write(stream,*) "Got coarse vector!"
505        !call MPI_BARRIER(MPI_COMM_WORLD,i)
506      end if
507
508      if (sctls%levels>1) then
509        if (isFirstIter) then
510          write (stream,*) &
511               'factorising coarse matrix of size',CoarseMtx_%nrows, &
512               ' and nnz:',CoarseMtx_%nnz
513          call free_spmtx_subsolves(CoarseMtx_)
514          allocate(CoarseMtx_%subsolve_ids(1))
515          CoarseMtx_%subsolve_ids=0
516          CoarseMtx_%nsubsolves=1
517        end if
518
519        ! Coarse solve
520        call sparse_singlesolve(CoarseMtx_%subsolve_ids(1),csol,crhs,&
521             nfreds=CoarseMtx_%nrows, &
522             nnz=CoarseMtx_%nnz,        &
523             indi=CoarseMtx_%indi,      &
524             indj=CoarseMtx_%indj,      &
525             val=CoarseMtx_%val)
526        if (bugtrack)write(stream,*) "(f) Coarse SOL is:",csol
527
528        if (cdat_vec%active) then
529          call Vect_remap(csol,clrhs,cdat_vec%gl_cfmap,dozero=.true.)
530          call SpMtx_Ax(tmpsol,Restrict,clrhs,dozero=.true.,transp=.true.)
531        elseif (cdat%active) then
532          call Vect_remap(csol,clrhs,cdat%gl_cfmap,dozero=.true.)
533          call SpMtx_Ax(tmpsol,Restrict,clrhs,dozero=.true.,transp=.true.)
534        else
535          if (sctls%smoothers==-1) then
536            call SpMtx_Ax(tmpsol2,Restrict,csol,dozero=.true.,transp=.true.) ! interpolation
537            tmpsol=0.0_rk
538            call exact_sparse_multismoother(tmpsol,A,tmpsol2)
539          else
540            call SpMtx_Ax(tmpsol,Restrict,csol,dozero=.true.,transp=.true.) ! interpolation
541          endif
542        endif
543
544        if (.not.(cdat_vec%active)) then
545          if (sctls%smoothers==-1) then
546            tmpsol2=0.0_rk
547            call exact_sparse_multismoother(tmpsol2,A,rhs)
548            call SpMtx_Ax(crhs,Restrict,tmpsol2,dozero=.true.) ! restriction
549          else
550            call SpMtx_Ax(crhs,Restrict,rhs,dozero=.true.) ! restriction
551          endif
552        endif
553      endif !}
554      if (sctls%method==1) then
555!call Print_Glob_Vect(sol,M,'sol===',chk_endind=M%ninner)
556        sol(1:A%nrows)=sol(1:A%nrows)+tmpsol(1:A%nrows)
557      elseif (sctls%method==3) then ! fully multiplicative Schwarz
558        sol(1:A%nrows)=sol(1:A%nrows)+tmpsol(1:A%nrows)
559        ! calculate the residual:
560        call SpMtx_Ax(res,A,sol,dozero=.true.) !
561        res=rhs-res
562      elseif (sctls%method==2.and.ol>0) then ! fully multiplicative Schwarz
563        sol(1:A%nrows)=tmpsol(1:A%nrows)
564        ! calculate the residual:
565        call SpMtx_Ax(res,A,sol,dozero=.true.) !
566        res=rhs-res
567      endif
568      if (sctls%method>1) then ! For multiplicative Schwarz method...:
569        refactor_=.false.
570        if (sctls%input_type==DCTL_INPUT_TYPE_ELEMENTAL) then
571            call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
572                          A_interf_=A_interf_, &
573                          refactor=refactor_) !fine solves
574        else
575            call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
576                          A_interf_=A_interf_,AC=CoarseMtx_, &
577                          refactor=refactor_,Restrict=Restrict) !fine solves
578        endif
579      endif
580      if (((ol==0.and.sctls%method==2).or.sctls%method==5).and.sctls%levels>1) then
581        ! multiplicative on fine level, additive with coarse level:
582        sol(1:A%nrows)=sol(1:A%nrows)+tmpsol(1:A%nrows)
583      endif
584    else !}{
585       call preconditioner_1level(sol,A,rhs,M,res,A_interf_,refactor_)
586    endif !}
587    if (sctls%method>1) then ! For multiplicative Schwarz method...:
588      if (associated(res)) deallocate(res)
589    endif
590
591    time_preconditioner = time_preconditioner + MPI_WTime()-t1
592
593  end subroutine preconditioner
594
595  !--------------------------
596  ! Solve
597  !--------------------------
598  subroutine msolve (m, r, z)
599    implicit none
600    real(kind=rk), dimension(:),     intent(in) :: m, r
601    real(kind=rk), dimension(:), intent(in out) :: z
602    integer :: i
603
604    do i = 1, size(r)
605       z(i) = r(i) !/ m(i)
606    end do
607  end subroutine msolve
608
609  !--------------------------
610  !> Preconditioned conjugent gradient method with eigenvalues
611  !--------------------------
612  subroutine pcg_weigs (A,b,x,Msh,it,cond_num,A_interf_,tol_,maxit_, &
613       x0_,solinf,resvects_,CoarseMtx_,Restrict,refactor_)
614    use CoarseAllgathers
615
616    implicit none
617   
618    type(SpMtx),intent(in out)          :: A !< System matrix (sparse)
619    float(kind=rk),dimension(:),pointer :: b !< right hand side
620    float(kind=rk),dimension(:),pointer :: x !< Solution
621    type(Mesh),intent(in)               :: Msh !< Mesh - aux data for Ax operation
622
623    integer,intent(out) :: it
624    real(kind=rk),intent(out) :: cond_num
625   
626    ! Optional arguments
627    type(SpMtx),intent(in out),optional                :: A_interf_ !< matr@interf.
628    real(kind=rk), intent(in), optional                :: tol_ !< Tolerance of the method
629    integer,       intent(in), optional                :: maxit_ !< Max number of iterations
630    float(kind=rk), dimension(:), intent(in), optional :: x0_ !< Initial guess
631    type(ConvInf),            intent(in out), optional :: solinf !< Solution statistics
632    logical,                      intent(in), optional :: resvects_ !< Fill in the 'resvect' or not
633    type(SpMtx),optional                               :: CoarseMtx_ !< Coarse matrix
634    type(SpMtx),optional                               :: Restrict !< Restriction mtx
635    logical,intent(in),optional                        :: refactor_
636   
637    ! Local variables
638    logical :: refactor
639
640    real(kind=rk) :: tol   ! Tolerance
641    integer       :: maxit ! Max number of iterations
642    logical       :: resvects=.false.
643
644    real(kind=rk),dimension(:),pointer,save :: r, p, q, m, z, b_k
645    real(kind=rk),dimension(:),pointer,save :: ztmp !TODO remove me
646    real(kind=rk) :: rho_curr, rho_prev, init_norm, res_norm
647    real(kind=rk) :: tmp,ratio_norm
648    integer :: nfreds
649    integer :: ierr
650    real(kind=rk),dimension(:),pointer,save :: dd,ee,alpha,beta
651    !logical,parameter :: bugtrack=.true.
652    logical,parameter :: bugtrack=.false.
653    ! profiling
654    real(8) :: time_iterations, t1
655
656    write(stream,'(/a)') 'Preconditioned conjugate gradient:'
657    if (present(refactor_).and.refactor_) then
658      refactor=.true.
659    else
660      refactor=.false.
661    endif
662
663    if (size(b) /= size(x)) &
664         call DOUG_abort('[pcg] : SEVERE : size(b) /= size(x)',-1)
665
666    nfreds=size(b)
667    if (.not.associated(r)) then
668      allocate(r(nfreds))
669      allocate(p(nfreds))
670      allocate(q(nfreds))
671      allocate(m(nfreds))
672      allocate(z(nfreds))
673      allocate(b_k(nfreds))
674    endif
675
676    tol = sctls%solve_tolerance
677    if (present(tol_)) &
678         tol = tol_
679
680    maxit = sctls%solve_maxiters
681    if (maxit==-1) maxit=100000
682    if (present(maxit_)) &
683         maxit = maxit_
684
685    if (present(resvects_).and.(.not.present(solinf))) &
686         call DOUG_abort('[pcg] : SEVERE : "resvects_" must be given along'//&
687         ' with "solinf".',-1)
688
689    ! Allocate convergence info structure
690    if (present(solinf).and.(.not.present(resvects_))) then
691       call ConvInf_Init(solinf)
692    else if (present(solinf).and.&
693         (present(resvects_).and.(resvects_.eqv.(.true.)))) then
694       call ConvInf_Init(solinf, maxit)
695       resvects = .true.
696    end if
697
698
699    ! Initialise auxiliary data structures
700    ! to assist with pmvm
701    call pmvmCommStructs_init(A, Msh)
702
703
704if (bugtrack)call Print_Glob_Vect(x,Msh,'global x===')
705    call SpMtx_pmvm(r,A,x,Msh)
706
707    r = b - r
708    init_norm = Vect_dot_product(r,r)
709    if (init_norm == 0.0) &
710         init_norm = 1.0
711
712
713    write(stream,'(a,i6,a,e8.2)') 'maxit = ',maxit,', tol = ',tol
714    write(stream,'(a,e22.15)') 'init_norm = ', dsqrt(init_norm)
715
716    it = 0
717    ratio_norm = 1.0_rk
718    ! arrays for eigenvalue calculation:
719    if (.not.associated(dd)) then
720      allocate(dd(maxit))
721      allocate(ee(maxit))
722      allocate(alpha(maxit))
723      allocate(beta(maxit))
724    endif
725
726    time_iterations = 0.
727    ! iterations
728    do while((ratio_norm > tol*tol).and.(it < maxit))
729      it = it + 1
730      t1 = MPI_WTime()
731
732if (bugtrack)call Print_Glob_Vect(r,Msh,'global r===',chk_endind=Msh%ninner)
733      call preconditioner(sol=z,          &
734                            A=A,          &
735                          rhs=r,          &
736                            M=Msh,        &
737                    A_interf_=A_interf_,  &
738                   CoarseMtx_=CoarseMtx_, &
739                    Restrict=Restrict,    &
740                    refactor_=refactor,   &
741                    bugtrack_=bugtrack)
742
743      refactor=.false.
744if (bugtrack)call Print_Glob_Vect(z,Msh,'global bef comm z===',chk_endind=Msh%ninner)
745      if (sctls%method/=0) then
746        call Add_common_interf(z,A,Msh)
747      endif
748!call Print_Glob_Vect(z,Msh,'global aft comm z===',chk_endind=Msh%ninner)
749!call Print_Glob_Vect(z,Msh,'global z===')
750      ! compute current rho
751      rho_curr = Vect_dot_product(r,z)
752
753      if (it == 1) then
754         p = z
755      else
756         beta(it) = rho_curr / rho_prev
757         p = z + beta(it) * p
758      end if
759      call SpMtx_pmvm(q,A, p, Msh)
760if (bugtrack)call Print_Glob_Vect(q,Msh,'global q===')
761      ! compute alpha
762      alpha(it) = rho_curr / Vect_dot_product(p,q)
763      x = x + alpha(it) * p
764if (bugtrack)call Print_Glob_Vect(x,Msh,'global x===')
765      r = r - alpha(it) * q
766if (bugtrack)call Print_Glob_Vect(r,Msh,'global r===')
767      rho_prev = rho_curr
768      ! check
769      res_norm = Vect_dot_product(r,r)
770      !write (stream,*) "Norm is", res_norm
771      ratio_norm = res_norm / init_norm
772
773      time_iterations = time_iterations + MPI_WTime()-t1
774
775      if (ismaster()) &
776           write(stream, '(i5,a,e22.15)') it,': res_norm=',dsqrt(res_norm)
777    end do
778
779    if(pstream/=0) then
780       write(pstream, "(I0,':pcg iterations:',I0)") myrank, it
781       write(pstream, "(I0,':iterations time:',F0.3)") myrank, time_iterations
782       write(pstream, "(I0,':preconditioner time:',F0.3)") myrank, time_preconditioner
783    end if
784
785    if (ismaster()) then
786      call CalculateEigenvalues(it,dd,ee,alpha,beta,maxit)
787      cond_num=dd(it)/dd(1)
788      write(stream,'(a,i3,a,e10.4)') '#it:',it,' Cond#: ',cond_num
789    endif
790    !deallocate(beta)
791    !deallocate(alpha)
792    !deallocate(ee)
793    !deallocate(dd)
794    ! Deallocate auxiliary data structures
795    ! helped to assist with pmvm
796 !  call pmvmCommStructs_destroy()
797
798    !deallocate(b_k)
799    !deallocate(z)
800    !deallocate(m)
801    !deallocate(q)
802    !deallocate(p)
803    !deallocate(r)
804  end subroutine pcg_weigs
805
806
807  subroutine CalculateEigenvalues(it,dd,ee,alpha,beta,MaxIt)
808    ! Finds eigenvalues using Lanczos connection
809    implicit none
810    integer :: it,MaxIt
811    real(kind=rk),dimension(:),pointer :: dd,ee,alpha,beta
812    integer :: i,j,ierr
813
814    dd(1) = 1.0_rk/alpha(1)
815    ee(1) = 0.0_rk
816    do i=2,it
817      if (alpha(i)==0.or.beta(i)<0) then
818        print *,'Unable to compute eigenvalue number',i,' and '
819        do j=1,i
820          print *,j,' alpha:',alpha(j),' beta:',beta
821        enddo
822        return
823      else
824        dd(i) = 1.0_rk/alpha(i) + beta(i)/alpha(i-1)
825        ee(i) = -dsqrt(beta(i))/alpha(i-1)
826      endif
827    enddo
828    !call tql1(it,dd,ee,MaxIt,ierr)
829    call tql1(it,dd,ee,ierr)
830    if (ierr > 0) then
831      print *,'Cancelled after EigMaxIt while calculating eig',ierr
832    endif
833  end subroutine CalculateEigenvalues
834
835  subroutine tql1(n,d,e,ierr)
836    implicit none
837    integer :: i,j,l,m,n,ii,l1,l2,mml,ierr
838    real(kind=rk) :: d(n),e(n)
839    real(kind=rk) :: c,c2,c3,dl1,el1,f,g,h,p,r,s,s2,tst1,tst2 !,pythag
840!
841!       this subroutine is a translation of the algol procedure tql1,
842!       num. math. 11, 293-306(1968) by bowdler, martin, reinsch, and
843!       wilkinson.
844!       handbook for auto. comp., vol.ii-linear algebra, 227-240(1971).
845!
846!       this subroutine finds the eigenvalues of a symmetric
847!       tridiagonal matrix by the ql method.
848!
849!       on input
850!          n is the order of the matrix.
851!          d contains the diagonal elements of the input matrix.
852!          e contains the subdiagonal elements of the input matrix
853!            in its last n-1 positions.  e(1) is arbitrary.
854!
855!        on output
856!          d contains the eigenvalues in ascending order.  if an
857!            error exit is made, the eigenvalues are correct and
858!            ordered for indices 1,2,...ierr-1, but may not be
859!            the smallest eigenvalues.
860!          e has been destroyed.
861!          ierr is set to
862!            zero       for normal return,
863!            j          if the j-th eigenvalue has not been
864!                       determined after 30 iterations.
865!       calls pythag for  sqrt(a*a + b*b) .
866!
867!       questions and comments should be directed to burton s. garbow,
868!       mathematics and computer science div, argonne national laboratory
869!
870!       this version dated august 1983.
871!
872!       ------------------------------------------------------------------
873!
874        ierr = 0
875        if (n .eq. 1) go to 1001
876!
877        do 100 i = 2, n
878    100 e(i-1) = e(i)
879!
880        f = 0.0e0
881        tst1 = 0.0e0
882        e(n) = 0.0e0
883
884        do 290 l = 1, n
885           j = 0
886           h = abs(d(l)) + abs(e(l))
887           if (tst1 .lt. h) tst1 = h
888!       .......... look for small sub-diagonal element ..........
889           do 110 m = l, n
890              tst2 = tst1 + abs(e(m))
891              if (tst2 .eq. tst1) go to 120
892!       .......... e(n) is always zero, so there is no exit
893!                  throu     gh the bottom of the loop ..........
894    110    continue
895
896    120    if (m .eq. l)      go to 210
897    130    if (j .eq. 30     ) go to 1000
898           j = j + 1
899!       .......... form shift ..........
900           l1 = l + 1
901           l2 = l1 + 1
902           g = d(l)
903           p = (d(l1) - g) / (2.0e0 * e(l))
904           r = pythag(p,1.0e0_rk)
905           d(l) = e(l) / (p + sign(r,p))
906           d(l1) = e(l) * (p + sign(r,p))
907           dl1 = d(l1)
908           h = g - d(l)
909           if (l2 .gt. n) go to 145
910
911           do 140 i = l2, n
912    140    d(i) = d(i) - h
913
914    145    f = f + h
915!       .......... ql transformation ..........
916           p = d(m)
917           c = 1.0e0
918           c2 = c
919           el1 = e(l1)
920           s = 0.0e0
921           mml = m - l
922!       .......... for i=m-1 step -1 until l do -- ..........
923           do 200 ii = 1, mml
924              c3 = c2
925              c2 = c
926              s2 = s
927              i = m - ii
928              g = c * e(i)
929              h = c * p
930              r = pythag(p,e(i))
931              e(i+1) = s * r
932              s = e(i) / r
933              c = p / r
934              p = c * d(i) - s * g
935              d(i+1) = h + s * (c * g + s * d(i))
936    200    continue
937
938           p = -s * s2 * c3 * el1 * e(l) / dl1
939           e(l) = s * p
940           d(l) = c * p
941           tst2 = tst1 + abs(e(l))
942           if (tst2 .gt. tst1) go to 130
943    210    p = d(l) + f
944!       .......... order eigenvalues ..........
945           if (l .eq. 1) go to 250
946!       .......... for i=l step -1 until 2 do -- ..........
947           do 230 ii = 2, l
948              i = l + 2 - ii
949              if (p .ge. d(i-1)) go to 270
950              d(i) = d(i-1)
951    230    continue
952
953    250    i = 1
954    270    d(i) = p
955    290 continue
956
957        go to 1001
958!       .......... set error -- no convergence to an
959!                  eigenvalue after 30 iterations ..........
960   1000 ierr = l
961   1001 return
962  end subroutine tql1
963
964  function pythag(a,b) result(pyt)
965    ! finds dsqrt(a**2+b**2) without overflow or destructive underflow
966    implicit none
967    real(kind=rk) :: a,b,pyt
968    real(kind=rk) :: p,r,s,t,u
969
970    p = dmax1(dabs(a),dabs(b))
971    if (p .eq. 0.0d0) go to 20
972    r = (dmin1(dabs(a),dabs(b))/p)**2
973 10 continue
974    t = 4.0d0 + r
975    if (t .eq. 4.0d0) go to 20
976    s = r/t
977    u = 1.0d0 + 2.0d0*s
978    p = u*p
979    r = (s/u)**2 * r
980    go to 10
981 20 pyt = p
982    return
983  end function pythag
984
985
986end module pcg_mod
987
988
Note: See TracBrowser for help on using the repository browser.