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