source: src/solvers/pcg.F90 @ 8a2eebe

Revision 8a2eebe, 33.4 KB checked in by Oleg Batrashev <ogbash@…>, 2 years ago (diff)

Fix some tests with 1 node two-level.

  • 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
432        ! factorise coarse matrix
433        write (stream,*) &
434             'factorising coarse matrix of size',CoarseMtx_%nrows, &
435             ' and nnz:',CoarseMtx_%nnz
436        call free_spmtx_subsolves(CoarseMtx_)
437        allocate(CoarseMtx_%subsolve_ids(1))
438        CoarseMtx_%subsolve_ids=0
439        CoarseMtx_%nsubsolves=1
440
441        ! Factorise the matrix
442        if (sctls%verbose>2) then
443          write(stream,*)'Global coarse matrix is:---------'
444          do i=1,CoarseMtx_%nnz
445             write(stream,*) CoarseMtx_%indi(i),&
446                  CoarseMtx_%indj(i),CoarseMtx_%val(i)
447          enddo
448          write(stream,*)'---------------------------------'
449        endif
450        call factorise(CoarseMtx_%subsolve_ids(1), &
451             nfreds=CoarseMtx_%nrows,   &
452             nnz=CoarseMtx_%nnz,        &
453             indi=CoarseMtx_%indi,      &
454             indj=CoarseMtx_%indj,      &
455             val=CoarseMtx_%val)
456        CoarseMtx_%indi=CoarseMtx_%indi+1
457        CoarseMtx_%indj=CoarseMtx_%indj+1
458     end if
459
460     ! after coarse matrix is received allocate coarse vectors
461     if (isFirstIter.and.sctls%levels>1) then
462       if (associated(crhs)) then
463         if (size(crhs)/=CoarseMtx_%ncols) then
464           deallocate(csol)
465           deallocate(crhs)
466           allocate(crhs(CoarseMtx_%ncols))
467           allocate(csol(CoarseMtx_%nrows))
468         endif
469       else
470         allocate(crhs(CoarseMtx_%ncols))
471         allocate(csol(CoarseMtx_%nrows))
472       endif
473       allocate(clrhs(Restrict%nrows)) ! allocate memory for vector
474     end if
475     if (.not.associated(tmpsol)) then
476       !allocate(tmpsol(A%nrows))
477       allocate(tmpsol(size(rhs)))
478     endif
479
480     if (cdat%active) then
481        ! Send coarse vector
482        call SpMtx_Ax(clrhs,Restrict,rhs,dozero=.true.) ! restrict <RA>
483        if (cdat_vec%active) then
484          call AllSendCoarseVector(clrhs,cdat_vec%nprocs,cdat_vec%cdisps,&
485               cdat_vec%send,useprev=.not.isFirstIter)
486        else
487          call AllSendCoarseVector(clrhs,cdat%nprocs,cdat%cdisps,&
488               cdat%send,useprev=.not.isFirstIter)
489        endif
490      end if ! cdat%active
491
492      ! first level prec
493      if (sctls%input_type==DCTL_INPUT_TYPE_ELEMENTAL.or.numprocs>1) then
494        call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
495                        A_interf_=A_interf_, &
496                        refactor=refactor_) !fine solves
497      else
498        call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
499                        A_interf_=A_interf_,AC=CoarseMtx_, &
500                        refactor=refactor_,Restrict=Restrict) !fine solves
501      endif
502
503      if (sctls%levels>1) then
504        ol=max(sctls%overlap,sctls%smoothers)
505        if (isFirstIter) then
506          if (.not.cdat%active) then ! 1 processor case
507            if (sctls%smoothers==-1) then
508              allocate(tmpsol2(A%nrows))
509              tmpsol2=0.0_rk
510              call exact_sparse_multismoother(tmpsol2,A,rhs)
511              call SpMtx_Ax(crhs,Restrict,tmpsol2,dozero=.true.) ! restriction
512            else
513              if (sctls%method>1.and.sctls%method/=5) then ! multiplicative Schwarz
514                call SpMtx_Ax(crhs,Restrict,res,dozero=.true.) ! restriction
515              else
516                call SpMtx_Ax(crhs,Restrict,rhs,dozero=.true.) ! restriction
517              endif
518            endif
519          endif
520        end if
521      end if
522
523        !csol=0.0_rk
524      if (cdat%active) then
525        ! Recieve the vector for solve
526        if (cdat_vec%active) then
527          call AllRecvCoarseVector(crhs,cdat_vec%nprocs,&
528               cdat_vec%cdisps,cdat_vec%glg_cfmap,cdat_vec%send)
529        else
530          call AllRecvCoarseVector(crhs,cdat%nprocs,&
531               cdat%cdisps,cdat%glg_cfmap,cdat%send)
532        endif
533        !write(stream,*) "Got coarse vector!"
534        !call MPI_BARRIER(MPI_COMM_WORLD,i)
535
536        call sparse_singlesolve(CoarseMtx_%subsolve_ids(1),csol,crhs,nfreds=CoarseMtx_%nrows)
537        write (stream,*) 'coarse factorisation done!',CoarseMtx_%subsolve_ids(1)
538        if (bugtrack)write(stream,*) "(f) Coarse SOL is:",csol
539      end if
540
541      if (sctls%levels>1) then
542        if (cdat_vec%active) then
543          call Vect_remap(csol,clrhs,cdat_vec%gl_cfmap,dozero=.true.)
544          call SpMtx_Ax(tmpsol,Restrict,clrhs,dozero=.true.,transp=.true.)
545        elseif (cdat%active) then
546          call Vect_remap(csol,clrhs,cdat%gl_cfmap,dozero=.true.)
547          call SpMtx_Ax(tmpsol,Restrict,clrhs,dozero=.true.,transp=.true.)
548        else
549          if (sctls%smoothers==-1) then
550            call SpMtx_Ax(tmpsol2,Restrict,csol,dozero=.true.,transp=.true.) ! interpolation
551            tmpsol=0.0_rk
552            call exact_sparse_multismoother(tmpsol,A,tmpsol2)
553          else
554            call SpMtx_Ax(tmpsol,Restrict,csol,dozero=.true.,transp=.true.) ! interpolation
555          endif
556        endif
557
558        if (.not.(cdat_vec%active)) then
559          if (sctls%smoothers==-1) then
560            tmpsol2=0.0_rk
561            call exact_sparse_multismoother(tmpsol2,A,rhs)
562            call SpMtx_Ax(crhs,Restrict,tmpsol2,dozero=.true.) ! restriction
563          else
564            call SpMtx_Ax(crhs,Restrict,rhs,dozero=.true.) ! restriction
565          endif
566        endif
567      endif !}
568      if (sctls%method==1) then
569!call Print_Glob_Vect(sol,M,'sol===',chk_endind=M%ninner)
570        sol(1:A%nrows)=sol(1:A%nrows)+tmpsol(1:A%nrows)
571      elseif (sctls%method==3) then ! fully multiplicative Schwarz
572        sol(1:A%nrows)=sol(1:A%nrows)+tmpsol(1:A%nrows)
573        ! calculate the residual:
574        call SpMtx_Ax(res,A,sol,dozero=.true.) !
575        res=rhs-res
576      elseif (sctls%method==2.and.ol>0) then ! fully multiplicative Schwarz
577        sol(1:A%nrows)=tmpsol(1:A%nrows)
578        ! calculate the residual:
579        call SpMtx_Ax(res,A,sol,dozero=.true.) !
580        res=rhs-res
581      endif
582      if (sctls%method>1) then ! For multiplicative Schwarz method...:
583        refactor_=.false.
584        if (sctls%input_type==DCTL_INPUT_TYPE_ELEMENTAL) then
585            call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
586                          A_interf_=A_interf_, &
587                          refactor=refactor_) !fine solves
588        else
589            call sparse_multisolve(sol=sol,A=A,M=M,rhs=rhs,res=res, &
590                          A_interf_=A_interf_,AC=CoarseMtx_, &
591                          refactor=refactor_,Restrict=Restrict) !fine solves
592        endif
593      endif
594      if (((ol==0.and.sctls%method==2).or.sctls%method==5).and.sctls%levels>1) then
595        ! multiplicative on fine level, additive with coarse level:
596        sol(1:A%nrows)=sol(1:A%nrows)+tmpsol(1:A%nrows)
597      endif
598    else !}{
599       call preconditioner_1level(sol,A,rhs,M,res,A_interf_,refactor_)
600    endif !}
601    if (sctls%method>1) then ! For multiplicative Schwarz method...:
602      if (associated(res)) deallocate(res)
603    endif
604
605    time_preconditioner = time_preconditioner + MPI_WTime()-t1
606
607  end subroutine preconditioner
608
609  !--------------------------
610  ! Solve
611  !--------------------------
612  subroutine msolve (m, r, z)
613    implicit none
614    real(kind=rk), dimension(:),     intent(in) :: m, r
615    real(kind=rk), dimension(:), intent(in out) :: z
616    integer :: i
617
618    do i = 1, size(r)
619       z(i) = r(i) !/ m(i)
620    end do
621  end subroutine msolve
622
623  !--------------------------
624  !> Preconditioned conjugent gradient method with eigenvalues
625  !--------------------------
626  subroutine pcg_weigs (A,b,x,Msh,it,cond_num,A_interf_,tol_,maxit_, &
627       x0_,solinf,resvects_,CoarseMtx_,Restrict,refactor_)
628    use CoarseAllgathers
629
630    implicit none
631   
632    type(SpMtx),intent(in out)          :: A !< System matrix (sparse)
633    float(kind=rk),dimension(:),pointer :: b !< right hand side
634    float(kind=rk),dimension(:),pointer :: x !< Solution
635    type(Mesh),intent(in)               :: Msh !< Mesh - aux data for Ax operation
636
637    integer,intent(out) :: it
638    real(kind=rk),intent(out) :: cond_num
639   
640    ! Optional arguments
641    type(SpMtx),intent(in out),optional                :: A_interf_ !< matr@interf.
642    real(kind=rk), intent(in), optional                :: tol_ !< Tolerance of the method
643    integer,       intent(in), optional                :: maxit_ !< Max number of iterations
644    float(kind=rk), dimension(:), intent(in), optional :: x0_ !< Initial guess
645    type(ConvInf),            intent(in out), optional :: solinf !< Solution statistics
646    logical,                      intent(in), optional :: resvects_ !< Fill in the 'resvect' or not
647    type(SpMtx),optional                               :: CoarseMtx_ !< Coarse matrix
648    type(SpMtx),optional                               :: Restrict !< Restriction mtx
649    logical,intent(in),optional                        :: refactor_
650   
651    ! Local variables
652    logical :: refactor
653
654    real(kind=rk) :: tol   ! Tolerance
655    integer       :: maxit ! Max number of iterations
656    logical       :: resvects=.false.
657
658    real(kind=rk),dimension(:),pointer,save :: r, p, q, m, z, b_k
659    real(kind=rk),dimension(:),pointer,save :: ztmp !TODO remove me
660    real(kind=rk) :: rho_curr, rho_prev, init_norm, res_norm
661    real(kind=rk) :: tmp,ratio_norm
662    integer :: nfreds
663    integer :: ierr
664    real(kind=rk),dimension(:),pointer,save :: dd,ee,alpha,beta
665    !logical,parameter :: bugtrack=.true.
666    logical,parameter :: bugtrack=.false.
667    ! profiling
668    real(8) :: time_iterations, t1
669
670    write(stream,'(/a)') 'Preconditioned conjugate gradient:'
671    if (present(refactor_).and.refactor_) then
672      refactor=.true.
673    else
674      refactor=.false.
675    endif
676
677    if (size(b) /= size(x)) &
678         call DOUG_abort('[pcg] : SEVERE : size(b) /= size(x)',-1)
679
680    nfreds=size(b)
681    if (.not.associated(r)) then
682      allocate(r(nfreds))
683      allocate(p(nfreds))
684      allocate(q(nfreds))
685      allocate(m(nfreds))
686      allocate(z(nfreds))
687      allocate(b_k(nfreds))
688    endif
689
690    tol = sctls%solve_tolerance
691    if (present(tol_)) &
692         tol = tol_
693
694    maxit = sctls%solve_maxiters
695    if (maxit==-1) maxit=100000
696    if (present(maxit_)) &
697         maxit = maxit_
698
699    if (present(resvects_).and.(.not.present(solinf))) &
700         call DOUG_abort('[pcg] : SEVERE : "resvects_" must be given along'//&
701         ' with "solinf".',-1)
702
703    ! Allocate convergence info structure
704    if (present(solinf).and.(.not.present(resvects_))) then
705       call ConvInf_Init(solinf)
706    else if (present(solinf).and.&
707         (present(resvects_).and.(resvects_.eqv.(.true.)))) then
708       call ConvInf_Init(solinf, maxit)
709       resvects = .true.
710    end if
711
712
713    ! Initialise auxiliary data structures
714    ! to assist with pmvm
715    call pmvmCommStructs_init(A, Msh)
716
717
718if (bugtrack)call Print_Glob_Vect(x,Msh,'global x===')
719    call SpMtx_pmvm(r,A,x,Msh)
720
721    r = b - r
722    init_norm = Vect_dot_product(r,r)
723    if (init_norm == 0.0) &
724         init_norm = 1.0
725
726
727    write(stream,'(a,i6,a,e8.2)') 'maxit = ',maxit,', tol = ',tol
728    write(stream,'(a,e22.15)') 'init_norm = ', dsqrt(init_norm)
729
730    it = 0
731    ratio_norm = 1.0_rk
732    ! arrays for eigenvalue calculation:
733    if (.not.associated(dd)) then
734      allocate(dd(maxit))
735      allocate(ee(maxit))
736      allocate(alpha(maxit))
737      allocate(beta(maxit))
738    endif
739
740    time_iterations = 0.
741    ! iterations
742    do while((ratio_norm > tol*tol).and.(it < maxit))
743      it = it + 1
744      t1 = MPI_WTime()
745
746if (bugtrack)call Print_Glob_Vect(r,Msh,'global r===',chk_endind=Msh%ninner)
747      call preconditioner(sol=z,          &
748                            A=A,          &
749                          rhs=r,          &
750                            M=Msh,        &
751                    A_interf_=A_interf_,  &
752                   CoarseMtx_=CoarseMtx_, &
753                    Restrict=Restrict,    &
754                    refactor_=refactor,   &
755                    bugtrack_=bugtrack)
756
757      refactor=.false.
758if (bugtrack)call Print_Glob_Vect(z,Msh,'global bef comm z===',chk_endind=Msh%ninner)
759      if (sctls%method/=0) then
760        call Add_common_interf(z,A,Msh)
761      endif
762!call Print_Glob_Vect(z,Msh,'global aft comm z===',chk_endind=Msh%ninner)
763!call Print_Glob_Vect(z,Msh,'global z===')
764      ! compute current rho
765      rho_curr = Vect_dot_product(r,z)
766
767      if (it == 1) then
768         p = z
769      else
770         beta(it) = rho_curr / rho_prev
771         p = z + beta(it) * p
772      end if
773      call SpMtx_pmvm(q,A, p, Msh)
774if (bugtrack)call Print_Glob_Vect(q,Msh,'global q===')
775      ! compute alpha
776      alpha(it) = rho_curr / Vect_dot_product(p,q)
777      x = x + alpha(it) * p
778if (bugtrack)call Print_Glob_Vect(x,Msh,'global x===')
779      r = r - alpha(it) * q
780if (bugtrack)call Print_Glob_Vect(r,Msh,'global r===')
781      rho_prev = rho_curr
782      ! check
783      res_norm = Vect_dot_product(r,r)
784      !write (stream,*) "Norm is", res_norm
785      ratio_norm = res_norm / init_norm
786
787      time_iterations = time_iterations + MPI_WTime()-t1
788
789      if (ismaster()) &
790           write(stream, '(i5,a,e22.15)') it,': res_norm=',dsqrt(res_norm)
791    end do
792
793    if(pstream/=0) then
794       write(pstream, "(I0,':pcg iterations:',I0)") myrank, it
795       write(pstream, "(I0,':iterations time:',F0.3)") myrank, time_iterations
796       write(pstream, "(I0,':preconditioner time:',F0.3)") myrank, time_preconditioner
797    end if
798
799    if (ismaster()) then
800      call CalculateEigenvalues(it,dd,ee,alpha,beta,maxit)
801      cond_num=dd(it)/dd(1)
802      write(stream,'(a,i3,a,e10.4)') '#it:',it,' Cond#: ',cond_num
803    endif
804    !deallocate(beta)
805    !deallocate(alpha)
806    !deallocate(ee)
807    !deallocate(dd)
808    ! Deallocate auxiliary data structures
809    ! helped to assist with pmvm
810 !  call pmvmCommStructs_destroy()
811
812    !deallocate(b_k)
813    !deallocate(z)
814    !deallocate(m)
815    !deallocate(q)
816    !deallocate(p)
817    !deallocate(r)
818  end subroutine pcg_weigs
819
820
821  subroutine CalculateEigenvalues(it,dd,ee,alpha,beta,MaxIt)
822    ! Finds eigenvalues using Lanczos connection
823    implicit none
824    integer :: it,MaxIt
825    real(kind=rk),dimension(:),pointer :: dd,ee,alpha,beta
826    integer :: i,j,ierr
827
828    dd(1) = 1.0_rk/alpha(1)
829    ee(1) = 0.0_rk
830    do i=2,it
831      if (alpha(i)==0.or.beta(i)<0) then
832        print *,'Unable to compute eigenvalue number',i,' and '
833        do j=1,i
834          print *,j,' alpha:',alpha(j),' beta:',beta
835        enddo
836        return
837      else
838        dd(i) = 1.0_rk/alpha(i) + beta(i)/alpha(i-1)
839        ee(i) = -dsqrt(beta(i))/alpha(i-1)
840      endif
841    enddo
842    !call tql1(it,dd,ee,MaxIt,ierr)
843    call tql1(it,dd,ee,ierr)
844    if (ierr > 0) then
845      print *,'Cancelled after EigMaxIt while calculating eig',ierr
846    endif
847  end subroutine CalculateEigenvalues
848
849  subroutine tql1(n,d,e,ierr)
850    implicit none
851    integer :: i,j,l,m,n,ii,l1,l2,mml,ierr
852    real(kind=rk) :: d(n),e(n)
853    real(kind=rk) :: c,c2,c3,dl1,el1,f,g,h,p,r,s,s2,tst1,tst2 !,pythag
854!
855!       this subroutine is a translation of the algol procedure tql1,
856!       num. math. 11, 293-306(1968) by bowdler, martin, reinsch, and
857!       wilkinson.
858!       handbook for auto. comp., vol.ii-linear algebra, 227-240(1971).
859!
860!       this subroutine finds the eigenvalues of a symmetric
861!       tridiagonal matrix by the ql method.
862!
863!       on input
864!          n is the order of the matrix.
865!          d contains the diagonal elements of the input matrix.
866!          e contains the subdiagonal elements of the input matrix
867!            in its last n-1 positions.  e(1) is arbitrary.
868!
869!        on output
870!          d contains the eigenvalues in ascending order.  if an
871!            error exit is made, the eigenvalues are correct and
872!            ordered for indices 1,2,...ierr-1, but may not be
873!            the smallest eigenvalues.
874!          e has been destroyed.
875!          ierr is set to
876!            zero       for normal return,
877!            j          if the j-th eigenvalue has not been
878!                       determined after 30 iterations.
879!       calls pythag for  sqrt(a*a + b*b) .
880!
881!       questions and comments should be directed to burton s. garbow,
882!       mathematics and computer science div, argonne national laboratory
883!
884!       this version dated august 1983.
885!
886!       ------------------------------------------------------------------
887!
888        ierr = 0
889        if (n .eq. 1) go to 1001
890!
891        do 100 i = 2, n
892    100 e(i-1) = e(i)
893!
894        f = 0.0e0
895        tst1 = 0.0e0
896        e(n) = 0.0e0
897
898        do 290 l = 1, n
899           j = 0
900           h = abs(d(l)) + abs(e(l))
901           if (tst1 .lt. h) tst1 = h
902!       .......... look for small sub-diagonal element ..........
903           do 110 m = l, n
904              tst2 = tst1 + abs(e(m))
905              if (tst2 .eq. tst1) go to 120
906!       .......... e(n) is always zero, so there is no exit
907!                  throu     gh the bottom of the loop ..........
908    110    continue
909
910    120    if (m .eq. l)      go to 210
911    130    if (j .eq. 30     ) go to 1000
912           j = j + 1
913!       .......... form shift ..........
914           l1 = l + 1
915           l2 = l1 + 1
916           g = d(l)
917           p = (d(l1) - g) / (2.0e0 * e(l))
918           r = pythag(p,1.0e0_rk)
919           d(l) = e(l) / (p + sign(r,p))
920           d(l1) = e(l) * (p + sign(r,p))
921           dl1 = d(l1)
922           h = g - d(l)
923           if (l2 .gt. n) go to 145
924
925           do 140 i = l2, n
926    140    d(i) = d(i) - h
927
928    145    f = f + h
929!       .......... ql transformation ..........
930           p = d(m)
931           c = 1.0e0
932           c2 = c
933           el1 = e(l1)
934           s = 0.0e0
935           mml = m - l
936!       .......... for i=m-1 step -1 until l do -- ..........
937           do 200 ii = 1, mml
938              c3 = c2
939              c2 = c
940              s2 = s
941              i = m - ii
942              g = c * e(i)
943              h = c * p
944              r = pythag(p,e(i))
945              e(i+1) = s * r
946              s = e(i) / r
947              c = p / r
948              p = c * d(i) - s * g
949              d(i+1) = h + s * (c * g + s * d(i))
950    200    continue
951
952           p = -s * s2 * c3 * el1 * e(l) / dl1
953           e(l) = s * p
954           d(l) = c * p
955           tst2 = tst1 + abs(e(l))
956           if (tst2 .gt. tst1) go to 130
957    210    p = d(l) + f
958!       .......... order eigenvalues ..........
959           if (l .eq. 1) go to 250
960!       .......... for i=l step -1 until 2 do -- ..........
961           do 230 ii = 2, l
962              i = l + 2 - ii
963              if (p .ge. d(i-1)) go to 270
964              d(i) = d(i-1)
965    230    continue
966
967    250    i = 1
968    270    d(i) = p
969    290 continue
970
971        go to 1001
972!       .......... set error -- no convergence to an
973!                  eigenvalue after 30 iterations ..........
974   1000 ierr = l
975   1001 return
976  end subroutine tql1
977
978  function pythag(a,b) result(pyt)
979    ! finds dsqrt(a**2+b**2) without overflow or destructive underflow
980    implicit none
981    real(kind=rk) :: a,b,pyt
982    real(kind=rk) :: p,r,s,t,u
983
984    p = dmax1(dabs(a),dabs(b))
985    if (p .eq. 0.0d0) go to 20
986    r = (dmin1(dabs(a),dabs(b))/p)**2
987 10 continue
988    t = 4.0d0 + r
989    if (t .eq. 4.0d0) go to 20
990    s = r/t
991    u = 1.0d0 + 2.0d0*s
992    p = u*p
993    r = (s/u)**2 * r
994    go to 10
995 20 pyt = p
996    return
997  end function pythag
998
999
1000end module pcg_mod
1001
1002
Note: See TracBrowser for help on using the repository browser.