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