xref: /phasta/phSolver/AMG/ramg_ITAI.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen
2*59599516SKenneth E. Jansen!***********************************************************
3*59599516SKenneth E. Jansen!      mtx_colm,mtx_rowp,mtx_mtx
4*59599516SKenneth E. Jansen!**********************************************************
5*59599516SKenneth E. Jansen      subroutine mtx_colm(a_colm,a_rowp,
6*59599516SKenneth E. Jansen     &                    b_colm,b_rowp,
7*59599516SKenneth E. Jansen     &                    c_colm,m,l,n)
8*59599516SKenneth E. Jansen      integer,intent(in)                :: m,l,n
9*59599516SKenneth E. Jansen      integer,intent(in) :: a_colm(m+1),b_colm(l+1)
10*59599516SKenneth E. Jansen      integer,intent(inout) :: c_colm(m+1)
11*59599516SKenneth E. Jansen      integer,intent(in) :: a_rowp(a_colm(m+1)-1)
12*59599516SKenneth E. Jansen      integer,intent(in) :: b_rowp(b_colm(l+1)-1)
13*59599516SKenneth E. Jansen
14*59599516SKenneth E. Jansen      integer      :: i,j,k,p,q,mnnz
15*59599516SKenneth E. Jansen      integer      :: itemp(max(m,l))
16*59599516SKenneth E. Jansen
17*59599516SKenneth E. Jansen      mnnz = 1
18*59599516SKenneth E. Jansen      itemp = 0
19*59599516SKenneth E. Jansen      do i=1,m
20*59599516SKenneth E. Jansen         c_colm(i) = mnnz
21*59599516SKenneth E. Jansen         do k=a_colm(i),a_colm(i+1)-1
22*59599516SKenneth E. Jansen            j=a_rowp(k)
23*59599516SKenneth E. Jansen            do p=b_colm(j),b_colm(j+1)-1
24*59599516SKenneth E. Jansen               q=b_rowp(p)
25*59599516SKenneth E. Jansen               if (itemp(q).ne.i) then
26*59599516SKenneth E. Jansen                   itemp(q)=i
27*59599516SKenneth E. Jansen                   mnnz = mnnz + 1
28*59599516SKenneth E. Jansen               end if
29*59599516SKenneth E. Jansen            enddo
30*59599516SKenneth E. Jansen         enddo
31*59599516SKenneth E. Jansen      enddo
32*59599516SKenneth E. Jansen      c_colm(m+1) = mnnz
33*59599516SKenneth E. Jansen
34*59599516SKenneth E. Jansen      end subroutine!mtx_colm
35*59599516SKenneth E. Jansen
36*59599516SKenneth E. Jansen      subroutine mtx_rowp(a_colm,a_rowp,b_colm,b_rowp,
37*59599516SKenneth E. Jansen     &                    c_colm,c_rowp,m,l,n,diagflag)
38*59599516SKenneth E. Jansen      integer,intent(in)                :: m,l,n
39*59599516SKenneth E. Jansen      integer,intent(in) :: a_colm(m+1),b_colm(l+1)
40*59599516SKenneth E. Jansen      integer,intent(in) :: c_colm(m+1)
41*59599516SKenneth E. Jansen      integer,intent(inout) :: c_rowp(c_colm(m+1)-1)
42*59599516SKenneth E. Jansen      integer,intent(in) :: a_rowp(a_colm(m+1)-1)
43*59599516SKenneth E. Jansen      integer,intent(in) :: b_rowp(b_colm(l+1)-1)
44*59599516SKenneth E. Jansen      logical,intent(in) :: diagflag
45*59599516SKenneth E. Jansen
46*59599516SKenneth E. Jansen      integer      :: i,j,k,p,q,mnnz
47*59599516SKenneth E. Jansen      integer      :: itemp(max(m,l))
48*59599516SKenneth E. Jansen
49*59599516SKenneth E. Jansen      mnnz = 0
50*59599516SKenneth E. Jansen      itemp = 0
51*59599516SKenneth E. Jansen      do i=1,m
52*59599516SKenneth E. Jansen         if (diagflag) then
53*59599516SKenneth E. Jansen         mnnz = mnnz + 1
54*59599516SKenneth E. Jansen         c_rowp(mnnz) = i
55*59599516SKenneth E. Jansen         itemp(i)=i
56*59599516SKenneth E. Jansen         endif
57*59599516SKenneth E. Jansen         do k=a_colm(i),a_colm(i+1)-1
58*59599516SKenneth E. Jansen            j=a_rowp(k)
59*59599516SKenneth E. Jansen            do p=b_colm(j),b_colm(j+1)-1
60*59599516SKenneth E. Jansen               q=b_rowp(p)
61*59599516SKenneth E. Jansen               if (itemp(q).ne.i) then
62*59599516SKenneth E. Jansen                   itemp(q)=i
63*59599516SKenneth E. Jansen                   mnnz = mnnz + 1
64*59599516SKenneth E. Jansen                   c_rowp(mnnz) = q
65*59599516SKenneth E. Jansen               end if
66*59599516SKenneth E. Jansen            enddo
67*59599516SKenneth E. Jansen         enddo
68*59599516SKenneth E. Jansen      enddo
69*59599516SKenneth E. Jansen
70*59599516SKenneth E. Jansen      end subroutine !mtx_rowp
71*59599516SKenneth E. Jansen
72*59599516SKenneth E. Jansen      subroutine mtx_mtx(a_colm,a_rowp,a_mtx,
73*59599516SKenneth E. Jansen     &                   b_colm,b_rowp,b_mtx,
74*59599516SKenneth E. Jansen     &                   c_colm,c_rowp,c_mtx,
75*59599516SKenneth E. Jansen     &                   m,l,n)
76*59599516SKenneth E. Jansen      integer,intent(in)                :: m,l,n
77*59599516SKenneth E. Jansen      integer,intent(in) :: a_colm(m+1),b_colm(l+1)
78*59599516SKenneth E. Jansen      integer,intent(in) :: c_colm(m+1)
79*59599516SKenneth E. Jansen      integer,intent(in) :: c_rowp(c_colm(m+1)-1)
80*59599516SKenneth E. Jansen      integer,intent(in) :: a_rowp(a_colm(m+1)-1)
81*59599516SKenneth E. Jansen      integer,intent(in) :: b_rowp(b_colm(l+1)-1)
82*59599516SKenneth E. Jansen
83*59599516SKenneth E. Jansen      real(kind=8),intent(in) :: a_mtx(a_colm(m+1)-1)
84*59599516SKenneth E. Jansen      real(kind=8),intent(in) :: b_mtx(b_colm(l+1)-1)
85*59599516SKenneth E. Jansen      real(kind=8),intent(inout) :: c_mtx(c_colm(m+1)-1)
86*59599516SKenneth E. Jansen
87*59599516SKenneth E. Jansen      integer      :: i,j,k,p,q
88*59599516SKenneth E. Jansen      integer      :: itemp(max(m,l))
89*59599516SKenneth E. Jansen      real(kind=8) :: rtemp
90*59599516SKenneth E. Jansen
91*59599516SKenneth E. Jansen      c_mtx = 0
92*59599516SKenneth E. Jansen      do i=1,m
93*59599516SKenneth E. Jansen         do k=c_colm(i),c_colm(i+1)-1
94*59599516SKenneth E. Jansen            itemp(c_rowp(k)) = k
95*59599516SKenneth E. Jansen         enddo
96*59599516SKenneth E. Jansen         do k=a_colm(i),a_colm(i+1)-1
97*59599516SKenneth E. Jansen            j = a_rowp(k)
98*59599516SKenneth E. Jansen            rtemp = a_mtx(k)
99*59599516SKenneth E. Jansen            do p = b_colm(j),b_colm(j+1)-1
100*59599516SKenneth E. Jansen               q = b_rowp(p)
101*59599516SKenneth E. Jansen               c_mtx(itemp(q)) = c_mtx(itemp(q)) + rtemp*b_mtx(p)
102*59599516SKenneth E. Jansen            enddo
103*59599516SKenneth E. Jansen         enddo
104*59599516SKenneth E. Jansen      enddo
105*59599516SKenneth E. Jansen
106*59599516SKenneth E. Jansen      end subroutine !mtx_mtx
107*59599516SKenneth E. Jansen
108*59599516SKenneth E. Jansen!************************************************************
109*59599516SKenneth E. Jansen!     ramg_calcIAI
110*59599516SKenneth E. Jansen!        do matrix multiplication for I^T A I
111*59599516SKenneth E. Jansen!        actual operation is: A2_ij = I^T_ip A_pq I_qj
112*59599516SKenneth E. Jansen!        where p,q are repeated indicies
113*59599516SKenneth E. Jansen!
114*59599516SKenneth E. Jansen!************************************************************
115*59599516SKenneth E. Jansen      subroutine ramg_calcITAI(level1,level2,maxstopsign)
116*59599516SKenneth E. Jansen      use ramg_data
117*59599516SKenneth E. Jansen      implicit none
118*59599516SKenneth E. Jansen
119*59599516SKenneth E. Jansen      integer,intent(in)             :: level1,level2
120*59599516SKenneth E. Jansen      logical,intent(inout)          :: maxstopsign
121*59599516SKenneth E. Jansen
122*59599516SKenneth E. Jansen      integer                        :: i,j,k,m,p,q,r,n,s,p0
123*59599516SKenneth E. Jansen      integer,dimension(amg_nshg(level1))  :: itemp
124*59599516SKenneth E. Jansen      integer,dimension(amg_nshg(level2))  :: itemp2
125*59599516SKenneth E. Jansen      integer                              :: mnnz
126*59599516SKenneth E. Jansen
127*59599516SKenneth E. Jansen      integer,dimension(amg_nshg(level2)+1)   :: ITA_colm
128*59599516SKenneth E. Jansen      integer,dimension(:),allocatable        :: ITA_rowp
129*59599516SKenneth E. Jansen      real(kind=8),dimension(:),allocatable   :: ITA
130*59599516SKenneth E. Jansen
131*59599516SKenneth E. Jansen      integer                              :: mem_err, mem_err_s
132*59599516SKenneth E. Jansen
133*59599516SKenneth E. Jansen      real(kind=8)                             :: rtp
134*59599516SKenneth E. Jansen      character                            :: fname*80
135*59599516SKenneth E. Jansen      real,dimension(10)                   :: cpusec
136*59599516SKenneth E. Jansen
137*59599516SKenneth E. Jansen      call cpu_time(cpusec(1))
138*59599516SKenneth E. Jansen
139*59599516SKenneth E. Jansen      IF (level1.eq.0) THEN
140*59599516SKenneth E. Jansen          ! reduced serial test here, DO NOT USE IT unless RECOVER
141*59599516SKenneth E. Jansen          ! from extract.
142*59599516SKenneth E. Jansen      call mtx_colm(I_fc_colm(level1)%p,I_fc_rowp(level1)%p,
143*59599516SKenneth E. Jansen!     &              amg_A_colm(level1)%p,amg_A_rowp(level1)%p,
144*59599516SKenneth E. Jansen     &              reducecolm,reducerowp,
145*59599516SKenneth E. Jansen     &              ITA_colm,amg_nshg(level2),amg_nshg(level1),
146*59599516SKenneth E. Jansen     &              amg_nshg(level1))
147*59599516SKenneth E. Jansen      mnnz = ITA_colm(amg_nshg(level2)+1)-1
148*59599516SKenneth E. Jansen      allocate(ITA_rowp(mnnz),stat=mem_err)
149*59599516SKenneth E. Jansen      allocate(ITA(mnnz),stat=mem_err)
150*59599516SKenneth E. Jansen      ! IT * A
151*59599516SKenneth E. Jansen      ! rowp
152*59599516SKenneth E. Jansen      call mtx_rowp(I_fc_colm(level1)%p,I_fc_rowp(level1)%p,
153*59599516SKenneth E. Jansen!     &              amg_A_colm(level1)%p,amg_A_rowp(level1)%p,
154*59599516SKenneth E. Jansen     &              reducecolm,reducerowp,
155*59599516SKenneth E. Jansen     &              ITA_colm,ITA_rowp,
156*59599516SKenneth E. Jansen     &              amg_nshg(level2),amg_nshg(level1),
157*59599516SKenneth E. Jansen     &              amg_nshg(level1),.false.)
158*59599516SKenneth E. Jansen      ! IT * A
159*59599516SKenneth E. Jansen      ! value
160*59599516SKenneth E. Jansen      call mtx_mtx(I_fc_colm(level1)%p,I_fc_rowp(level1)%p,
161*59599516SKenneth E. Jansen     &              I_fc(level1)%p,
162*59599516SKenneth E. Jansen!     &              amg_A_colm(level1)%p,amg_A_rowp(level1)%p,
163*59599516SKenneth E. Jansen!     &              amg_A_lhs(level1)%p,
164*59599516SKenneth E. Jansen     &              reducecolm,reducerowp,reducelhs,
165*59599516SKenneth E. Jansen     &              ITA_colm,ITA_rowp,ITA,
166*59599516SKenneth E. Jansen     &              amg_nshg(level2),amg_nshg(level1),
167*59599516SKenneth E. Jansen     &              amg_nshg(level1))
168*59599516SKenneth E. Jansen      call cpu_time(cpusec(2))
169*59599516SKenneth E. Jansen
170*59599516SKenneth E. Jansen      ELSE
171*59599516SKenneth E. Jansen
172*59599516SKenneth E. Jansen      ! IT * A
173*59599516SKenneth E. Jansen      ! colm and nnz
174*59599516SKenneth E. Jansen      call mtx_colm(I_fc_colm(level1)%p,I_fc_rowp(level1)%p,
175*59599516SKenneth E. Jansen     &              amg_A_colm(level1)%p,amg_A_rowp(level1)%p,
176*59599516SKenneth E. Jansen     &              ITA_colm,amg_nshg(level2),amg_nshg(level1),
177*59599516SKenneth E. Jansen     &              amg_nshg(level1))
178*59599516SKenneth E. Jansen      mnnz = ITA_colm(amg_nshg(level2)+1)-1
179*59599516SKenneth E. Jansen      allocate(ITA_rowp(mnnz),stat=mem_err)
180*59599516SKenneth E. Jansen      allocate(ITA(mnnz),stat=mem_err)
181*59599516SKenneth E. Jansen      ! IT * A
182*59599516SKenneth E. Jansen      ! rowp
183*59599516SKenneth E. Jansen      call mtx_rowp(I_fc_colm(level1)%p,I_fc_rowp(level1)%p,
184*59599516SKenneth E. Jansen     &              amg_A_colm(level1)%p,amg_A_rowp(level1)%p,
185*59599516SKenneth E. Jansen     &              ITA_colm,ITA_rowp,
186*59599516SKenneth E. Jansen     &              amg_nshg(level2),amg_nshg(level1),
187*59599516SKenneth E. Jansen     &              amg_nshg(level1),.false.)
188*59599516SKenneth E. Jansen      ! IT * A
189*59599516SKenneth E. Jansen      ! value
190*59599516SKenneth E. Jansen      call mtx_mtx(I_fc_colm(level1)%p,I_fc_rowp(level1)%p,
191*59599516SKenneth E. Jansen     &              I_fc(level1)%p,
192*59599516SKenneth E. Jansen     &              amg_A_colm(level1)%p,amg_A_rowp(level1)%p,
193*59599516SKenneth E. Jansen     &              amg_A_lhs(level1)%p,
194*59599516SKenneth E. Jansen     &              ITA_colm,ITA_rowp,ITA,
195*59599516SKenneth E. Jansen     &              amg_nshg(level2),amg_nshg(level1),
196*59599516SKenneth E. Jansen     &              amg_nshg(level1))
197*59599516SKenneth E. Jansen      call cpu_time(cpusec(2))
198*59599516SKenneth E. Jansen
199*59599516SKenneth E. Jansen      ENDIF
200*59599516SKenneth E. Jansen
201*59599516SKenneth E. Jansen      ! ( IT*A ) * I
202*59599516SKenneth E. Jansen      ! colm and nnz
203*59599516SKenneth E. Jansen      call ramg_allocate(level2,amg_nshg(level2),0,1)
204*59599516SKenneth E. Jansen      call mtx_colm(ITA_colm,ITA_rowp,I_cf_colm(level1)%p,
205*59599516SKenneth E. Jansen     &              I_cf_rowp(level1)%p,amg_A_colm(level2)%p,
206*59599516SKenneth E. Jansen     &              amg_nshg(level2),amg_nshg(level1),
207*59599516SKenneth E. Jansen     &              amg_nshg(level2))
208*59599516SKenneth E. Jansen
209*59599516SKenneth E. Jansen      call cpu_time(cpusec(3))
210*59599516SKenneth E. Jansen      mnnz = amg_A_colm(level2)%p(amg_nshg(level2)+1)-1
211*59599516SKenneth E. Jansen      call ramg_allocate(level2,0,mnnz,1)
212*59599516SKenneth E. Jansen      ! ( IT*A ) * I
213*59599516SKenneth E. Jansen      ! rowp
214*59599516SKenneth E. Jansen      call mtx_rowp(ITA_colm,ITA_rowp,I_cf_colm(level1)%p,
215*59599516SKenneth E. Jansen     &              I_cf_rowp(level1)%p,amg_A_colm(level2)%p,
216*59599516SKenneth E. Jansen     &              amg_A_rowp(level2)%p,
217*59599516SKenneth E. Jansen     &              amg_nshg(level2),amg_nshg(level1),
218*59599516SKenneth E. Jansen     &              amg_nshg(level2),.true.)
219*59599516SKenneth E. Jansen      call cpu_time(cpusec(4))
220*59599516SKenneth E. Jansen      ! ( IT*A ) * I
221*59599516SKenneth E. Jansen      ! value
222*59599516SKenneth E. Jansen      call mtx_mtx(ITA_colm,ITA_rowp,ITA,I_cf_colm(level1)%p,
223*59599516SKenneth E. Jansen     &              I_cf_rowp(level1)%p,I_cf(level1)%p,
224*59599516SKenneth E. Jansen     &              amg_A_colm(level2)%p,
225*59599516SKenneth E. Jansen     &              amg_A_rowp(level2)%p,amg_A_lhs(level2)%p,
226*59599516SKenneth E. Jansen     &              amg_nshg(level2),amg_nshg(level1),
227*59599516SKenneth E. Jansen     &              amg_nshg(level2))
228*59599516SKenneth E. Jansen
229*59599516SKenneth E. Jansen      deallocate(ITA_rowp,stat=mem_err)
230*59599516SKenneth E. Jansen      deallocate(ITA,stat=mem_err)
231*59599516SKenneth E. Jansen
232*59599516SKenneth E. Jansen      call cpu_time(cpusec(5))
233*59599516SKenneth E. Jansen
234*59599516SKenneth E. Jansen      amg_nshg(level2+1) = 0
235*59599516SKenneth E. Jansen
236*59599516SKenneth E. Jansen      ! put on right hand side
237*59599516SKenneth E. Jansen      call ramg_calcIvFC(level1,level2,amg_A_rhs(level1)%p,
238*59599516SKenneth E. Jansen     &               amg_A_rhs(level2)%p)
239*59599516SKenneth E. Jansen
240*59599516SKenneth E. Jansen      call cpu_time(cpusec(10))
241*59599516SKenneth E. Jansen
242*59599516SKenneth E. Jansen      do i=1,amg_nshg(level2)
243*59599516SKenneth E. Jansen         amg_ppeDiag(level2)%p(i) =
244*59599516SKenneth E. Jansen     &   1.0/sqrt(amg_A_lhs(level2)%p(amg_A_colm(level2)%p(i),1))
245*59599516SKenneth E. Jansen      enddo
246*59599516SKenneth E. Jansen      do i=1,amg_nshg(level2)
247*59599516SKenneth E. Jansen         do m=amg_A_colm(level2)%p(i),amg_A_colm(level2)%p(i+1)-1
248*59599516SKenneth E. Jansen            j = amg_A_rowp(level2)%p(m)
249*59599516SKenneth E. Jansen            amg_A_lhs(level2)%p(m,1)=amg_A_lhs(level2)%p(m,1)*
250*59599516SKenneth E. Jansen     &        amg_ppeDiag(level2)%p(i)*amg_ppeDiag(level2)%p(j)
251*59599516SKenneth E. Jansen         enddo
252*59599516SKenneth E. Jansen      enddo
253*59599516SKenneth E. Jansen
254*59599516SKenneth E. Jansen!      do i=1,amg_nshg(level1)
255*59599516SKenneth E. Jansen!         do j=I_cf_colm(level1)%p(i),I_cf_colm(level1)%p(i+1)-1
256*59599516SKenneth E. Jansen!            k = I_cf_rowp(level1)%p(j)
257*59599516SKenneth E. Jansen!            I_cf(level1)%p(j) =
258*59599516SKenneth E. Jansen!     &      I_cf(level1)%p(j)*amg_ppeDiag(level2)%p(k)
259*59599516SKenneth E. Jansen!         enddo
260*59599516SKenneth E. Jansen!      enddo
261*59599516SKenneth E. Jansen
262*59599516SKenneth E. Jansen!      do i=1,amg_nshg(level2)
263*59599516SKenneth E. Jansen!         do j=I_fc_colm(level1)%p(i),I_fc_colm(level1)%p(i+1)-1
264*59599516SKenneth E. Jansen!         I_cf(level1)%p(j)=
265*59599516SKenneth E. Jansen!     &   I_cf(level1)%p(j)*amg_ppeDiag(level2)%p(i)
266*59599516SKenneth E. Jansen!         enddo
267*59599516SKenneth E. Jansen!      enddo
268*59599516SKenneth E. Jansen
269*59599516SKenneth E. Jansen      end subroutine ! ramg_calcIAI
270*59599516SKenneth E. Jansen
271*59599516SKenneth E. Jansen
272