127d8b95cSPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/ 2c8a8475eSBarry Smith 3c8a8475eSBarry Smith typedef struct { 47cfd0b8cSBarry Smith Mat A; 5a48507d3SJose E. Roman Mat D; /* local submatrix for diagonal part */ 627d8b95cSPierre Jolivet Vec w; 7c8a8475eSBarry Smith } Mat_Normal; 8c8a8475eSBarry Smith 966976f2fSJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Normal(Mat A, PetscInt is_max, IS is[], PetscInt ov) 10d71ae5a4SJacob Faibussowitsch { 1127d8b95cSPierre Jolivet Mat_Normal *a; 12fa80d070SPierre Jolivet Mat pattern; 13fa80d070SPierre Jolivet 14fa80d070SPierre Jolivet PetscFunctionBegin; 1508401ef6SPierre Jolivet PetscCheck(ov >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified"); 1627d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 179566063dSJacob Faibussowitsch PetscCall(MatProductCreate(a->A, a->A, NULL, &pattern)); 189566063dSJacob Faibussowitsch PetscCall(MatProductSetType(pattern, MATPRODUCT_AtB)); 199566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(pattern)); 209566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(pattern)); 219566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(pattern, is_max, is, ov)); 229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pattern)); 233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24fa80d070SPierre Jolivet } 25fa80d070SPierre Jolivet 2666976f2fSJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Normal(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 27d71ae5a4SJacob Faibussowitsch { 2827d8b95cSPierre Jolivet Mat_Normal *a; 2927d8b95cSPierre Jolivet Mat B, *suba; 30fa80d070SPierre Jolivet IS *row; 31*b9c875b8SPierre Jolivet PetscScalar shift, scale; 32fa80d070SPierre Jolivet PetscInt M; 33fa80d070SPierre Jolivet 34fa80d070SPierre Jolivet PetscFunctionBegin; 3527d8b95cSPierre Jolivet PetscCheck(irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented"); 36*b9c875b8SPierre Jolivet PetscCall(MatShellGetScalingShifts(mat, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 3727d8b95cSPierre Jolivet PetscCall(MatShellGetContext(mat, &a)); 3827d8b95cSPierre Jolivet B = a->A; 3948a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat)); 409566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, &M, NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &row)); 429566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0])); 439566063dSJacob Faibussowitsch PetscCall(ISSetIdentity(row[0])); 44fa80d070SPierre Jolivet for (M = 1; M < n; ++M) row[M] = row[0]; 459566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba)); 46fa80d070SPierre Jolivet for (M = 0; M < n; ++M) { 479566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(suba[M], *submat + M)); 48*b9c875b8SPierre Jolivet PetscCall(MatShift((*submat)[M], shift)); 49*b9c875b8SPierre Jolivet PetscCall(MatScale((*submat)[M], scale)); 50fa80d070SPierre Jolivet } 519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row[0])); 529566063dSJacob Faibussowitsch PetscCall(PetscFree(row)); 539566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(n, &suba)); 543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55fa80d070SPierre Jolivet } 56fa80d070SPierre Jolivet 5766976f2fSJacob Faibussowitsch static PetscErrorCode MatPermute_Normal(Mat A, IS rowp, IS colp, Mat *B) 58d71ae5a4SJacob Faibussowitsch { 5927d8b95cSPierre Jolivet Mat_Normal *a; 6027d8b95cSPierre Jolivet Mat C, Aa; 61fa80d070SPierre Jolivet IS row; 62*b9c875b8SPierre Jolivet PetscScalar shift, scale; 63fa80d070SPierre Jolivet 64fa80d070SPierre Jolivet PetscFunctionBegin; 6508401ef6SPierre Jolivet PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same"); 66*b9c875b8SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 6727d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 6827d8b95cSPierre Jolivet Aa = a->A; 699566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row)); 709566063dSJacob Faibussowitsch PetscCall(ISSetIdentity(row)); 719566063dSJacob Faibussowitsch PetscCall(MatPermute(Aa, row, colp, &C)); 729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 739566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(C, B)); 749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 75*b9c875b8SPierre Jolivet PetscCall(MatShift(*B, shift)); 76*b9c875b8SPierre Jolivet PetscCall(MatScale(*B, scale)); 773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78fa80d070SPierre Jolivet } 79fa80d070SPierre Jolivet 8066976f2fSJacob Faibussowitsch static PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B) 81d71ae5a4SJacob Faibussowitsch { 8227d8b95cSPierre Jolivet Mat_Normal *a; 83fa80d070SPierre Jolivet Mat C; 84fa80d070SPierre Jolivet 85fa80d070SPierre Jolivet PetscFunctionBegin; 8627d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 879566063dSJacob Faibussowitsch PetscCall(MatDuplicate(a->A, op, &C)); 889566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(C, B)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 9027d8b95cSPierre Jolivet if (op == MAT_COPY_VALUES) PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN)); 913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92fa80d070SPierre Jolivet } 93fa80d070SPierre Jolivet 9466976f2fSJacob Faibussowitsch static PetscErrorCode MatCopy_Normal(Mat A, Mat B, MatStructure str) 95d71ae5a4SJacob Faibussowitsch { 9627d8b95cSPierre Jolivet Mat_Normal *a, *b; 97fa80d070SPierre Jolivet 98fa80d070SPierre Jolivet PetscFunctionBegin; 9927d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 10027d8b95cSPierre Jolivet PetscCall(MatShellGetContext(B, &b)); 1019566063dSJacob Faibussowitsch PetscCall(MatCopy(a->A, b->A, str)); 1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103fa80d070SPierre Jolivet } 104fa80d070SPierre Jolivet 10566976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_Normal(Mat N, Vec x, Vec y) 106d71ae5a4SJacob Faibussowitsch { 10727d8b95cSPierre Jolivet Mat_Normal *Na; 108c8a8475eSBarry Smith 109c8a8475eSBarry Smith PetscFunctionBegin; 11027d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 11127d8b95cSPierre Jolivet PetscCall(MatMult(Na->A, x, Na->w)); 1129566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A, Na->w, y)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114db090513SMatthew Knepley } 115db090513SMatthew Knepley 11666976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_Normal(Mat N) 117d71ae5a4SJacob Faibussowitsch { 11827d8b95cSPierre Jolivet Mat_Normal *Na; 119c8a8475eSBarry Smith 120c8a8475eSBarry Smith PetscFunctionBegin; 12127d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 123a48507d3SJose E. Roman PetscCall(MatDestroy(&Na->D)); 1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->w)); 12527d8b95cSPierre Jolivet PetscCall(PetscFree(Na)); 1269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL)); 1279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_seqaij_C", NULL)); 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_mpiaij_C", NULL)); 1292f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE) 1302f0a5c5aSJose E. Roman PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_hypre_C", NULL)); 1312f0a5c5aSJose E. Roman #endif 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_seqdense_C", NULL)); 1339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_mpidense_C", NULL)); 13427d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatShellSetContext_C", NULL)); 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136c8a8475eSBarry Smith } 137c8a8475eSBarry Smith 13817a6fd94SBarry Smith /* 13917a6fd94SBarry Smith Slow, nonscalable version 14017a6fd94SBarry Smith */ 14166976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Normal(Mat N, Vec v) 142d71ae5a4SJacob Faibussowitsch { 14327d8b95cSPierre Jolivet Mat_Normal *Na; 14427d8b95cSPierre Jolivet Mat A; 145b24ad042SBarry Smith PetscInt i, j, rstart, rend, nnz; 146b24ad042SBarry Smith const PetscInt *cols; 14717a6fd94SBarry Smith PetscScalar *diag, *work, *values; 148b3cc6726SBarry Smith const PetscScalar *mvalues; 14917a6fd94SBarry Smith 15017a6fd94SBarry Smith PetscFunctionBegin; 15127d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 15227d8b95cSPierre Jolivet A = Na->A; 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work)); 1549566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(work, A->cmap->N)); 1559566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 15617a6fd94SBarry Smith for (i = rstart; i < rend; i++) { 1579566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues)); 158ad540459SPierre Jolivet for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * mvalues[j]; 1599566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues)); 16017a6fd94SBarry Smith } 1611c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N))); 162d0f46423SBarry Smith rstart = N->cmap->rstart; 163d0f46423SBarry Smith rend = N->cmap->rend; 1649566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &values)); 1659566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart)); 1669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &values)); 1679566063dSJacob Faibussowitsch PetscCall(PetscFree2(diag, work)); 1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16917a6fd94SBarry Smith } 170c8a8475eSBarry Smith 17166976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Normal(Mat N, Mat *D) 172d71ae5a4SJacob Faibussowitsch { 17327d8b95cSPierre Jolivet Mat_Normal *Na; 17427d8b95cSPierre Jolivet Mat M, A; 175a48507d3SJose E. Roman 176a48507d3SJose E. Roman PetscFunctionBegin; 17727d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)N->data)->zrows && !((Mat_Shell *)N->data)->zcols, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME 17827d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)N->data)->axpy, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatAXPY() has been called on the input Mat"); // TODO FIXME 17927d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)N->data)->left && !((Mat_Shell *)N->data)->right, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME 18027d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)N->data)->dshift, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME 18127d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 18227d8b95cSPierre Jolivet A = Na->A; 183a48507d3SJose E. Roman PetscCall(MatGetDiagonalBlock(A, &M)); 184a48507d3SJose E. Roman PetscCall(MatCreateNormal(M, &Na->D)); 185a48507d3SJose E. Roman *D = Na->D; 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 187a48507d3SJose E. Roman } 188a48507d3SJose E. Roman 18966976f2fSJacob Faibussowitsch static PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M) 190d71ae5a4SJacob Faibussowitsch { 19127d8b95cSPierre Jolivet Mat_Normal *Aa; 192fa80d070SPierre Jolivet 193fa80d070SPierre Jolivet PetscFunctionBegin; 19427d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &Aa)); 195fa80d070SPierre Jolivet *M = Aa->A; 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197fa80d070SPierre Jolivet } 198fa80d070SPierre Jolivet 199fa80d070SPierre Jolivet /*@ 20011a5261eSBarry Smith MatNormalGetMat - Gets the `Mat` object stored inside a `MATNORMAL` 201fa80d070SPierre Jolivet 2022ef1f0ffSBarry Smith Logically Collective 203fa80d070SPierre Jolivet 204fa80d070SPierre Jolivet Input Parameter: 20511a5261eSBarry Smith . A - the `MATNORMAL` matrix 206fa80d070SPierre Jolivet 207fa80d070SPierre Jolivet Output Parameter: 2082ef1f0ffSBarry Smith . M - the matrix object stored inside `A` 209fa80d070SPierre Jolivet 210fa80d070SPierre Jolivet Level: intermediate 211fa80d070SPierre Jolivet 2121cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatCreateNormal()` 213fa80d070SPierre Jolivet @*/ 214d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat(Mat A, Mat *M) 215d71ae5a4SJacob Faibussowitsch { 216fa80d070SPierre Jolivet PetscFunctionBegin; 217fa80d070SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 218fa80d070SPierre Jolivet PetscValidType(A, 1); 2194f572ea9SToby Isaac PetscAssertPointer(M, 2); 220cac4c232SBarry Smith PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M)); 2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 222fa80d070SPierre Jolivet } 223fa80d070SPierre Jolivet 22466976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 225d71ae5a4SJacob Faibussowitsch { 22627d8b95cSPierre Jolivet Mat_Normal *Aa; 227fa80d070SPierre Jolivet Mat B; 228fa80d070SPierre Jolivet PetscInt m, n, M, N; 229fa80d070SPierre Jolivet 230fa80d070SPierre Jolivet PetscFunctionBegin; 23127d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &Aa)); 2329566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 2339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 234fa80d070SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 235fa80d070SPierre Jolivet B = *newmat; 2369566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B)); 237fa80d070SPierre Jolivet } else { 2389566063dSJacob Faibussowitsch PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B)); 2399566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B, MATPRODUCT_AtB)); 2409566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 2419566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 2429566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE)); 243fa80d070SPierre Jolivet } 2449566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 245fa80d070SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 2469566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 247fa80d070SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 2489566063dSJacob Faibussowitsch PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat)); 2493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 250fa80d070SPierre Jolivet } 251fa80d070SPierre Jolivet 2522f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE) 25366976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 2542f0a5c5aSJose E. Roman { 2552f0a5c5aSJose E. Roman PetscFunctionBegin; 2562f0a5c5aSJose E. Roman if (reuse == MAT_INITIAL_MATRIX) { 2572f0a5c5aSJose E. Roman PetscCall(MatConvert(A, MATAIJ, reuse, B)); 2582f0a5c5aSJose E. Roman PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B)); 2592f0a5c5aSJose E. Roman } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */ 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2612f0a5c5aSJose E. Roman } 2622f0a5c5aSJose E. Roman #endif 2632f0a5c5aSJose E. Roman 264fa80d070SPierre Jolivet typedef struct { 265fa80d070SPierre Jolivet Mat work[2]; 266fa80d070SPierre Jolivet } Normal_Dense; 267fa80d070SPierre Jolivet 26866976f2fSJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Normal_Dense(Mat C) 269d71ae5a4SJacob Faibussowitsch { 270fa80d070SPierre Jolivet Mat A, B; 271fa80d070SPierre Jolivet Normal_Dense *contents; 272fa80d070SPierre Jolivet Mat_Normal *a; 273*b9c875b8SPierre Jolivet Vec right; 274*b9c875b8SPierre Jolivet PetscScalar *array, scale; 275fa80d070SPierre Jolivet 276fa80d070SPierre Jolivet PetscFunctionBegin; 2778f83a4d9SJacob Faibussowitsch MatCheckProduct(C, 1); 278fa80d070SPierre Jolivet A = C->product->A; 279fa80d070SPierre Jolivet B = C->product->B; 28027d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 281fa80d070SPierre Jolivet contents = (Normal_Dense *)C->product->data; 28228b400f6SJacob Faibussowitsch PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty"); 283*b9c875b8SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, &right, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 284*b9c875b8SPierre Jolivet if (right) { 2859566063dSJacob Faibussowitsch PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN)); 286*b9c875b8SPierre Jolivet PetscCall(MatDiagonalScale(C, right, NULL)); 287fa80d070SPierre Jolivet } 2889566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(contents->work[0])); 2899566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 2909566063dSJacob Faibussowitsch PetscCall(MatDensePlaceArray(contents->work[1], array)); 2919566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(contents->work[1])); 2929566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 2939566063dSJacob Faibussowitsch PetscCall(MatDenseResetArray(contents->work[1])); 2949566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 2959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 2969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 297*b9c875b8SPierre Jolivet PetscCall(MatScale(C, scale)); 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 299fa80d070SPierre Jolivet } 300fa80d070SPierre Jolivet 30166976f2fSJacob Faibussowitsch static PetscErrorCode MatNormal_DenseDestroy(void *ctx) 302d71ae5a4SJacob Faibussowitsch { 303fa80d070SPierre Jolivet Normal_Dense *contents = (Normal_Dense *)ctx; 304fa80d070SPierre Jolivet 305fa80d070SPierre Jolivet PetscFunctionBegin; 3069566063dSJacob Faibussowitsch PetscCall(MatDestroy(contents->work)); 3079566063dSJacob Faibussowitsch PetscCall(MatDestroy(contents->work + 1)); 3089566063dSJacob Faibussowitsch PetscCall(PetscFree(contents)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310fa80d070SPierre Jolivet } 311fa80d070SPierre Jolivet 31266976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C) 313d71ae5a4SJacob Faibussowitsch { 314fa80d070SPierre Jolivet Mat A, B; 315fa80d070SPierre Jolivet Normal_Dense *contents = NULL; 316fa80d070SPierre Jolivet Mat_Normal *a; 317*b9c875b8SPierre Jolivet Vec right; 318*b9c875b8SPierre Jolivet PetscScalar *array, scale; 319fa80d070SPierre Jolivet PetscInt n, N, m, M; 320fa80d070SPierre Jolivet 321fa80d070SPierre Jolivet PetscFunctionBegin; 3228f83a4d9SJacob Faibussowitsch MatCheckProduct(C, 1); 32328b400f6SJacob Faibussowitsch PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty"); 324fa80d070SPierre Jolivet A = C->product->A; 325fa80d070SPierre Jolivet B = C->product->B; 326*b9c875b8SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, &right, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 32727d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 3289566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &m, &n)); 3299566063dSJacob Faibussowitsch PetscCall(MatGetSize(C, &M, &N)); 330fa80d070SPierre Jolivet if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) { 3319566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, NULL, &n)); 3329566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, NULL, &N)); 3339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 3349566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL)); 3359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m, n, M, N)); 336fa80d070SPierre Jolivet } 3379566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)B)->type_name)); 3389566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 3399566063dSJacob Faibussowitsch PetscCall(PetscNew(&contents)); 340fa80d070SPierre Jolivet C->product->data = contents; 341fa80d070SPierre Jolivet C->product->destroy = MatNormal_DenseDestroy; 342*b9c875b8SPierre Jolivet if (right) PetscCall(MatProductCreate(a->A, C, NULL, contents->work)); 343*b9c875b8SPierre Jolivet else PetscCall(MatProductCreate(a->A, B, NULL, contents->work)); 3449566063dSJacob Faibussowitsch PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB)); 3459566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(contents->work[0])); 3469566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(contents->work[0])); 3479566063dSJacob Faibussowitsch PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1)); 3489566063dSJacob Faibussowitsch PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB)); 3499566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(contents->work[1])); 3509566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(contents->work[1])); 3519566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 3529566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array)); 3539566063dSJacob Faibussowitsch PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array)); 3549566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 355fa80d070SPierre Jolivet C->ops->productnumeric = MatProductNumeric_Normal_Dense; 3563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 357fa80d070SPierre Jolivet } 358fa80d070SPierre Jolivet 35966976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C) 360d71ae5a4SJacob Faibussowitsch { 361fa80d070SPierre Jolivet PetscFunctionBegin; 362fa80d070SPierre Jolivet C->ops->productsymbolic = MatProductSymbolic_Normal_Dense; 3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 364fa80d070SPierre Jolivet } 365fa80d070SPierre Jolivet 36666976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C) 367d71ae5a4SJacob Faibussowitsch { 368fa80d070SPierre Jolivet Mat_Product *product = C->product; 369fa80d070SPierre Jolivet 370fa80d070SPierre Jolivet PetscFunctionBegin; 37148a46eb9SPierre Jolivet if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Normal_Dense_AB(C)); 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373fa80d070SPierre Jolivet } 374fa80d070SPierre Jolivet 3752ef1f0ffSBarry Smith /*MC 3762ef1f0ffSBarry Smith MATNORMAL - a matrix that behaves like A'*A for `MatMult()` while only containing A 3772ef1f0ffSBarry Smith 3782ef1f0ffSBarry Smith Level: intermediate 3792ef1f0ffSBarry Smith 38027d8b95cSPierre Jolivet Developer Notes: 38127d8b95cSPierre Jolivet This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code 38227d8b95cSPierre Jolivet 38327d8b95cSPierre Jolivet Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage 38427d8b95cSPierre Jolivet 3851cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 3862ef1f0ffSBarry Smith M*/ 3872ef1f0ffSBarry Smith 388c8a8475eSBarry Smith /*@ 38911a5261eSBarry Smith MatCreateNormal - Creates a new `MATNORMAL` matrix object that behaves like A'*A. 390c8a8475eSBarry Smith 391c3339decSBarry Smith Collective 392c8a8475eSBarry Smith 393c8a8475eSBarry Smith Input Parameter: 394c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 395c8a8475eSBarry Smith 396c8a8475eSBarry Smith Output Parameter: 397c8a8475eSBarry Smith . N - the matrix that represents A'*A 398c8a8475eSBarry Smith 399c30e7cdbSSatish Balay Level: intermediate 400c30e7cdbSSatish Balay 40195452b02SPatrick Sanan Notes: 40295452b02SPatrick Sanan The product A'*A is NOT actually formed! Rather the new matrix 40311a5261eSBarry Smith object performs the matrix-vector product, `MatMult()`, by first multiplying by 404c8a8475eSBarry Smith A and then A' 40511a5261eSBarry Smith 4061cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 407c8a8475eSBarry Smith @*/ 408d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormal(Mat A, Mat *N) 409d71ae5a4SJacob Faibussowitsch { 410c8a8475eSBarry Smith Mat_Normal *Na; 4119ca0e1b6SStefano Zampini VecType vtype; 412c8a8475eSBarry Smith 413c8a8475eSBarry Smith PetscFunctionBegin; 4149566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N)); 4159566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap)); 4169566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap)); 41727d8b95cSPierre Jolivet PetscCall(MatSetType(*N, MATSHELL)); 4184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na)); 41927d8b95cSPierre Jolivet PetscCall(MatShellSetContext(*N, Na)); 4209566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 421c3122656SLisandro Dalcin Na->A = A; 4229566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &Na->w)); 4232205254eSKarl Rupp 42427d8b95cSPierre Jolivet PetscCall(MatSetBlockSizes(*N, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs))); 42527d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (void (*)(void))MatDestroy_Normal)); 42627d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_MULT, (void (*)(void))MatMult_Normal)); 42727d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Normal)); 42827d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_Normal)); 42927d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Normal)); 43027d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_COPY, (void (*)(void))MatCopy_Normal)); 431a48507d3SJose E. Roman (*N)->ops->getdiagonalblock = MatGetDiagonalBlock_Normal; 432fa80d070SPierre Jolivet (*N)->ops->increaseoverlap = MatIncreaseOverlap_Normal; 433fa80d070SPierre Jolivet (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal; 434fa80d070SPierre Jolivet (*N)->ops->permute = MatPermute_Normal; 4359ca0e1b6SStefano Zampini 43627d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalGetMat_Normal)); 43727d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ)); 43827d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ)); 4392f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE) 44027d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_hypre_C", MatConvert_Normal_HYPRE)); 4412f0a5c5aSJose E. Roman #endif 44227d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense)); 44327d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense)); 44427d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable)); 44527d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 44627d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 4479566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE)); 4489566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A, &vtype)); 4499566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N, vtype)); 4509ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 4519566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N, A->boundtocpu)); 4529ca0e1b6SStefano Zampini #endif 45327d8b95cSPierre Jolivet PetscCall(MatSetUp(*N)); 45427d8b95cSPierre Jolivet PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL)); 4553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 456c8a8475eSBarry Smith } 457