Lines Matching +full:- +full:r

3           C = R * A * R^T
15 PetscCall(MatTransposeColoringDestroy(&rart->matcoloring));
16 PetscCall(MatDestroy(&rart->Rt));
17 PetscCall(MatDestroy(&rart->RARt));
18 PetscCall(MatDestroy(&rart->ARt));
19 PetscCall(PetscFree(rart->work));
20 if (rart->destroy) PetscCall((*rart->destroy)(&rart->data));
25 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A, Mat R, PetscReal fill, Mat C)
36 PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
38 PetscCall(MatTransposeSymbolic(R, &P));
42 PetscCall(MatSetBlockSizes(C, R->rmap->bs, R->rmap->bs));
43 C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart;
47 C->product->data = rart;
48 C->product->destroy = MatProductCtxDestroy_SeqAIJ_RARt;
54 /* Create MatTransposeColoring from symbolic C=R*A*R^T */
63 rart->matcoloring = matcoloring;
68 PetscCall(MatSetSizes(Rt_dense, A->cmap->n, matcoloring->ncolors, A->cmap->n, matcoloring->ncolors));
72 Rt_dense->assembled = PETSC_TRUE;
73 rart->Rt = Rt_dense;
75 /* Create RARt_dense = R*A*Rt_dense */
77 PetscCall(MatSetSizes(RARt_dense, C->rmap->n, matcoloring->ncolors, C->rmap->n, matcoloring->ncolors));
81 rart->RARt = RARt_dense;
83 /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
84 PetscCall(PetscMalloc1(A->rmap->n * 4, &rart->work));
91 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
92 PetscReal density = (PetscReal)c->nz / (RARt_dense->rmap->n * RARt_dense->cmap->n);
93 PetscCall(PetscInfo(C, "C=R*(A*Rt) via coloring C - use sparse-dense inner products\n"));
94 PetscCall(PetscInfo(C, "RARt_den %" PetscInt_FMT " %" PetscInt_FMT "; Rt %" PetscInt_FMT " %" PetscInt_FMT " (RARt->nz %" PetscInt_FMT ")/(m*ncolors)=%g\n", RARt_dense->rmap->n, RARt_dense->cmap->n, R->cmap->n, R->rmap->n, c->nz, (double)density));
101 RAB = R * A * B, R and A in seqaij format, B in dense format;
103 static PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R, Mat A, Mat B, Mat RAB, PetscScalar *work)
105 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *r = (Mat_SeqAIJ *)R->data;
109 PetscInt cn = B->cmap->n, bm = B->rmap->n, col, i, j, n, *ai = a->i, *aj, am = A->rmap->n;
112 PetscInt *rj, rm = R->rmap->n, dm = RAB->rmap->n, dn = RAB->cmap->n;
117 PetscCheck(bm == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number columns in A %" PetscInt_FMT " not equal rows in B %" PetscInt_FMT, A->cmap->n, bm);
118 PetscCheck(am == R->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number columns in R %" PetscInt_FMT " not equal rows in A %" PetscInt_FMT, R->cmap->n, am);
119 PetscCheck(R->rmap->n == RAB->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number rows in RAB %" PetscInt_FMT " not equal rows in R %" PetscInt_FMT, RAB->rmap->n, R->rmap->n);
120 PetscCheck(B->cmap->n == RAB->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number columns in RAB %" PetscInt_FMT " not equal columns in B %" PetscInt_FMT, RAB->cmap->n, B->cmap->n);
124 AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/tutorials/ex56.c
127 PetscCall(PetscOptionsGetBool(NULL, NULL, "-matrart_via_matmatmult", &via_matmatmult, NULL));
133 PetscCall(MatMatMultNumeric_SeqAIJ_SeqDense(R, AB_den, RAB));
149 for (col = 0; col < cn - 4; col += 4) { /* over columns of C */
152 n = ai[i + 1] - ai[i];
153 aj = a->j + ai[i];
154 aa = a->a + ai[i];
171 /* RAB[:,col] = R*C[:,col] */
173 for (i = 0; i < rm; i++) { /* over rows of R in those columns */
175 n = r->i[i + 1] - r->i[i];
176 rj = r->j + r->i[i];
177 ra = r->a + r->i[i];
193 n = a->i[i + 1] - a->i[i];
194 aj = a->j + a->i[i];
195 aa = a->a + a->i[i];
201 for (i = 0; i < rm; i++) { /* over rows of R in those columns */
203 n = r->i[i + 1] - r->i[i];
204 rj = r->j + r->i[i];
205 ra = r->a + r->i[i];
210 PetscCall(PetscLogFlops(cn * 2.0 * (a->nz + r->nz)));
219 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A, Mat R, Mat C)
227 PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
228 rart = (MatProductCtx_RARt *)C->product->data;
230 /* Get dense Rt by Apply MatTransposeColoring to R */
231 matcoloring = rart->matcoloring;
232 Rt = rart->Rt;
233 PetscCall(MatTransColoringApplySpToDen(matcoloring, R, Rt));
235 /* Get dense RARt = R*A*Rt -- dominates! */
236 RARt = rart->RARt;
237 PetscCall(MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R, A, Rt, RARt, rart->work));
244 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A, Mat R, PetscReal fill, Mat C)
252 PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
253 /* create symbolic ARt = A*R^T */
254 PetscCall(MatProductCreate(A, R, NULL, &ARt));
261 /* compute symbolic C = R*ARt */
262 /* set algorithm for C = R*ARt */
263 PetscCall(PetscStrallocpy(C->product->alg, &alg));
265 PetscCall(MatMatMultSymbolic_SeqAIJ_SeqAIJ(R, ARt, fill, C));
270 C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult;
273 rart->ARt = ARt;
274 C->product->data = rart;
275 C->product->destroy = MatProductCtxDestroy_SeqAIJ_RARt;
276 PetscCall(PetscInfo(C, "Use ARt=A*R^T, C=R*ARt via MatMatTransposeMult(). Coloring can be applied to A*R^T.\n"));
280 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A, Mat R, Mat C)
286 PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
287 rart = (MatProductCtx_RARt *)C->product->data;
288 PetscCall(MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A, R, rart->ARt)); /* dominate! */
289 PetscCall(MatMatMultNumeric_SeqAIJ_SeqAIJ(R, rart->ARt, C));
293 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A, Mat R, PetscReal fill, Mat C)
300 PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
301 PetscCall(MatTranspose(R, MAT_INITIAL_MATRIX, &Rt));
302 PetscCall(MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R, A, Rt, fill, C));
305 rart->data = C->product->data;
306 rart->destroy = C->product->destroy;
307 rart->Rt = Rt;
308 C->product->data = rart;
309 C->product->destroy = MatProductCtxDestroy_SeqAIJ_RARt;
310 C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ;
311 PetscCall(PetscInfo(C, "Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n"));
315 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A, Mat R, Mat C)
321 PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
322 rart = (MatProductCtx_RARt *)C->product->data;
323 PetscCall(MatTranspose(R, MAT_REUSE_MATRIX, &rart->Rt));
325 C->product->data = rart->data;
326 PetscCall(MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R, A, rart->Rt, C));
327 C->product->data = rart;
331 PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A, Mat R, MatReuse scall, PetscReal fill, Mat *C)
338 PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "MatRARt", "Mat");
339 PetscCall(PetscOptionsEList("-matrart_via", "Algorithmic approach", "MatRARt", algTypes, 3, algTypes[0], &alg, NULL));
342 PetscCall(PetscLogEventBegin(MAT_RARtSymbolic, A, R, 0, 0));
346 /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */
347 PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A, R, fill, *C));
350 /* via coloring_rart: apply coloring C = R*A*R^T */
351 PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A, R, fill, *C));
354 /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */
355 PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ(A, R, fill, *C));
358 PetscCall(PetscLogEventEnd(MAT_RARtSymbolic, A, R, 0, 0));
361 PetscCall(PetscLogEventBegin(MAT_RARtNumeric, A, R, 0, 0));
362 PetscCall(((*C)->ops->rartnumeric)(A, R, *C));
363 PetscCall(PetscLogEventEnd(MAT_RARtNumeric, A, R, 0, 0));
369 Mat_Product *product = C->product;
370 Mat A = product->A, R = product->B;
371 MatProductAlgorithm alg = product->alg;
372 PetscReal fill = product->fill;
376 PetscCall(PetscStrcmp(alg, "r*a*rt", &flg));
378 PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ(A, R, fill, C));
382 PetscCall(PetscStrcmp(alg, "r*art", &flg));
384 PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A, R, fill, C));
390 PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A, R, fill, C));
397 C->ops->productnumeric = MatProductNumeric_RARt;