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