xref: /phasta/phSolver/AMG/ramg_cheby.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen!***********************************************************
2*59599516SKenneth E. Jansen!      ramg_cheby_apply: u = Cheby(v,(lhs,r))
3*59599516SKenneth E. Jansen!      v : approx. solution before
4*59599516SKenneth E. Jansen!      u : approx. solution (after)
5*59599516SKenneth E. Jansen!      r : residual vector, r=Ax-b
6*59599516SKenneth E. Jansen!***********************************************************
7*59599516SKenneth E. Jansen
8*59599516SKenneth E. Jansen
9*59599516SKenneth E. Jansen      subroutine ramg_cheby_apply(u,v,r,level,colm,rowp,lhsK,lhsP,
10*59599516SKenneth E. Jansen     &                          ilwork,BC,iBC,iper,fwdbck)
11*59599516SKenneth E. Jansen      use ramg_data
12*59599516SKenneth E. Jansen      include "common.h"
13*59599516SKenneth E. Jansen      include "mpif.h"
14*59599516SKenneth E. Jansen      include "auxmpi.h"
15*59599516SKenneth E. Jansen
16*59599516SKenneth E. Jansen      integer,intent(in),dimension(nlwork)            :: ilwork
17*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg)              :: iBC,iper
18*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg,ndofBC)  :: BC
19*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg+1)       :: colm
20*59599516SKenneth E. Jansen      integer,intent(in),dimension(nnz_tot)      :: rowp
21*59599516SKenneth E. Jansen
22*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(9,nnz_tot)    :: lhsK
23*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(4,nnz_tot)    :: lhsP
24*59599516SKenneth E. Jansen
25*59599516SKenneth E. Jansen
26*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(amg_nshg(level)) :: v,r
27*59599516SKenneth E. Jansen      real(kind=8),intent(inout),dimension(amg_nshg(level)) :: u
28*59599516SKenneth E. Jansen      integer,intent(in) :: fwdbck
29*59599516SKenneth E. Jansen
30*59599516SKenneth E. Jansen      ! Local
31*59599516SKenneth E. Jansen      real(kind=8),dimension(amg_nshg(level)) :: pAux
32*59599516SKenneth E. Jansen      real(kind=8),dimension(amg_nshg(level)) :: epdk,epx
33*59599516SKenneth E. Jansen
34*59599516SKenneth E. Jansen      real(kind=8),dimension(amg_nshg(level)) :: cy1,cz1
35*59599516SKenneth E. Jansen
36*59599516SKenneth E. Jansen      real(kind=8) :: cf,beta,alpha,delta,theta,s1,rhok,invtheta
37*59599516SKenneth E. Jansen      real(kind=8) :: rhokp1,dtemp1,dtemp2,twoodelta
38*59599516SKenneth E. Jansen      real(kind=8) :: chebyniu,chebyomega,chebyck1,chebyck2,chebygamma
39*59599516SKenneth E. Jansen      integer :: i,j,k
40*59599516SKenneth E. Jansen
41*59599516SKenneth E. Jansen      beta = 1.1*mlsCf(level,1) ! 1.1 * Ev_max
42*59599516SKenneth E. Jansen      alpha = mlsCf(level,1)*ramg_chebyratio
43*59599516SKenneth E. Jansen
44*59599516SKenneth E. Jansen      IF (.FALSE.) THEN
45*59599516SKenneth E. Jansen      ! mine ( matrix computing )
46*59599516SKenneth E. Jansen      chebyniu = 1.0+2.0*(1.0-beta)/(beta-alpha)
47*59599516SKenneth E. Jansen      chebygamma = 2.0/(2.0-alpha-beta)
48*59599516SKenneth E. Jansen      chebyomega = chebyniu/(2.0*chebyniu*chebyniu-1.0)
49*59599516SKenneth E. Jansen      chebyomega = chebyomega*2.0*(2.0-beta-alpha)/(beta-alpha)
50*59599516SKenneth E. Jansen      if (fwdbck.eq.2) then ! initial guess x0 = 0
51*59599516SKenneth E. Jansen          cy1 = r
52*59599516SKenneth E. Jansen      else
53*59599516SKenneth E. Jansen          call ramg_calcAv_g(level,pAux,v,colm,rowp,lhsK,lhsP,
54*59599516SKenneth E. Jansen     &         ilwork,BC,iBC,iper,1)
55*59599516SKenneth E. Jansen          cy1 = v-1.1*3.0/beta*(pAux-r)
56*59599516SKenneth E. Jansen      endif
57*59599516SKenneth E. Jansen
58*59599516SKenneth E. Jansen      call ramg_calcAv_g(level,pAux,cy1,colm,rowp,lhsK,lhsP,
59*59599516SKenneth E. Jansen     &                   ilwork,BC,iBC,iper,1)
60*59599516SKenneth E. Jansen      cz1 = r-pAux
61*59599516SKenneth E. Jansen
62*59599516SKenneth E. Jansen      if (fwdbck.eq.2) then
63*59599516SKenneth E. Jansen          u = chebyomega*(cy1+chebygamma*cz1)
64*59599516SKenneth E. Jansen      else
65*59599516SKenneth E. Jansen          u = chebyomega*(cy1-v+chebygamma*cz1)+v
66*59599516SKenneth E. Jansen      endif
67*59599516SKenneth E. Jansen
68*59599516SKenneth E. Jansen      ELSE
69*59599516SKenneth E. Jansen      ! trilinos
70*59599516SKenneth E. Jansen      delta = (beta-alpha)/2.0
71*59599516SKenneth E. Jansen      twoodelta = 4.0/(beta-alpha)
72*59599516SKenneth E. Jansen      theta = (beta+alpha)/2.0
73*59599516SKenneth E. Jansen      invtheta = 2.0/(beta+alpha)
74*59599516SKenneth E. Jansen      s1 = 2.0*theta/delta
75*59599516SKenneth E. Jansen      rhok = delta/theta
76*59599516SKenneth E. Jansen
77*59599516SKenneth E. Jansen      !write(*,*)delta,twoodelta,2.0/delta
78*59599516SKenneth E. Jansen      if (fwdbck.eq.2) then ! pre-smoothing, initial guess v=0
79*59599516SKenneth E. Jansen          epdk = r*invtheta!/theta
80*59599516SKenneth E. Jansen          epx  = epdk
81*59599516SKenneth E. Jansen      else
82*59599516SKenneth E. Jansen          call ramg_calcAv_g(level,pAux,v,colm,rowp,lhsK,lhsP,
83*59599516SKenneth E. Jansen     &         ilwork,BC,iBC,iper,1)
84*59599516SKenneth E. Jansen          epdk = r-pAux
85*59599516SKenneth E. Jansen          epdk = epdk*invtheta!/theta
86*59599516SKenneth E. Jansen          epx = v+epdk
87*59599516SKenneth E. Jansen      endif
88*59599516SKenneth E. Jansen
89*59599516SKenneth E. Jansen      do k=1,mlsDeg-1
90*59599516SKenneth E. Jansen          call ramg_calcAv_g(level,pAux,epx,colm,rowp,lhsK,lhsP,
91*59599516SKenneth E. Jansen     &         ilwork,BC,iBC,iper,1)
92*59599516SKenneth E. Jansen          rhokp1 = 1.0/(s1-rhok)
93*59599516SKenneth E. Jansen          dtemp1 = rhokp1*rhok
94*59599516SKenneth E. Jansen          dtemp2 = rhokp1*twoodelta
95*59599516SKenneth E. Jansen          rhok = rhokp1
96*59599516SKenneth E. Jansen          epdk = dtemp1*epdk+dtemp2*(r-pAux)
97*59599516SKenneth E. Jansen          epx = epx + epdk
98*59599516SKenneth E. Jansen      enddo
99*59599516SKenneth E. Jansen      u = epx
100*59599516SKenneth E. Jansen
101*59599516SKenneth E. Jansen      ENDIF
102*59599516SKenneth E. Jansen
103*59599516SKenneth E. Jansen      end subroutine ! ramg_cheby_apply
104*59599516SKenneth E. Jansen
105*59599516SKenneth E. Jansen
106*59599516SKenneth E. Jansen      subroutine ramg_cheby_setup(colm,rowp,lhsK,lhsP,
107*59599516SKenneth E. Jansen     &                          ilwork,BC,iBC,iper)
108*59599516SKenneth E. Jansen      use ramg_data
109*59599516SKenneth E. Jansen      include "common.h"
110*59599516SKenneth E. Jansen      include "mpif.h"
111*59599516SKenneth E. Jansen      include "auxmpi.h"
112*59599516SKenneth E. Jansen
113*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg+1)            :: colm
114*59599516SKenneth E. Jansen      integer,intent(in),dimension(nnz_tot)           :: rowp
115*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(9,nnz_tot)    :: lhsK
116*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(4,nnz_tot)    :: lhsP
117*59599516SKenneth E. Jansen      integer,intent(in),dimension(nlwork)            :: ilwork
118*59599516SKenneth E. Jansen      integer,intent(in),dimension(nshg)              :: iBC,iper
119*59599516SKenneth E. Jansen      real(kind=8),intent(in),dimension(nshg,ndofBC)  :: BC
120*59599516SKenneth E. Jansen
121*59599516SKenneth E. Jansen      integer :: i,j,k,lvl,ideg
122*59599516SKenneth E. Jansen
123*59599516SKenneth E. Jansen      if (ramg_setup_flag .ne. 0 ) return
124*59599516SKenneth E. Jansen
125*59599516SKenneth E. Jansen      if (myrank.eq.master) then
126*59599516SKenneth E. Jansen         write(*,*)'Start Chebyshev setup'
127*59599516SKenneth E. Jansen      endif
128*59599516SKenneth E. Jansen      do lvl=1,ramg_levelx
129*59599516SKenneth E. Jansen      call ramg_mls_eigen(mlsCf(lvl,1),lvl,1,colm,rowp,lhsK,lhsP,
130*59599516SKenneth E. Jansen     &                    ilwork,BC,iBC,iper)
131*59599516SKenneth E. Jansen!      call ramg_mls_min_eigen(mlsCf(lvl,2),lvl,1,colm,rowp,lhsK,lhsP,
132*59599516SKenneth E. Jansen!     &                    ilwork,BC,iBC,iper)
133*59599516SKenneth E. Jansen!      mlsCf(lvl,2) = mlsCf(lvl,1)
134*59599516SKenneth E. Jansen      if (myrank.eq.master) then
135*59599516SKenneth E. Jansen          write(*,7000)lvl,mlsCf(lvl,1)
136*59599516SKenneth E. Jansen7000      format('Level: ',I2,T15,'Ev_max: ',E10.4)
137*59599516SKenneth E. Jansen      endif
138*59599516SKenneth E. Jansen      enddo
139*59599516SKenneth E. Jansen      if (myrank.eq.master) then
140*59599516SKenneth E. Jansen          write(*,*)'End Chebyshev setup'
141*59599516SKenneth E. Jansen      endif
142*59599516SKenneth E. Jansen
143*59599516SKenneth E. Jansen      end subroutine !ramg_cheby_setup
144