xref: /petsc/src/mat/impls/aij/mpi/crl/mcrl.c (revision 183aafc7325fac25f4e3f5d49e5dca787cc2dd63)
1 #define PETSCMAT_DLL
2 
3 /*
4   Defines a matrix-vector product for the MATMPIAIJCRL matrix class.
5   This class is derived from the MATMPIAIJ class and retains the
6   compressed row storage (aka Yale sparse matrix format) but augments
7   it with a column oriented storage that is more efficient for
8   matrix vector products on Vector machines.
9 
10   CRL stands for constant row length (that is the same number of columns
11   is kept (padded with zeros) for each row of the sparse matrix.
12 
13    See src/mat/impls/aij/seq/crl/crl.c for the sequential version
14 */
15 
16 #include "src/mat/impls/aij/mpi/mpiaij.h"
17 #include "src/mat/impls/aij/seq/crl/crl.h"
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "MatDestroy_MPICRL"
21 PetscErrorCode MatDestroy_MPICRL(Mat A)
22 {
23   PetscErrorCode ierr;
24   Mat_CRL        *crl = (Mat_CRL *) A->spptr;
25 
26   /* We are going to convert A back into a MPIAIJ matrix, since we are
27    * eventually going to use MatDestroy_MPIAIJ() to destroy everything
28    * that is not specific to CRL.
29    * In preparation for this, reset the operations pointers in A to
30    * their MPIAIJ versions. */
31   A->ops->assemblyend = crl->AssemblyEnd;
32   A->ops->destroy     = crl->MatDestroy;
33   A->ops->duplicate   = crl->MatDuplicate;
34 
35   /* Free everything in the Mat_CRL data structure. */
36   if (crl->icols) {
37     ierr = PetscFree2(crl->acols,crl->icols);CHKERRQ(ierr);
38   }
39   /* Free the Mat_CRL struct itself. */
40   ierr = PetscFree(crl);CHKERRQ(ierr);
41 
42   /* Change the type of A back to MPIAIJ and use MatDestroy_MPIAIJ()
43    * to destroy everything that remains. */
44   ierr = PetscObjectChangeTypeName( (PetscObject)A, MATMPIAIJ);CHKERRQ(ierr);
45   /* Note that I don't call MatSetType().  I believe this is because that
46    * is only to be called when *building* a matrix. */
47   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
48   PetscFunctionReturn(0);
49 }
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "MPICRL_create_crl"
53 PetscErrorCode MPICRL_create_crl(Mat A)
54 {
55   Mat_MPIAIJ     *a = (Mat_MPIAIJ *)(A)->data;
56   Mat_SeqAIJ     *Aij = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->A->data);
57   Mat_CRL        *crl = (Mat_CRL*) A->spptr;
58   PetscInt       m = A->m;  /* Number of rows in the matrix. */
59   PetscInt       *aj = Aij->j,*bj = Bij->j;  /* From the CSR representation; points to the beginning  of each row. */
60   PetscInt       i, j,rmax = 0,*icols, *ailen = Aij->ilen, *bilen = Bij->ilen;
61   PetscScalar    *aa = Aij->a,*ba = Bij->a,*acols,*array;
62   PetscErrorCode ierr;
63 
64   PetscFunctionBegin;
65   /* determine the row with the most columns */
66   for (i=0; i<m; i++) {
67     rmax = PetscMax(rmax,ailen[i]+bilen[i]);
68   }
69   crl->nz   = Aij->nz+Bij->nz;
70   crl->m    = A->m;
71   crl->rmax = rmax;
72   ierr  = PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);CHKERRQ(ierr);
73   acols = crl->acols;
74   icols = crl->icols;
75   for (i=0; i<m; i++) {
76     for (j=0; j<ailen[i]; j++) {
77       acols[j*m+i] = *aa++;
78       icols[j*m+i] = *aj++;
79     }
80     for (;j<ailen[i]+bilen[i]; j++) {
81       acols[j*m+i] = *ba++;
82       icols[j*m+i] = *bj++;
83     }
84     for (;j<rmax; j++) { /* empty column entries */
85       acols[j*m+i] = 0.0;
86       icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0;  /* handle case where row is EMPTY */
87     }
88   }
89   ierr = PetscLogInfo((A,"MPICRL_create_crl: Percentage of 0's introduced for vectorized multiply %g\n",1.0-((double)(Bij->nz + Bij->nz))/((double)(rmax*m))));
90 
91   ierr = VecCreateSeqWithArray(A->comm,Aij->rmax,array,&crl->xwork);CHKERRQ(ierr);
92 
93   PetscFunctionReturn(0);
94 }
95 
96 #undef __FUNCT__
97 #define __FUNCT__ "MatAssemblyEnd_MPICRL"
98 PetscErrorCode MatAssemblyEnd_MPICRL(Mat A, MatAssemblyType mode)
99 {
100   PetscErrorCode ierr;
101   Mat_CRL        *crl = (Mat_CRL*) A->spptr;
102   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
103   Mat_SeqAIJ     *Aij = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->A->data);
104 
105   PetscFunctionBegin;
106   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
107 
108   /* Since a MATMPICRL matrix is really just a MATMPIAIJ with some
109    * extra information, call the AssemblyEnd routine for a MATMPIAIJ.
110    * I'm not sure if this is the best way to do this, but it avoids
111    * a lot of code duplication.
112    * I also note that currently MATMPICRL doesn't know anything about
113    * the Mat_CompressedRow data structure that MPIAIJ now uses when there
114    * are many zero rows.  If the MPIAIJ assembly end routine decides to use
115    * this, this may break things.  (Don't know... haven't looked at it.) */
116   Aij->inode.use = PETSC_FALSE;
117   Bij->inode.use = PETSC_FALSE;
118   (*crl->AssemblyEnd)(A, mode);
119 
120   /* Now calculate the permutation and grouping information. */
121   ierr = MPICRL_create_crl(A);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 extern PetscErrorCode MatMult_CRL(Mat,Vec,Vec);
126 extern PetscErrorCode MatDuplicate_CRL(Mat,MatDuplicateOption,Mat*);
127 
128 /* MatConvert_MPIAIJ_MPICRL converts a MPIAIJ matrix into a
129  * MPICRL matrix.  This routine is called by the MatCreate_MPICRL()
130  * routine, but can also be used to convert an assembled MPIAIJ matrix
131  * into a MPICRL one. */
132 EXTERN_C_BEGIN
133 #undef __FUNCT__
134 #define __FUNCT__ "MatConvert_MPIAIJ_MPICRL"
135 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MPIAIJ_MPICRL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
136 {
137   /* This routine is only called to convert to MATMPICRL
138    * from MATMPIAIJ, so we can ignore 'MatType Type'. */
139   PetscErrorCode ierr;
140   Mat            B = *newmat;
141   Mat_CRL        *crl;
142 
143   PetscFunctionBegin;
144   if (reuse == MAT_INITIAL_MATRIX) {
145     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
146   }
147 
148   ierr = PetscNew(Mat_CRL,&crl);CHKERRQ(ierr);
149   B->spptr = (void *) crl;
150 
151   /* Save a pointer to the original MPIAIJ assembly end routine, because we
152    * will want to use it later in the CRL assembly end routine.
153    * Also, save a pointer to the original MPIAIJ Destroy routine, because we
154    * will want to use it in the CRL destroy routine. */
155   crl->AssemblyEnd  = A->ops->assemblyend;
156   crl->MatDestroy   = A->ops->destroy;
157   crl->MatDuplicate = A->ops->duplicate;
158 
159   /* Set function pointers for methods that we inherit from AIJ but
160    * override. */
161   B->ops->duplicate   = MatDuplicate_CRL;
162   B->ops->assemblyend = MatAssemblyEnd_MPICRL;
163   B->ops->destroy     = MatDestroy_MPICRL;
164   B->ops->mult        = MatMult_CRL;
165 
166   /* If A has already been assembled, compute the permutation. */
167   if (A->assembled == PETSC_TRUE) {
168     ierr = MPICRL_create_crl(B);CHKERRQ(ierr);
169   }
170   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPICRL);CHKERRQ(ierr);
171   *newmat = B;
172   PetscFunctionReturn(0);
173 }
174 EXTERN_C_END
175 
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "MatCreateMPICRL"
179 /*@C
180    MatCreateMPICRL - Creates a sparse matrix of type MPICRL.
181    This type inherits from AIJ, but stores some additional
182    information that is used to allow better vectorization of
183    the matrix-vector product. At the cost of increased storage, the AIJ formatted
184    matrix can be copied to a format in which pieces of the matrix are
185    stored in ELLPACK format, allowing the vectorized matrix multiply
186    routine to use stride-1 memory accesses.  As with the AIJ type, it is
187    important to preallocate matrix storage in order to get good assembly
188    performance.
189 
190    Collective on MPI_Comm
191 
192    Input Parameters:
193 +  comm - MPI communicator, set to PETSC_COMM_SELF
194 .  m - number of rows
195 .  n - number of columns
196 .  nz - number of nonzeros per row (same for all rows)
197 -  nnz - array containing the number of nonzeros in the various rows
198          (possibly different for each row) or PETSC_NULL
199 
200    Output Parameter:
201 .  A - the matrix
202 
203    Notes:
204    If nnz is given then nz is ignored
205 
206    Level: intermediate
207 
208 .keywords: matrix, cray, sparse, parallel
209 
210 .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues()
211 @*/
212 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[],Mat *A)
213 {
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   ierr = MatCreate(comm,A);CHKERRQ(ierr);
218   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
219   ierr = MatSetType(*A,MATMPICRL);CHKERRQ(ierr);
220   ierr = MatMPIAIJSetPreallocation_MPIAIJ(*A,nz,(PetscInt*)nnz,onz,(PetscInt*)onnz);CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 
225 EXTERN_C_BEGIN
226 #undef __FUNCT__
227 #define __FUNCT__ "MatCreate_MPICRL"
228 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPICRL(Mat A)
229 {
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   /* Change the type name before calling MatSetType() to force proper construction of MPIAIJ
234      and MATMPICRL types. */
235   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMPICRL);CHKERRQ(ierr);
236   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
237   ierr = MatConvert_MPIAIJ_MPICRL(A,MATMPICRL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
238   PetscFunctionReturn(0);
239 }
240 EXTERN_C_END
241 
242