Actual source code: mmaij.c
1: /*$Id: mmaij.c,v 1.59 2001/08/07 03:02:49 balay Exp $*/
3: /*
4: Support for the parallel AIJ matrix vector multiply
5: */
6: #include src/mat/impls/aij/mpi/mpiaij.h
7: #include src/vec/vecimpl.h
9: int MatSetUpMultiply_MPIAIJ(Mat mat)
10: {
11: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data;
12: Mat_SeqAIJ *B = (Mat_SeqAIJ*)(aij->B->data);
13: int N = mat->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
14: int shift = B->indexshift;
15: IS from,to;
16: Vec gvec;
17: #if defined (PETSC_USE_CTABLE)
18: PetscTable gid1_lid1;
19: PetscTablePosition tpos;
20: int gid,lid;
21: #endif
25: #if defined (PETSC_USE_CTABLE)
26: /* use a table - Mark Adams (this has not been tested with "shift") */
27: PetscTableCreate(aij->B->m,&gid1_lid1);
28: for (i=0; i<aij->B->m; i++) {
29: for (j=0; j<B->ilen[i]; j++) {
30: int data,gid1 = aj[B->i[i] + shift + j] + 1 + shift;
31: PetscTableFind(gid1_lid1,gid1,&data);
32: if (!data) {
33: /* one based table */
34: PetscTableAdd(gid1_lid1,gid1,++ec);
35: }
36: }
37: }
38: /* form array of columns we need */
39: PetscMalloc((ec+1)*sizeof(int),&garray);
40: PetscTableGetHeadPosition(gid1_lid1,&tpos);
41: while (tpos) {
42: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
43: gid--;
44: lid--;
45: garray[lid] = gid;
46: }
47: PetscSortInt(ec,garray); /* sort, and rebuild */
48: PetscTableRemoveAll(gid1_lid1);
49: for (i=0; i<ec; i++) {
50: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
51: }
52: /* compact out the extra columns in B */
53: for (i=0; i<aij->B->m; i++) {
54: for (j=0; j<B->ilen[i]; j++) {
55: int gid1 = aj[B->i[i] + shift + j] + 1 + shift;
56: PetscTableFind(gid1_lid1,gid1,&lid);
57: lid --;
58: aj[B->i[i] + shift + j] = lid - shift;
59: }
60: }
61: aij->B->n = aij->B->N = ec;
62: PetscTableDelete(gid1_lid1);
63: /* Mark Adams */
64: #else
65: /* For the first stab we make an array as long as the number of columns */
66: /* mark those columns that are in aij->B */
67: PetscMalloc((N+1)*sizeof(int),&indices);
68: PetscMemzero(indices,N*sizeof(int));
69: for (i=0; i<aij->B->m; i++) {
70: for (j=0; j<B->ilen[i]; j++) {
71: if (!indices[aj[B->i[i] +shift + j] + shift]) ec++;
72: indices[aj[B->i[i] + shift + j] + shift] = 1;
73: }
74: }
76: /* form array of columns we need */
77: PetscMalloc((ec+1)*sizeof(int),&garray);
78: ec = 0;
79: for (i=0; i<N; i++) {
80: if (indices[i]) garray[ec++] = i;
81: }
83: /* make indices now point into garray */
84: for (i=0; i<ec; i++) {
85: indices[garray[i]] = i-shift;
86: }
88: /* compact out the extra columns in B */
89: for (i=0; i<aij->B->m; i++) {
90: for (j=0; j<B->ilen[i]; j++) {
91: aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift];
92: }
93: }
94: aij->B->n = aij->B->N = ec;
95: PetscFree(indices);
96: #endif
97: /* create local vector that is used to scatter into */
98: VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);
100: /* create two temporary Index sets for build scatter gather */
101: ISCreateGeneral(mat->comm,ec,garray,&from);
102: ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);
104: /* create temporary global vector to generate scatter context */
105: /* this is inefficient, but otherwise we must do either
106: 1) save garray until the first actual scatter when the vector is known or
107: 2) have another way of generating a scatter context without a vector.*/
108: VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);
110: /* generate the scatter context */
111: VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);
112: PetscLogObjectParent(mat,aij->Mvctx);
113: PetscLogObjectParent(mat,aij->lvec);
114: PetscLogObjectParent(mat,from);
115: PetscLogObjectParent(mat,to);
116: aij->garray = garray;
117: PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
118: ISDestroy(from);
119: ISDestroy(to);
120: VecDestroy(gvec);
121: return(0);
122: }
125: /*
126: Takes the local part of an already assembled MPIAIJ matrix
127: and disassembles it. This is to allow new nonzeros into the matrix
128: that require more communication in the matrix vector multiply.
129: Thus certain data-structures must be rebuilt.
131: Kind of slow! But that's what application programmers get when
132: they are sloppy.
133: */
134: int DisAssemble_MPIAIJ(Mat A)
135: {
136: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data;
137: Mat B = aij->B,Bnew;
138: Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data;
139: int ierr,i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray;
140: int *nz,ec,shift = Baij->indexshift;
141: PetscScalar v;
144: /* free stuff related to matrix-vec multiply */
145: VecGetSize(aij->lvec,&ec); /* needed for PetscLogObjectMemory below */
146: VecDestroy(aij->lvec); aij->lvec = 0;
147: VecScatterDestroy(aij->Mvctx); aij->Mvctx = 0;
148: if (aij->colmap) {
149: #if defined (PETSC_USE_CTABLE)
150: PetscTableDelete(aij->colmap);
151: aij->colmap = 0;
152: #else
153: PetscFree(aij->colmap);
154: aij->colmap = 0;
155: PetscLogObjectMemory(A,-aij->B->n*sizeof(int));
156: #endif
157: }
159: /* make sure that B is assembled so we can access its values */
160: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
161: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
163: /* invent new B and copy stuff over */
164: PetscMalloc((m+1)*sizeof(int),&nz);
165: for (i=0; i<m; i++) {
166: nz[i] = Baij->i[i+1] - Baij->i[i];
167: }
168: MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew);
169: PetscFree(nz);
170: for (i=0; i<m; i++) {
171: for (j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++) {
172: col = garray[Baij->j[ct]+shift];
173: v = Baij->a[ct++];
174: MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);
175: }
176: }
177: PetscFree(aij->garray);
178: aij->garray = 0;
179: PetscLogObjectMemory(A,-ec*sizeof(int));
180: MatDestroy(B);
181: PetscLogObjectParent(A,Bnew);
182: aij->B = Bnew;
183: A->was_assembled = PETSC_FALSE;
184: return(0);
185: }
187: /* ugly stuff added for Glenn someday we should fix this up */
189: static int *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal"
190: parts of the local matrix */
191: static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */
194: int MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
195: {
196: Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
197: int ierr,i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
198: int *r_rmapd,*r_rmapo;
199:
201: MatGetOwnershipRange(inA,&cstart,&cend);
202: MatGetSize(ina->A,PETSC_NULL,&n);
203: PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);
204: PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));
205: nt = 0;
206: for (i=0; i<inA->mapping->n; i++) {
207: if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) {
208: nt++;
209: r_rmapd[i] = inA->mapping->indices[i] + 1;
210: }
211: }
212: if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n);
213: PetscMalloc((n+1)*sizeof(int),&auglyrmapd);
214: for (i=0; i<inA->mapping->n; i++) {
215: if (r_rmapd[i]){
216: auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
217: }
218: }
219: PetscFree(r_rmapd);
220: VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);
222: PetscMalloc((inA->N+1)*sizeof(int),&lindices);
223: PetscMemzero(lindices,inA->N*sizeof(int));
224: for (i=0; i<ina->B->n; i++) {
225: lindices[garray[i]] = i+1;
226: }
227: no = inA->mapping->n - nt;
228: PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);
229: PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));
230: nt = 0;
231: for (i=0; i<inA->mapping->n; i++) {
232: if (lindices[inA->mapping->indices[i]]) {
233: nt++;
234: r_rmapo[i] = lindices[inA->mapping->indices[i]];
235: }
236: }
237: if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n);
238: PetscFree(lindices);
239: PetscMalloc((nt+1)*sizeof(int),&auglyrmapo);
240: for (i=0; i<inA->mapping->n; i++) {
241: if (r_rmapo[i]){
242: auglyrmapo[(r_rmapo[i]-1)] = i;
243: }
244: }
245: PetscFree(r_rmapo);
246: VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);
248: return(0);
249: }
251: int MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
252: {
253: Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
254: int ierr,n,i;
255: PetscScalar *d,*o,*s;
256:
258: if (!auglyrmapd) {
259: MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
260: }
262: VecGetArray(scale,&s);
263:
264: VecGetLocalSize(auglydd,&n);
265: VecGetArray(auglydd,&d);
266: for (i=0; i<n; i++) {
267: d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
268: }
269: VecRestoreArray(auglydd,&d);
270: /* column scale "diagonal" portion of local matrix */
271: MatDiagonalScale(a->A,PETSC_NULL,auglydd);
273: VecGetLocalSize(auglyoo,&n);
274: VecGetArray(auglyoo,&o);
275: for (i=0; i<n; i++) {
276: o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
277: }
278: VecRestoreArray(scale,&s);
279: VecRestoreArray(auglyoo,&o);
280: /* column scale "off-diagonal" portion of local matrix */
281: MatDiagonalScale(a->B,PETSC_NULL,auglyoo);
283: return(0);
284: }