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