Actual source code: mmsell.c
petsc-3.9.4 2018-09-11
1: /*
2: Support for the parallel SELL matrix vector multiply
3: */
4: #include <../src/mat/impls/sell/mpi/mpisell.h>
5: #include <petsc/private/isimpl.h>
7: /*
8: Takes the local part of an already assembled MPISELL matrix
9: and disassembles it. This is to allow new nonzeros into the matrix
10: that require more communication in the matrix vector multiply.
11: Thus certain data-structures must be rebuilt.
13: Kind of slow! But that's what application programmers get when
14: they are sloppy.
15: */
16: PetscErrorCode MatDisAssemble_MPISELL(Mat A)
17: {
18: Mat_MPISELL *sell=(Mat_MPISELL*)A->data;
19: Mat B=sell->B,Bnew;
20: Mat_SeqSELL *Bsell=(Mat_SeqSELL*)B->data;
21: PetscInt i,j,totalslices,N=A->cmap->N,ec,row;
22: PetscBool isnonzero;
26: /* free stuff related to matrix-vec multiply */
27: VecGetSize(sell->lvec,&ec); /* needed for PetscLogObjectMemory below */
28: VecDestroy(&sell->lvec);
29: VecScatterDestroy(&sell->Mvctx);
30: if (sell->colmap) {
31: #if defined(PETSC_USE_CTABLE)
32: PetscTableDestroy(&sell->colmap);
33: #else
34: PetscFree(sell->colmap);
35: PetscLogObjectMemory((PetscObject)A,-sell->B->cmap->n*sizeof(PetscInt));
36: #endif
37: }
39: /* make sure that B is assembled so we can access its values */
40: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
41: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
43: /* invent new B and copy stuff over */
44: MatCreate(PETSC_COMM_SELF,&Bnew);
45: MatSetSizes(Bnew,B->rmap->n,N,B->rmap->n,N);
46: MatSetBlockSizesFromMats(Bnew,A,A);
47: MatSetType(Bnew,((PetscObject)B)->type_name);
48: MatSeqSELLSetPreallocation(Bnew,0,Bsell->rlen);
49: if (Bsell->nonew >= 0) { /* Inherit insertion error options (if positive). */
50: ((Mat_SeqSELL*)Bnew->data)->nonew = Bsell->nonew;
51: }
53: /*
54: Ensure that B's nonzerostate is monotonically increasing.
55: Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
56: */
57: Bnew->nonzerostate = B->nonzerostate;
59: totalslices = B->rmap->n/8+((B->rmap->n & 0x07)?1:0); /* floor(n/8) */
60: for (i=0; i<totalslices; i++) { /* loop over slices */
61: for (j=Bsell->sliidx[i],row=0; j<Bsell->sliidx[i+1]; j++,row=((row+1)&0x07)) {
62: isnonzero = (PetscBool)((j-Bsell->sliidx[i])/8 < Bsell->rlen[8*i+row]);
63: if (isnonzero) {
64: MatSetValue(Bnew,8*i+row,sell->garray[Bsell->colidx[j]],Bsell->val[j],B->insertmode);
65: }
66: }
67: }
69: PetscFree(sell->garray);
70: PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
71: MatDestroy(&B);
72: PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);
74: sell->B = Bnew;
75: A->was_assembled = PETSC_FALSE;
76: A->assembled = PETSC_FALSE;
77: return(0);
78: }
80: PetscErrorCode MatSetUpMultiply_MPISELL(Mat mat)
81: {
82: Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
83: Mat_SeqSELL *B=(Mat_SeqSELL*)(sell->B->data);
85: PetscInt i,j,*bcolidx=B->colidx,ec=0,*garray,totalslices;
86: IS from,to;
87: Vec gvec;
88: PetscBool isnonzero;
89: #if defined(PETSC_USE_CTABLE)
90: PetscTable gid1_lid1;
91: PetscTablePosition tpos;
92: PetscInt gid,lid;
93: #else
94: PetscInt N = mat->cmap->N,*indices;
95: #endif
98: totalslices = sell->B->rmap->n/8+((sell->B->rmap->n & 0x07)?1:0); /* floor(n/8) */
100: /* ec counts the number of columns that contain nonzeros */
101: #if defined(PETSC_USE_CTABLE)
102: /* use a table */
103: PetscTableCreate(sell->B->rmap->n,mat->cmap->N+1,&gid1_lid1);
104: for (i=0; i<totalslices; i++) { /* loop over slices */
105: for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
106: isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
107: if (isnonzero) { /* check the mask bit */
108: PetscInt data,gid1 = bcolidx[j] + 1;
109: PetscTableFind(gid1_lid1,gid1,&data);
110: if (!data) {
111: /* one based table */
112: PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
113: }
114: }
115: }
116: }
118: /* form array of columns we need */
119: PetscMalloc1(ec+1,&garray);
120: PetscTableGetHeadPosition(gid1_lid1,&tpos);
121: while (tpos) {
122: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
123: gid--;
124: lid--;
125: garray[lid] = gid;
126: }
127: PetscSortInt(ec,garray); /* sort, and rebuild */
128: PetscTableRemoveAll(gid1_lid1);
129: for (i=0; i<ec; i++) {
130: PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
131: }
133: /* compact out the extra columns in B */
134: for (i=0; i<totalslices; i++) { /* loop over slices */
135: for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
136: isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
137: if (isnonzero) {
138: PetscInt gid1 = bcolidx[j] + 1;
139: PetscTableFind(gid1_lid1,gid1,&lid);
140: lid--;
141: bcolidx[j] = lid;
142: }
143: }
144: }
145: sell->B->cmap->n = sell->B->cmap->N = ec;
146: sell->B->cmap->bs = 1;
148: PetscLayoutSetUp((sell->B->cmap));
149: PetscTableDestroy(&gid1_lid1);
150: #else
151: /* Make an array as long as the number of columns */
152: PetscCalloc1(N+1,&indices);
153: /* mark those columns that are in sell->B */
154: for (i=0; i<totalslices; i++) { /* loop over slices */
155: for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
156: isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
157: if (isnonzero) {
158: if (!indices[bcolidx[j]]) ec++;
159: indices[bcolidx[j]] = 1;
160: }
161: }
162: }
164: /* form array of columns we need */
165: PetscMalloc1(ec+1,&garray);
166: ec = 0;
167: for (i=0; i<N; i++) {
168: if (indices[i]) garray[ec++] = i;
169: }
171: /* make indices now point into garray */
172: for (i=0; i<ec; i++) {
173: indices[garray[i]] = i;
174: }
176: /* compact out the extra columns in B */
177: for (i=0; i<totalslices; i++) { /* loop over slices */
178: for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
179: isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
180: if (isnonzero) bcolidx[j] = indices[bcolidx[j]];
181: }
182: }
183: sell->B->cmap->n = sell->B->cmap->N = ec; /* number of columns that are not all zeros */
184: sell->B->cmap->bs = 1;
186: PetscLayoutSetUp((sell->B->cmap));
187: PetscFree(indices);
188: #endif
189: /* create local vector that is used to scatter into */
190: VecCreateSeq(PETSC_COMM_SELF,ec,&sell->lvec);
191: /* create two temporary Index sets for build scatter gather */
192: ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);
193: ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);
195: /* create temporary global vector to generate scatter context */
196: /* This does not allocate the array's memory so is efficient */
197: VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
199: /* generate the scatter context */
200: VecScatterCreate(gvec,from,sell->lvec,to,&sell->Mvctx);
201: PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->Mvctx);
202: PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->lvec);
203: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
204: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
206: sell->garray = garray;
208: PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
209: ISDestroy(&from);
210: ISDestroy(&to);
211: VecDestroy(&gvec);
212: return(0);
213: }
215: /* ugly stuff added for Glenn someday we should fix this up */
216: static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
217: static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */
219: PetscErrorCode MatMPISELLDiagonalScaleLocalSetUp(Mat inA,Vec scale)
220: {
221: Mat_MPISELL *ina=(Mat_MPISELL*)inA->data; /*access private part of matrix */
223: PetscInt i,n,nt,cstart,cend,no,*garray=ina->garray,*lindices;
224: PetscInt *r_rmapd,*r_rmapo;
227: MatGetOwnershipRange(inA,&cstart,&cend);
228: MatGetSize(ina->A,NULL,&n);
229: PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);
230: nt = 0;
231: for (i=0; i<inA->rmap->mapping->n; i++) {
232: if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
233: nt++;
234: r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
235: }
236: }
237: if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
238: PetscMalloc1(n+1,&auglyrmapd);
239: for (i=0; i<inA->rmap->mapping->n; i++) {
240: if (r_rmapd[i]) {
241: auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
242: }
243: }
244: PetscFree(r_rmapd);
245: VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);
246: PetscCalloc1(inA->cmap->N+1,&lindices);
247: for (i=0; i<ina->B->cmap->n; 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+1,&auglyrmapo);
262: for (i=0; i<inA->rmap->mapping->n; i++) {
263: if (r_rmapo[i]) {
264: auglyrmapo[(r_rmapo[i]-1)] = i;
265: }
266: }
267: PetscFree(r_rmapo);
268: VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);
269: return(0);
270: }
272: PetscErrorCode MatDiagonalScaleLocal_MPISELL(Mat A,Vec scale)
273: {
274: Mat_MPISELL *a=(Mat_MPISELL*)A->data; /*access private part of matrix */
275: PetscErrorCode ierr;
276: PetscInt n,i;
277: PetscScalar *d,*o;
278: const PetscScalar *s;
281: if (!auglyrmapd) {
282: MatMPISELLDiagonalScaleLocalSetUp(A,scale);
283: }
284: VecGetArrayRead(scale,&s);
285: VecGetLocalSize(auglydd,&n);
286: VecGetArray(auglydd,&d);
287: for (i=0; i<n; i++) {
288: d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
289: }
290: VecRestoreArray(auglydd,&d);
291: /* column scale "diagonal" portion of local matrix */
292: MatDiagonalScale(a->A,NULL,auglydd);
293: VecGetLocalSize(auglyoo,&n);
294: VecGetArray(auglyoo,&o);
295: for (i=0; i<n; i++) {
296: o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
297: }
298: VecRestoreArrayRead(scale,&s);
299: VecRestoreArray(auglyoo,&o);
300: /* column scale "off-diagonal" portion of local matrix */
301: MatDiagonalScale(a->B,NULL,auglyoo);
302: return(0);
303: }