source: src/solvers/pcg.F90 @ 82752f1

Revision 82752f1, 33.3 KB checked in by Oleg Batrashev <ogbash@…>, 2 years ago (diff)

Another fix for coarse preconditioners.

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