source: src/solvers/pcg.F90 @ dc3a733

Revision dc3a733, 29.5 KB checked in by Oleg Batrashev <ogbash@…>, 2 years ago (diff)

Combine first level and second level preconditioner codes.

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