xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision dfbe8321245cdf9338262824a41580dace4e4143)
1 /*
2    Support for the parallel AIJ matrix vector multiply
3 */
4 #include "src/mat/impls/aij/mpi/mpiaij.h"
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "MatSetUpMultiply_MPIAIJ"
8 PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
9 {
10   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
11   Mat_SeqAIJ         *B = (Mat_SeqAIJ*)(aij->B->data);
12   int                i,j,*aj = B->j,ierr,ec = 0,*garray;
13   IS                 from,to;
14   Vec                gvec;
15 #if defined (PETSC_USE_CTABLE)
16   PetscTable         gid1_lid1;
17   PetscTablePosition tpos;
18   int                gid,lid;
19 #else
20   int                N = mat->N,*indices;
21 
22 #endif
23 
24   PetscFunctionBegin;
25 
26 #if defined (PETSC_USE_CTABLE)
27   /* use a table - Mark Adams (this has not been tested with "shift") */
28   ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr);
29   for (i=0; i<aij->B->m; i++) {
30     for (j=0; j<B->ilen[i]; j++) {
31       int data,gid1 = aj[B->i[i] + j] + 1;
32       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
33       if (!data) {
34         /* one based table */
35         ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr);
36       }
37     }
38   }
39   /* form array of columns we need */
40   ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr);
41   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
42   while (tpos) {
43     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
44     gid--;
45     lid--;
46     garray[lid] = gid;
47   }
48   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */
49   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
50   for (i=0; i<ec; i++) {
51     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr);
52   }
53   /* compact out the extra columns in B */
54   for (i=0; i<aij->B->m; i++) {
55     for (j=0; j<B->ilen[i]; j++) {
56       int gid1 = aj[B->i[i] + j] + 1;
57       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
58       lid --;
59       aj[B->i[i] + j]  = lid;
60     }
61   }
62   aij->B->n = aij->B->N = ec;
63   ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr);
64   /* Mark Adams */
65 #else
66   /* For the first stab we make an array as long as the number of columns */
67   /* mark those columns that are in aij->B */
68   ierr = PetscMalloc((N+1)*sizeof(int),&indices);CHKERRQ(ierr);
69   ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr);
70   for (i=0; i<aij->B->m; i++) {
71     for (j=0; j<B->ilen[i]; j++) {
72       if (!indices[aj[B->i[i] + j] ]) ec++;
73       indices[aj[B->i[i] + j] ] = 1;
74     }
75   }
76 
77   /* form array of columns we need */
78   ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr);
79   ec = 0;
80   for (i=0; i<N; i++) {
81     if (indices[i]) garray[ec++] = i;
82   }
83 
84   /* make indices now point into garray */
85   for (i=0; i<ec; i++) {
86     indices[garray[i]] = i;
87   }
88 
89   /* compact out the extra columns in B */
90   for (i=0; i<aij->B->m; i++) {
91     for (j=0; j<B->ilen[i]; j++) {
92       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
93     }
94   }
95   aij->B->n = aij->B->N = ec;
96   ierr = PetscFree(indices);CHKERRQ(ierr);
97 #endif
98   /* create local vector that is used to scatter into */
99   ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr);
100 
101   /* create two temporary Index sets for build scatter gather */
102   ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr);
103   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
104 
105   /* create temporary global vector to generate scatter context */
106   /* this is inefficient, but otherwise we must do either
107      1) save garray until the first actual scatter when the vector is known or
108      2) have another way of generating a scatter context without a vector.*/
109   ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr);
110 
111   /* generate the scatter context */
112   ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr);
113   PetscLogObjectParent(mat,aij->Mvctx);
114   PetscLogObjectParent(mat,aij->lvec);
115   PetscLogObjectParent(mat,from);
116   PetscLogObjectParent(mat,to);
117   aij->garray = garray;
118   PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
119   ierr = ISDestroy(from);CHKERRQ(ierr);
120   ierr = ISDestroy(to);CHKERRQ(ierr);
121   ierr = VecDestroy(gvec);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 
126 #undef __FUNCT__
127 #define __FUNCT__ "DisAssemble_MPIAIJ"
128 /*
129      Takes the local part of an already assembled MPIAIJ matrix
130    and disassembles it. This is to allow new nonzeros into the matrix
131    that require more communication in the matrix vector multiply.
132    Thus certain data-structures must be rebuilt.
133 
134    Kind of slow! But that's what application programmers get when
135    they are sloppy.
136 */
137 PetscErrorCode DisAssemble_MPIAIJ(Mat A)
138 {
139   Mat_MPIAIJ   *aij = (Mat_MPIAIJ*)A->data;
140   Mat          B = aij->B,Bnew;
141   Mat_SeqAIJ   *Baij = (Mat_SeqAIJ*)B->data;
142   PetscErrorCode ierr;
143   int          i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray;
144   int          *nz,ec;
145   PetscScalar  v;
146 
147   PetscFunctionBegin;
148   /* free stuff related to matrix-vec multiply */
149   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
150   ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0;
151   ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0;
152   if (aij->colmap) {
153 #if defined (PETSC_USE_CTABLE)
154     ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr);
155     aij->colmap = 0;
156 #else
157     ierr = PetscFree(aij->colmap);CHKERRQ(ierr);
158     aij->colmap = 0;
159     PetscLogObjectMemory(A,-aij->B->n*sizeof(int));
160 #endif
161   }
162 
163   /* make sure that B is assembled so we can access its values */
164   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166 
167   /* invent new B and copy stuff over */
168   ierr = PetscMalloc((m+1)*sizeof(int),&nz);CHKERRQ(ierr);
169   for (i=0; i<m; i++) {
170     nz[i] = Baij->i[i+1] - Baij->i[i];
171   }
172   ierr = MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);CHKERRQ(ierr);
173   ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr);
174   ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr);
175   ierr = PetscFree(nz);CHKERRQ(ierr);
176   for (i=0; i<m; i++) {
177     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
178       col  = garray[Baij->j[ct]];
179       v    = Baij->a[ct++];
180       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr);
181     }
182   }
183   ierr = PetscFree(aij->garray);CHKERRQ(ierr);
184   aij->garray = 0;
185   PetscLogObjectMemory(A,-ec*sizeof(int));
186   ierr = MatDestroy(B);CHKERRQ(ierr);
187   PetscLogObjectParent(A,Bnew);
188   aij->B = Bnew;
189   A->was_assembled = PETSC_FALSE;
190   PetscFunctionReturn(0);
191 }
192 
193 /*      ugly stuff added for Glenn someday we should fix this up */
194 
195 static int *auglyrmapd = 0,*auglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal"
196                                       parts of the local matrix */
197 static Vec auglydd = 0,auglyoo = 0;   /* work vectors used to scale the two parts of the local matrix */
198 
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp"
202 PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
203 {
204   Mat_MPIAIJ  *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
205   PetscErrorCode ierr;
206   int          i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
207   int          *r_rmapd,*r_rmapo;
208 
209   PetscFunctionBegin;
210   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
211   ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr);
212   ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);CHKERRQ(ierr);
213   ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));CHKERRQ(ierr);
214   nt   = 0;
215   for (i=0; i<inA->mapping->n; i++) {
216     if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) {
217       nt++;
218       r_rmapd[i] = inA->mapping->indices[i] + 1;
219     }
220   }
221   if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n);
222   ierr = PetscMalloc((n+1)*sizeof(int),&auglyrmapd);CHKERRQ(ierr);
223   for (i=0; i<inA->mapping->n; i++) {
224     if (r_rmapd[i]){
225       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
226     }
227   }
228   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
229   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr);
230 
231   ierr = PetscMalloc((inA->N+1)*sizeof(int),&lindices);CHKERRQ(ierr);
232   ierr = PetscMemzero(lindices,inA->N*sizeof(int));CHKERRQ(ierr);
233   for (i=0; i<ina->B->n; i++) {
234     lindices[garray[i]] = i+1;
235   }
236   no   = inA->mapping->n - nt;
237   ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);CHKERRQ(ierr);
238   ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));CHKERRQ(ierr);
239   nt   = 0;
240   for (i=0; i<inA->mapping->n; i++) {
241     if (lindices[inA->mapping->indices[i]]) {
242       nt++;
243       r_rmapo[i] = lindices[inA->mapping->indices[i]];
244     }
245   }
246   if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n);
247   ierr = PetscFree(lindices);CHKERRQ(ierr);
248   ierr = PetscMalloc((nt+1)*sizeof(int),&auglyrmapo);CHKERRQ(ierr);
249   for (i=0; i<inA->mapping->n; i++) {
250     if (r_rmapo[i]){
251       auglyrmapo[(r_rmapo[i]-1)] = i;
252     }
253   }
254   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
255   ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr);
256 
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal"
262 PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
263 {
264   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
265   PetscErrorCode ierr,(*f)(Mat,Vec);
266 
267   PetscFunctionBegin;
268   ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr);
269   if (f) {
270     ierr = (*f)(A,scale);CHKERRQ(ierr);
271   }
272   PetscFunctionReturn(0);
273 }
274 
275 EXTERN_C_BEGIN
276 #undef __FUNCT__
277 #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ"
278 PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
279 {
280   Mat_MPIAIJ  *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
281   PetscErrorCode ierr;
282   int         n,i;
283   PetscScalar *d,*o,*s;
284 
285   PetscFunctionBegin;
286   if (!auglyrmapd) {
287     ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
288   }
289 
290   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
291 
292   ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr);
293   ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr);
294   for (i=0; i<n; i++) {
295     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
296   }
297   ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr);
298   /* column scale "diagonal" portion of local matrix */
299   ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr);
300 
301   ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr);
302   ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr);
303   for (i=0; i<n; i++) {
304     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
305   }
306   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
307   ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr);
308   /* column scale "off-diagonal" portion of local matrix */
309   ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr);
310 
311   PetscFunctionReturn(0);
312 }
313 EXTERN_C_END
314 
315 
316