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