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