xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision 95d5f7c29374efcfd2ca44c2fe93981fbc2b4454)
1 /*$Id: mmaij.c,v 1.52 2000/04/09 04:36:04 bsmith Exp bsmith $*/
2 
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"
8 
9 #undef __FUNC__
10 #define __FUNC__ /*<a name=""></a>*/"MatSetUpMultiply_MPIAIJ"
11 int MatSetUpMultiply_MPIAIJ(Mat mat)
12 {
13   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
14   Mat_SeqAIJ         *B = (Mat_SeqAIJ*)(aij->B->data);
15   int                N = aij->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
16   int                shift = B->indexshift;
17   IS                 from,to;
18   Vec                gvec;
19 #if defined (PETSC_USE_CTABLE)
20   PetscTable         gid1_lid1;
21   PetscTablePosition tpos;
22   int                gid,lid;
23 #endif
24 
25   PetscFunctionBegin;
26 
27 #if defined (PETSC_USE_CTABLE)
28   /* use a table - Mark Adams (this has not been tested with "shift") */
29   PetscTableCreate(B->m,&gid1_lid1);
30   for (i=0; i<B->m; i++) {
31     for (j=0; j<B->ilen[i]; j++) {
32       int data,gid1 = aj[B->i[i] + shift + j] + 1 + shift;
33       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
34       if (!data) {
35         /* one based table */
36         ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr);
37       }
38     }
39   }
40   /* form array of columns we need */
41   garray = (int *)PetscMalloc((ec+1)*sizeof(int));CHKPTRQ(garray);
42   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
43   while (tpos) {
44     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
45     gid--; lid--;
46     garray[lid] = gid;
47   }
48   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */
49   /* qsort(garray, ec, sizeof(int), intcomparc); */
50   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
51   for (i=0; i<ec; i++) {
52     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr);
53   }
54   /* compact out the extra columns in B */
55   for (i=0; i<B->m; i++) {
56     for (j=0; j<B->ilen[i]; j++) {
57       int gid1 = aj[B->i[i] + shift + j] + 1 + shift;
58       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
59       lid --;
60       aj[B->i[i] + shift + j]  = lid - shift;
61     }
62   }
63   B->n = aij->B->n = aij->B->N = ec;
64   ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr);
65   /* Mark Adams */
66 #else
67   /* For the first stab we make an array as long as the number of columns */
68   /* mark those columns that are in aij->B */
69   indices = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(indices);
70   ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr);
71   for (i=0; i<B->m; i++) {
72     for (j=0; j<B->ilen[i]; j++) {
73       if (!indices[aj[B->i[i] +shift + j] + shift]) ec++;
74       indices[aj[B->i[i] + shift + j] + shift] = 1;
75     }
76   }
77 
78   /* form array of columns we need */
79   garray = (int*)PetscMalloc((ec+1)*sizeof(int));CHKPTRQ(garray);
80   ec = 0;
81   for (i=0; i<N; i++) {
82     if (indices[i]) garray[ec++] = i;
83   }
84 
85   /* make indices now point into garray */
86   for (i=0; i<ec; i++) {
87     indices[garray[i]] = i-shift;
88   }
89 
90   /* compact out the extra columns in B */
91   for (i=0; i<B->m; i++) {
92     for (j=0; j<B->ilen[i]; j++) {
93       aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift];
94     }
95   }
96   B->n = aij->B->n = aij->B->N = ec;
97   ierr = PetscFree(indices);CHKERRQ(ierr);
98 #endif
99   /* create local vector that is used to scatter into */
100   ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr);
101 
102   /* create two temporary Index sets for build scatter gather */
103   ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,&from);CHKERRQ(ierr);
104   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
105 
106   /* create temporary global vector to generate scatter context */
107   /* this is inefficient, but otherwise we must do either
108      1) save garray until the first actual scatter when the vector is known or
109      2) have another way of generating a scatter context without a vector.*/
110   ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec);CHKERRQ(ierr);
111 
112   /* generate the scatter context */
113   ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr);
114   PLogObjectParent(mat,aij->Mvctx);
115   PLogObjectParent(mat,aij->lvec);
116   PLogObjectParent(mat,from);
117   PLogObjectParent(mat,to);
118   aij->garray = garray;
119   PLogObjectMemory(mat,(ec+1)*sizeof(int));
120   ierr = ISDestroy(from);CHKERRQ(ierr);
121   ierr = ISDestroy(to);CHKERRQ(ierr);
122   ierr = VecDestroy(gvec);CHKERRQ(ierr);
123   PetscFunctionReturn(0);
124 }
125 
126 
127 #undef __FUNC__
128 #define __FUNC__ /*<a name=""></a>*/"DisAssemble_MPIAIJ"
129 /*
130      Takes the local part of an already assembled MPIAIJ matrix
131    and disassembles it. This is to allow new nonzeros into the matrix
132    that require more communication in the matrix vector multiply.
133    Thus certain data-structures must be rebuilt.
134 
135    Kind of slow! But that's what application programmers get when
136    they are sloppy.
137 */
138 int DisAssemble_MPIAIJ(Mat A)
139 {
140   Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data;
141   Mat        B = aij->B,Bnew;
142   Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data;
143   int        ierr,i,j,m = Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray;
144   int        *nz,ec,shift = Baij->indexshift;
145   Scalar     v;
146 
147   PetscFunctionBegin;
148   /* free stuff related to matrix-vec multiply */
149   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PLogObjectMemory below */
150   ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0;
151   ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0;
152   if (aij->colmap) {
153 #if defined (PETSC_USE_CTABLE)
154     ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr);
155     aij->colmap = 0;
156 #else
157     ierr = PetscFree(aij->colmap);CHKERRQ(ierr);
158     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   ierr = PetscFree(nz);CHKERRQ(ierr);
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   ierr = PetscFree(aij->garray);CHKERRQ(ierr);
182   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