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