xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision 4ccb76dd2a06ade17c2800eec40bc3f6b9abd66f)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: mmaij.c,v 1.46 1999/05/12 03:29:11 bsmith Exp balay $";
3 #endif
4 
5 
6 /*
7    Support for the parallel AIJ matrix vector multiply
8 */
9 #include "src/mat/impls/aij/mpi/mpiaij.h"
10 #include "src/vec/vecimpl.h"
11 
12 #undef __FUNC__
13 #define __FUNC__ "MatSetUpMultiply_MPIAIJ"
14 int MatSetUpMultiply_MPIAIJ(Mat mat)
15 {
16   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
17   Mat_SeqAIJ *B = (Mat_SeqAIJ *) (aij->B->data);
18   int        N = aij->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
19   int        shift = B->indexshift;
20   IS         from,to;
21   Vec        gvec;
22 #if defined (PETSC_USE_CTABLE)
23   Table      gid1_lid1;
24   CTablePos  tpos;
25   int        gid, lid;
26 #endif
27 
28   PetscFunctionBegin;
29 
30 #if defined (PETSC_USE_CTABLE)
31   /* use a table - Mark Adams (this has not been tested with "shift") */
32   TableCreate(B->m,&gid1_lid1);
33   for ( i=0; i<B->m; i++ ) {
34     for ( j=0; j<B->ilen[i]; j++ ) {
35       int data,gid1 = aj[B->i[i] + shift + j] + 1 + shift;
36       ierr = TableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
37       if (!data) {
38         /* one based table */
39         ierr = TableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr);
40       }
41     }
42   }
43   /* form array of columns we need */
44   garray = (int *)PetscMalloc((ec+1)*sizeof(int));CHKPTRQ(garray);
45   ierr = TableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
46   while (tpos) {
47     ierr = TableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
48     gid--; lid--;
49     garray[lid] = gid;
50   }
51   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */
52   /* qsort( garray, ec, sizeof(int), intcomparc ); */
53   TableRemoveAll(gid1_lid1);
54   for ( i=0; i<ec; i++ ) {
55     ierr = TableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr);
56   }
57   /* compact out the extra columns in B */
58   for ( i=0; i<B->m; i++ ) {
59     for ( j=0; j<B->ilen[i]; j++ ) {
60       int gid1 = aj[B->i[i] + shift + j] + 1 + shift;
61       ierr = TableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
62       lid --;
63       aj[B->i[i] + shift + j]  = lid - shift;
64     }
65   }
66   B->n = aij->B->n = aij->B->N = ec;
67   TableDelete(gid1_lid1);
68   /* Mark Adams */
69 #else
70   /* For the first stab we make an array as long as the number of columns */
71   /* mark those columns that are in aij->B */
72   indices = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(indices);
73   ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr);
74   for ( i=0; i<B->m; i++ ) {
75     for ( j=0; j<B->ilen[i]; j++ ) {
76       if (!indices[aj[B->i[i] +shift + j] + shift]) ec++;
77       indices[aj[B->i[i] + shift + j] + shift] = 1;
78     }
79   }
80 
81   /* form array of columns we need */
82   garray = (int *) PetscMalloc( (ec+1)*sizeof(int) );CHKPTRQ(garray);
83   ec = 0;
84   for ( i=0; i<N; i++ ) {
85     if (indices[i]) garray[ec++] = i;
86   }
87 
88   /* make indices now point into garray */
89   for ( i=0; i<ec; i++ ) {
90     indices[garray[i]] = i-shift;
91   }
92 
93   /* compact out the extra columns in B */
94   for ( i=0; i<B->m; i++ ) {
95     for ( j=0; j<B->ilen[i]; j++ ) {
96       aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift];
97     }
98   }
99   B->n = aij->B->n = aij->B->N = ec;
100   PetscFree(indices);
101 #endif
102   /* create local vector that is used to scatter into */
103   ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr);
104 
105   /* create two temporary Index sets for build scatter gather */
106   ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,&from);CHKERRQ(ierr);
107   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
108 
109   /* create temporary global vector to generate scatter context */
110   /* this is inefficient, but otherwise we must do either
111      1) save garray until the first actual scatter when the vector is known or
112      2) have another way of generating a scatter context without a vector.*/
113   ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec);CHKERRQ(ierr);
114 
115   /* generate the scatter context */
116   ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr);
117   PLogObjectParent(mat,aij->Mvctx);
118   PLogObjectParent(mat,aij->lvec);
119   PLogObjectParent(mat,from);
120   PLogObjectParent(mat,to);
121   aij->garray = garray;
122   PLogObjectMemory(mat,(ec+1)*sizeof(int));
123   ierr = ISDestroy(from);CHKERRQ(ierr);
124   ierr = ISDestroy(to);CHKERRQ(ierr);
125   ierr = VecDestroy(gvec);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 
130 #undef __FUNC__
131 #define __FUNC__ "DisAssemble_MPIAIJ"
132 /*
133      Takes the local part of an already assembled MPIAIJ matrix
134    and disassembles it. This is to allow new nonzeros into the matrix
135    that require more communication in the matrix vector multiply.
136    Thus certain data-structures must be rebuilt.
137 
138    Kind of slow! But that's what application programmers get when
139    they are sloppy.
140 */
141 int DisAssemble_MPIAIJ(Mat A)
142 {
143   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data;
144   Mat        B = aij->B,Bnew;
145   Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data;
146   int        ierr,i,j,m = Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray;
147   int        *nz,ec,shift = Baij->indexshift;
148   Scalar     v;
149 
150   PetscFunctionBegin;
151   /* free stuff related to matrix-vec multiply */
152   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PLogObjectMemory below */
153   ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0;
154   ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0;
155   if (aij->colmap) {
156 #if defined (PETSC_USE_CTABLE)
157     TableDelete(aij->colmap); aij->colmap = 0;
158 #else
159     PetscFree(aij->colmap); aij->colmap = 0;
160     PLogObjectMemory(A,-Baij->n*sizeof(int));
161 #endif
162   }
163 
164   /* make sure that B is assembled so we can access its values */
165   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167 
168   /* invent new B and copy stuff over */
169   nz = (int *) PetscMalloc( (m+1)*sizeof(int) );CHKPTRQ(nz);
170   for ( i=0; i<m; i++ ) {
171     nz[i] = Baij->i[i+1] - Baij->i[i];
172   }
173   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew);CHKERRQ(ierr);
174   PetscFree(nz);
175   for ( i=0; i<m; i++ ) {
176     for ( j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++ ) {
177       col = garray[Baij->j[ct]+shift];
178       v = Baij->a[ct++];
179       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr);
180     }
181   }
182   PetscFree(aij->garray); aij->garray = 0;
183   PLogObjectMemory(A,-ec*sizeof(int));
184   ierr = MatDestroy(B);CHKERRQ(ierr);
185   PLogObjectParent(A,Bnew);
186   aij->B = Bnew;
187   A->was_assembled = PETSC_FALSE;
188   PetscFunctionReturn(0);
189 }
190 
191 
192