Actual source code: mmbaij.c
2: /*
3: Support for the parallel BAIJ matrix vector multiply
4: */
5: #include <../src/mat/impls/baij/mpi/mpibaij.h>
6: #include <petsc/private/isimpl.h>
8: PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat)
9: {
10: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
11: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
12: PetscInt i,j,*aj = B->j,ec = 0,*garray;
13: PetscInt bs = mat->rmap->bs,*stmp;
14: IS from,to;
15: Vec gvec;
16: #if defined(PETSC_USE_CTABLE)
17: PetscTable gid1_lid1;
18: PetscTablePosition tpos;
19: PetscInt gid,lid;
20: #else
21: PetscInt Nbs = baij->Nbs,*indices;
22: #endif
24: #if defined(PETSC_USE_CTABLE)
25: /* use a table - Mark Adams */
26: PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
27: for (i=0; i<B->mbs; i++) {
28: for (j=0; j<B->ilen[i]; j++) {
29: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
30: PetscTableFind(gid1_lid1,gid1,&data);
31: if (!data) {
32: /* one based table */
33: PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
34: }
35: }
36: }
37: /* form array of columns we need */
38: PetscMalloc1(ec,&garray);
39: PetscTableGetHeadPosition(gid1_lid1,&tpos);
40: while (tpos) {
41: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
42: gid--; lid--;
43: garray[lid] = gid;
44: }
45: PetscSortInt(ec,garray);
46: PetscTableRemoveAll(gid1_lid1);
47: for (i=0; i<ec; i++) {
48: PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
49: }
50: /* compact out the extra columns in B */
51: for (i=0; i<B->mbs; i++) {
52: for (j=0; j<B->ilen[i]; j++) {
53: PetscInt gid1 = aj[B->i[i] + j] + 1;
54: PetscTableFind(gid1_lid1,gid1,&lid);
55: lid--;
56: aj[B->i[i]+j] = lid;
57: }
58: }
59: B->nbs = ec;
60: PetscLayoutDestroy(&baij->B->cmap);
61: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap);
62: PetscTableDestroy(&gid1_lid1);
63: #else
64: /* Make an array as long as the number of columns */
65: /* mark those columns that are in baij->B */
66: PetscCalloc1(Nbs,&indices);
67: for (i=0; i<B->mbs; i++) {
68: for (j=0; j<B->ilen[i]; j++) {
69: if (!indices[aj[B->i[i] + j]]) ec++;
70: indices[aj[B->i[i] + j]] = 1;
71: }
72: }
74: /* form array of columns we need */
75: PetscMalloc1(ec,&garray);
76: ec = 0;
77: for (i=0; i<Nbs; i++) {
78: if (indices[i]) {
79: garray[ec++] = i;
80: }
81: }
83: /* make indices now point into garray */
84: for (i=0; i<ec; i++) {
85: indices[garray[i]] = i;
86: }
88: /* compact out the extra columns in B */
89: for (i=0; i<B->mbs; i++) {
90: for (j=0; j<B->ilen[i]; j++) {
91: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
92: }
93: }
94: B->nbs = ec;
95: PetscLayoutDestroy(&baij->B->cmap);
96: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap);
97: PetscFree(indices);
98: #endif
100: /* create local vector that is used to scatter into */
101: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
103: /* create two temporary index sets for building scatter-gather */
104: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
106: PetscMalloc1(ec,&stmp);
107: for (i=0; i<ec; i++) stmp[i] = i;
108: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);
110: /* create temporary global vector to generate scatter context */
111: VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
113: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
114: VecScatterViewFromOptions(baij->Mvctx,(PetscObject)mat,"-matmult_vecscatter_view");
116: PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);
117: PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);
118: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
119: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
121: baij->garray = garray;
123: PetscLogObjectMemory((PetscObject)mat,ec*sizeof(PetscInt));
124: ISDestroy(&from);
125: ISDestroy(&to);
126: VecDestroy(&gvec);
127: return 0;
128: }
130: /*
131: Takes the local part of an already assembled MPIBAIJ matrix
132: and disassembles it. This is to allow new nonzeros into the matrix
133: that require more communication in the matrix vector multiply.
134: Thus certain data-structures must be rebuilt.
136: Kind of slow! But that's what application programmers get when
137: they are sloppy.
138: */
139: PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A)
140: {
141: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data;
142: Mat B = baij->B,Bnew;
143: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
144: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
145: PetscInt bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
146: MatScalar *a = Bbaij->a;
147: MatScalar *atmp;
149: /* free stuff related to matrix-vec multiply */
150: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
151: VecDestroy(&baij->lvec); baij->lvec = NULL;
152: VecScatterDestroy(&baij->Mvctx); baij->Mvctx = NULL;
153: if (baij->colmap) {
154: #if defined(PETSC_USE_CTABLE)
155: PetscTableDestroy(&baij->colmap);
156: #else
157: PetscFree(baij->colmap);
158: PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));
159: #endif
160: }
162: /* make sure that B is assembled so we can access its values */
163: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
164: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
166: /* invent new B and copy stuff over */
167: PetscMalloc1(mbs,&nz);
168: for (i=0; i<mbs; i++) {
169: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
170: }
171: MatCreate(PetscObjectComm((PetscObject)B),&Bnew);
172: MatSetSizes(Bnew,m,n,m,n);
173: MatSetType(Bnew,((PetscObject)B)->type_name);
174: MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
175: if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */
176: ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew;
177: }
179: MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);
180: /*
181: Ensure that B's nonzerostate is monotonically increasing.
182: Or should this follow the MatSetValuesBlocked() loop to preserve B's nonzerstate across a MatDisAssemble() call?
183: */
184: Bnew->nonzerostate = B->nonzerostate;
186: for (i=0; i<mbs; i++) {
187: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
188: col = garray[Bbaij->j[j]];
189: atmp = a + j*bs2;
190: MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);
191: }
192: }
193: MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);
195: PetscFree(nz);
196: PetscFree(baij->garray);
197: PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
198: MatDestroy(&B);
199: PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);
201: baij->B = Bnew;
202: A->was_assembled = PETSC_FALSE;
203: A->assembled = PETSC_FALSE;
204: return 0;
205: }
207: /* ugly stuff added for Glenn someday we should fix this up */
209: static PetscInt *uglyrmapd = NULL,*uglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
210: static Vec uglydd = NULL,uglyoo = NULL; /* work vectors used to scale the two parts of the local matrix */
212: PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
213: {
214: Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
215: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data;
216: PetscInt bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
217: PetscInt *r_rmapd,*r_rmapo;
219: MatGetOwnershipRange(inA,&cstart,&cend);
220: MatGetSize(ina->A,NULL,&n);
221: PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);
222: nt = 0;
223: for (i=0; i<inA->rmap->mapping->n; i++) {
224: if (inA->rmap->mapping->indices[i]*bs >= cstart && inA->rmap->mapping->indices[i]*bs < cend) {
225: nt++;
226: r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
227: }
228: }
230: PetscMalloc1(n+1,&uglyrmapd);
231: for (i=0; i<inA->rmap->mapping->n; i++) {
232: if (r_rmapd[i]) {
233: for (j=0; j<bs; j++) {
234: uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
235: }
236: }
237: }
238: PetscFree(r_rmapd);
239: VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);
241: PetscCalloc1(ina->Nbs+1,&lindices);
242: for (i=0; i<B->nbs; i++) {
243: lindices[garray[i]] = i+1;
244: }
245: no = inA->rmap->mapping->n - nt;
246: PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);
247: nt = 0;
248: for (i=0; i<inA->rmap->mapping->n; i++) {
249: if (lindices[inA->rmap->mapping->indices[i]]) {
250: nt++;
251: r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
252: }
253: }
255: PetscFree(lindices);
256: PetscMalloc1(nt*bs+1,&uglyrmapo);
257: for (i=0; i<inA->rmap->mapping->n; i++) {
258: if (r_rmapo[i]) {
259: for (j=0; j<bs; j++) {
260: uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
261: }
262: }
263: }
264: PetscFree(r_rmapo);
265: VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);
266: return 0;
267: }
269: PetscErrorCode MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
270: {
271: /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
273: PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
274: return 0;
275: }
277: PetscErrorCode MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
278: {
279: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
280: PetscInt n,i;
281: PetscScalar *d,*o;
282: const PetscScalar *s;
284: if (!uglyrmapd) {
285: MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);
286: }
288: VecGetArrayRead(scale,&s);
290: VecGetLocalSize(uglydd,&n);
291: VecGetArray(uglydd,&d);
292: for (i=0; i<n; i++) {
293: d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
294: }
295: VecRestoreArray(uglydd,&d);
296: /* column scale "diagonal" portion of local matrix */
297: MatDiagonalScale(a->A,NULL,uglydd);
299: VecGetLocalSize(uglyoo,&n);
300: VecGetArray(uglyoo,&o);
301: for (i=0; i<n; i++) {
302: o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
303: }
304: VecRestoreArrayRead(scale,&s);
305: VecRestoreArray(uglyoo,&o);
306: /* column scale "off-diagonal" portion of local matrix */
307: MatDiagonalScale(a->B,NULL,uglyoo);
308: return 0;
309: }