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