xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision 63b9fa5e258ef37d5bcff491ebacde33f30ed49b)
1 /*$Id: mmaij.c,v 1.59 2001/08/07 03:02:49 balay Exp $*/
2 
3 /*
4    Support for the parallel AIJ matrix vector multiply
5 */
6 #include "src/mat/impls/aij/mpi/mpiaij.h"
7 #include "src/vec/vecimpl.h"
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "MatSetUpMultiply_MPIAIJ"
11 int MatSetUpMultiply_MPIAIJ(Mat mat)
12 {
13   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
14   Mat_SeqAIJ         *B = (Mat_SeqAIJ*)(aij->B->data);
15   int                N = mat->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
16   IS                 from,to;
17   Vec                gvec;
18 #if defined (PETSC_USE_CTABLE)
19   PetscTable         gid1_lid1;
20   PetscTablePosition tpos;
21   int                gid,lid;
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 int 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   int          ierr,i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray;
143   int          *nz,ec;
144   PetscScalar  v;
145 
146   PetscFunctionBegin;
147   /* free stuff related to matrix-vec multiply */
148   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
149   ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0;
150   ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0;
151   if (aij->colmap) {
152 #if defined (PETSC_USE_CTABLE)
153     ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr);
154     aij->colmap = 0;
155 #else
156     ierr = PetscFree(aij->colmap);CHKERRQ(ierr);
157     aij->colmap = 0;
158     PetscLogObjectMemory(A,-aij->B->n*sizeof(int));
159 #endif
160   }
161 
162   /* make sure that B is assembled so we can access its values */
163   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165 
166   /* invent new B and copy stuff over */
167   ierr = PetscMalloc((m+1)*sizeof(int),&nz);CHKERRQ(ierr);
168   for (i=0; i<m; i++) {
169     nz[i] = Baij->i[i+1] - Baij->i[i];
170   }
171   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew);CHKERRQ(ierr);
172   ierr = PetscFree(nz);CHKERRQ(ierr);
173   for (i=0; i<m; i++) {
174     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
175       col  = garray[Baij->j[ct]];
176       v    = Baij->a[ct++];
177       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr);
178     }
179   }
180   ierr = PetscFree(aij->garray);CHKERRQ(ierr);
181   aij->garray = 0;
182   PetscLogObjectMemory(A,-ec*sizeof(int));
183   ierr = MatDestroy(B);CHKERRQ(ierr);
184   PetscLogObjectParent(A,Bnew);
185   aij->B = Bnew;
186   A->was_assembled = PETSC_FALSE;
187   PetscFunctionReturn(0);
188 }
189 
190 /*      ugly stuff added for Glenn someday we should fix this up */
191 
192 static int *auglyrmapd = 0,*auglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal"
193                                       parts of the local matrix */
194 static Vec auglydd = 0,auglyoo = 0;   /* work vectors used to scale the two parts of the local matrix */
195 
196 
197 #undef __FUNCT__
198 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp"
199 int MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
200 {
201   Mat_MPIAIJ  *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
202   int          ierr,i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
203   int          *r_rmapd,*r_rmapo;
204 
205   PetscFunctionBegin;
206   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
207   ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr);
208   ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);CHKERRQ(ierr);
209   ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));CHKERRQ(ierr);
210   nt   = 0;
211   for (i=0; i<inA->mapping->n; i++) {
212     if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) {
213       nt++;
214       r_rmapd[i] = inA->mapping->indices[i] + 1;
215     }
216   }
217   if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n);
218   ierr = PetscMalloc((n+1)*sizeof(int),&auglyrmapd);CHKERRQ(ierr);
219   for (i=0; i<inA->mapping->n; i++) {
220     if (r_rmapd[i]){
221       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
222     }
223   }
224   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
225   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr);
226 
227   ierr = PetscMalloc((inA->N+1)*sizeof(int),&lindices);CHKERRQ(ierr);
228   ierr = PetscMemzero(lindices,inA->N*sizeof(int));CHKERRQ(ierr);
229   for (i=0; i<ina->B->n; i++) {
230     lindices[garray[i]] = i+1;
231   }
232   no   = inA->mapping->n - nt;
233   ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);CHKERRQ(ierr);
234   ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));CHKERRQ(ierr);
235   nt   = 0;
236   for (i=0; i<inA->mapping->n; i++) {
237     if (lindices[inA->mapping->indices[i]]) {
238       nt++;
239       r_rmapo[i] = lindices[inA->mapping->indices[i]];
240     }
241   }
242   if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n);
243   ierr = PetscFree(lindices);CHKERRQ(ierr);
244   ierr = PetscMalloc((nt+1)*sizeof(int),&auglyrmapo);CHKERRQ(ierr);
245   for (i=0; i<inA->mapping->n; i++) {
246     if (r_rmapo[i]){
247       auglyrmapo[(r_rmapo[i]-1)] = i;
248     }
249   }
250   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
251   ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr);
252 
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal"
258 int MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
259 {
260   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
261   int ierr,(*f)(Mat,Vec);
262 
263   PetscFunctionBegin;
264   ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr);
265   if (f) {
266     ierr = (*f)(A,scale);CHKERRQ(ierr);
267   }
268   PetscFunctionReturn(0);
269 }
270 
271 EXTERN_C_BEGIN
272 #undef __FUNCT__
273 #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ"
274 int MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
275 {
276   Mat_MPIAIJ  *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
277   int         ierr,n,i;
278   PetscScalar *d,*o,*s;
279 
280   PetscFunctionBegin;
281   if (!auglyrmapd) {
282     ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
283   }
284 
285   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
286 
287   ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr);
288   ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr);
289   for (i=0; i<n; i++) {
290     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
291   }
292   ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr);
293   /* column scale "diagonal" portion of local matrix */
294   ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr);
295 
296   ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr);
297   ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr);
298   for (i=0; i<n; i++) {
299     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
300   }
301   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
302   ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr);
303   /* column scale "off-diagonal" portion of local matrix */
304   ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr);
305 
306   PetscFunctionReturn(0);
307 }
308 EXTERN_C_END
309 
310 
311