xref: /phasta/phSolver/compressible/asaugmr.f (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08)
1        subroutine AsAuGMR (ien,  EGmass,  uBrg, uBtmp )
2c
3c----------------------------------------------------------------------
4c This routine computes and assembles the Au product for the
5c GMRES solver.
6c
7c input:
8c     ien    (npro,nshl)       : nodal connectivity
9c     EGmass (npro,nedof,nedof)  : element mass matrix
10c     uBrg   (nshg,nflow)         : u_i before product
11c
12c output:
13c     uBrg   (nshg,nflow)         : result of product ( u_{i+1} )
14c
15c----------------------------------------------------------------------
16c
17        include "common.h"
18c
19        dimension ien(npro,nshl),     EGmass(npro,nedof,nedof),
20     &            uBrg(nshg,nflow),    uBtmp(nshg,nflow)
21c
22        dimension uBrgl(npro,nedof),  ubBgl(npro,nedof)
23c
24c.... localize this K-vector for the EBE product
25c
26        call localt (uBrg,  uBrgl,  abs(ien),  nflow,  'gather  ')
27
28        ubBgl = zero
29c
30c.... ----------------------->  Au product  <---------------------------
31c
32        do i = 1, nflow*nshl, nflow
33           do j = 1, nflow*nshl, nflow
34              ubBgl(:,i  ) = ubBgl(:,i  )
35     &                     + EGmass(:,i  ,j  ) * uBrgl(:,j  )
36     &                     + EGmass(:,i  ,j+1) * uBrgl(:,j+1)
37     &                     + EGmass(:,i  ,j+2) * uBrgl(:,j+2)
38     &                     + EGmass(:,i  ,j+3) * uBrgl(:,j+3)
39     &                     + EGmass(:,i  ,j+4) * uBrgl(:,j+4)
40c
41              ubBgl(:,i+1) = ubBgl(:,i+1)
42     &                     + EGmass(:,i+1,j  ) * uBrgl(:,j  )
43     &                     + EGmass(:,i+1,j+1) * uBrgl(:,j+1)
44     &                     + EGmass(:,i+1,j+2) * uBrgl(:,j+2)
45     &                     + EGmass(:,i+1,j+3) * uBrgl(:,j+3)
46     &                     + EGmass(:,i+1,j+4) * uBrgl(:,j+4)
47c
48              ubBgl(:,i+2) = ubBgl(:,i+2)
49     &                     + EGmass(:,i+2,j  ) * uBrgl(:,j  )
50     &                     + EGmass(:,i+2,j+1) * uBrgl(:,j+1)
51     &                     + EGmass(:,i+2,j+2) * uBrgl(:,j+2)
52     &                     + EGmass(:,i+2,j+3) * uBrgl(:,j+3)
53     &                     + EGmass(:,i+2,j+4) * uBrgl(:,j+4)
54c
55              ubBgl(:,i+3) = ubBgl(:,i+3)
56     &                     + EGmass(:,i+3,j  ) * uBrgl(:,j  )
57     &                     + EGmass(:,i+3,j+1) * uBrgl(:,j+1)
58     &                     + EGmass(:,i+3,j+2) * uBrgl(:,j+2)
59     &                     + EGmass(:,i+3,j+3) * uBrgl(:,j+3)
60     &                     + EGmass(:,i+3,j+4) * uBrgl(:,j+4)
61c
62              ubBgl(:,i+4) = ubBgl(:,i+4)
63     &                     + EGmass(:,i+4,j  ) * uBrgl(:,j  )
64     &                     + EGmass(:,i+4,j+1) * uBrgl(:,j+1)
65     &                     + EGmass(:,i+4,j+2) * uBrgl(:,j+2)
66     &                     + EGmass(:,i+4,j+3) * uBrgl(:,j+3)
67     &                     + EGmass(:,i+4,j+4) * uBrgl(:,j+4)
68c
69           enddo
70        enddo
71c
72c.... assemble the result of the product
73c
74        call localt (uBtmp,  ubBgl,  abs(ien),  nflow,  'scatter ')
75c
76c.... end
77c
78        return
79        end
80c
81c
82c
83        subroutine AsAuGMRSclr (ien,  EGmass,  uBrg, uBtmp )
84c
85c----------------------------------------------------------------------
86c This routine computes and assembles the Au product for the
87c GMRES solver.
88c
89c input:
90c     ien    (npro,nshl)       : nodal connectivity
91c     EGmass (npro,nen,nen)    : element mass matrix
92c     uBrg   (nshg)           : u_i before product
93c
94c output:
95c     uBtmp   (nshg)           : result of product ( u_{i+1} )
96c
97c----------------------------------------------------------------------
98c
99        include "common.h"
100        include "mpif.h"
101        include "auxmpi.h"
102c
103       dimension ien(npro,nshl),   EGmass(npro,nshape,nshape),
104     &           uBrg(nshg),       uBtmp(nshg)
105c
106       dimension uBrgl(npro,nshl), ubBgl(npro,nshl)
107c
108c.... localize this K-vector for the EBE product
109c
110        uBrgl = zero
111        call localtSclr(uBrg,  uBrgl,  ien,  'gather  ')
112
113        ubBgl = zero
114c
115c.... ----------------------->  Au product  <---------------------------
116c
117        do i = 1, nshl
118           do j = 1, nshl
119              ubBgl(:,i  ) = ubBgl(:,i  )
120     &                     + EGmass(:,i  ,j  ) * uBrgl(:,j  )
121c
122           enddo
123         enddo
124c
125c.... assemble the result of the product
126c
127        call localtSclr(uBtmp,  ubBgl,  ien,  'scatter ')
128c
129c.... end
130c
131        return
132        end
133
134
135
136
137
138
139