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