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; 31b9c875b8SPierre 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"); 36b9c875b8SPierre 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)); 48b9c875b8SPierre Jolivet PetscCall(MatShift((*submat)[M], shift)); 49b9c875b8SPierre 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; 62b9c875b8SPierre 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"); 66b9c875b8SPierre 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)); 75b9c875b8SPierre Jolivet PetscCall(MatShift(*B, shift)); 76b9c875b8SPierre 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; 147*e91c04dfSPierre Jolivet PetscScalar *work, *values; 148b3cc6726SBarry Smith const PetscScalar *mvalues; 14917a6fd94SBarry Smith 15017a6fd94SBarry Smith PetscFunctionBegin; 15127d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 15227d8b95cSPierre Jolivet A = Na->A; 153*e91c04dfSPierre Jolivet PetscCall(PetscMalloc1(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 } 161*e91c04dfSPierre Jolivet PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, work, 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)); 165*e91c04dfSPierre Jolivet PetscCall(PetscArraycpy(values, work + rstart, rend - rstart)); 1669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &values)); 167*e91c04dfSPierre Jolivet PetscCall(PetscFree(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 PetscCall(MatShellGetContext(N, &Na)); 17827d8b95cSPierre Jolivet A = Na->A; 179a48507d3SJose E. Roman PetscCall(MatGetDiagonalBlock(A, &M)); 180a48507d3SJose E. Roman PetscCall(MatCreateNormal(M, &Na->D)); 181a48507d3SJose E. Roman *D = Na->D; 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183a48507d3SJose E. Roman } 184a48507d3SJose E. Roman 18566976f2fSJacob Faibussowitsch static PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M) 186d71ae5a4SJacob Faibussowitsch { 18727d8b95cSPierre Jolivet Mat_Normal *Aa; 188fa80d070SPierre Jolivet 189fa80d070SPierre Jolivet PetscFunctionBegin; 19027d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &Aa)); 191fa80d070SPierre Jolivet *M = Aa->A; 1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 193fa80d070SPierre Jolivet } 194fa80d070SPierre Jolivet 195fa80d070SPierre Jolivet /*@ 19611a5261eSBarry Smith MatNormalGetMat - Gets the `Mat` object stored inside a `MATNORMAL` 197fa80d070SPierre Jolivet 1982ef1f0ffSBarry Smith Logically Collective 199fa80d070SPierre Jolivet 200fa80d070SPierre Jolivet Input Parameter: 20111a5261eSBarry Smith . A - the `MATNORMAL` matrix 202fa80d070SPierre Jolivet 203fa80d070SPierre Jolivet Output Parameter: 2042ef1f0ffSBarry Smith . M - the matrix object stored inside `A` 205fa80d070SPierre Jolivet 206fa80d070SPierre Jolivet Level: intermediate 207fa80d070SPierre Jolivet 2081cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatCreateNormal()` 209fa80d070SPierre Jolivet @*/ 210d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat(Mat A, Mat *M) 211d71ae5a4SJacob Faibussowitsch { 212fa80d070SPierre Jolivet PetscFunctionBegin; 213fa80d070SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 214fa80d070SPierre Jolivet PetscValidType(A, 1); 2154f572ea9SToby Isaac PetscAssertPointer(M, 2); 216cac4c232SBarry Smith PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M)); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218fa80d070SPierre Jolivet } 219fa80d070SPierre Jolivet 22066976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 221d71ae5a4SJacob Faibussowitsch { 22227d8b95cSPierre Jolivet Mat_Normal *Aa; 223fa80d070SPierre Jolivet Mat B; 224fa80d070SPierre Jolivet PetscInt m, n, M, N; 225fa80d070SPierre Jolivet 226fa80d070SPierre Jolivet PetscFunctionBegin; 22727d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &Aa)); 2289566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 2299566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 230fa80d070SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 231fa80d070SPierre Jolivet B = *newmat; 2329566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B)); 233fa80d070SPierre Jolivet } else { 2349566063dSJacob Faibussowitsch PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B)); 2359566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B, MATPRODUCT_AtB)); 2369566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 2379566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE)); 239fa80d070SPierre Jolivet } 2409566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 241fa80d070SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 2429566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 243fa80d070SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 2449566063dSJacob Faibussowitsch PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat)); 2453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 246fa80d070SPierre Jolivet } 247fa80d070SPierre Jolivet 2482f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE) 24966976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 2502f0a5c5aSJose E. Roman { 2512f0a5c5aSJose E. Roman PetscFunctionBegin; 2522f0a5c5aSJose E. Roman if (reuse == MAT_INITIAL_MATRIX) { 2532f0a5c5aSJose E. Roman PetscCall(MatConvert(A, MATAIJ, reuse, B)); 2542f0a5c5aSJose E. Roman PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B)); 2552f0a5c5aSJose E. Roman } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */ 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2572f0a5c5aSJose E. Roman } 2582f0a5c5aSJose E. Roman #endif 2592f0a5c5aSJose E. Roman 260fa80d070SPierre Jolivet typedef struct { 261fa80d070SPierre Jolivet Mat work[2]; 262fa80d070SPierre Jolivet } Normal_Dense; 263fa80d070SPierre Jolivet 26466976f2fSJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Normal_Dense(Mat C) 265d71ae5a4SJacob Faibussowitsch { 266fa80d070SPierre Jolivet Mat A, B; 267fa80d070SPierre Jolivet Normal_Dense *contents; 268fa80d070SPierre Jolivet Mat_Normal *a; 269b9c875b8SPierre Jolivet Vec right; 270b9c875b8SPierre Jolivet PetscScalar *array, scale; 271fa80d070SPierre Jolivet 272fa80d070SPierre Jolivet PetscFunctionBegin; 2738f83a4d9SJacob Faibussowitsch MatCheckProduct(C, 1); 274fa80d070SPierre Jolivet A = C->product->A; 275fa80d070SPierre Jolivet B = C->product->B; 27627d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 277fa80d070SPierre Jolivet contents = (Normal_Dense *)C->product->data; 27828b400f6SJacob Faibussowitsch PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty"); 27905694d27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 280b9c875b8SPierre Jolivet if (right) { 2819566063dSJacob Faibussowitsch PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN)); 282b9c875b8SPierre Jolivet PetscCall(MatDiagonalScale(C, right, NULL)); 283fa80d070SPierre Jolivet } 2849566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(contents->work[0])); 2859566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 2869566063dSJacob Faibussowitsch PetscCall(MatDensePlaceArray(contents->work[1], array)); 2879566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(contents->work[1])); 2889566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 2899566063dSJacob Faibussowitsch PetscCall(MatDenseResetArray(contents->work[1])); 2909566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 2919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 2929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 293b9c875b8SPierre Jolivet PetscCall(MatScale(C, scale)); 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 295fa80d070SPierre Jolivet } 296fa80d070SPierre Jolivet 29766976f2fSJacob Faibussowitsch static PetscErrorCode MatNormal_DenseDestroy(void *ctx) 298d71ae5a4SJacob Faibussowitsch { 299fa80d070SPierre Jolivet Normal_Dense *contents = (Normal_Dense *)ctx; 300fa80d070SPierre Jolivet 301fa80d070SPierre Jolivet PetscFunctionBegin; 3029566063dSJacob Faibussowitsch PetscCall(MatDestroy(contents->work)); 3039566063dSJacob Faibussowitsch PetscCall(MatDestroy(contents->work + 1)); 3049566063dSJacob Faibussowitsch PetscCall(PetscFree(contents)); 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 306fa80d070SPierre Jolivet } 307fa80d070SPierre Jolivet 30866976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C) 309d71ae5a4SJacob Faibussowitsch { 310fa80d070SPierre Jolivet Mat A, B; 311fa80d070SPierre Jolivet Normal_Dense *contents = NULL; 312fa80d070SPierre Jolivet Mat_Normal *a; 313b9c875b8SPierre Jolivet Vec right; 314b9c875b8SPierre Jolivet PetscScalar *array, scale; 315fa80d070SPierre Jolivet PetscInt n, N, m, M; 316fa80d070SPierre Jolivet 317fa80d070SPierre Jolivet PetscFunctionBegin; 3188f83a4d9SJacob Faibussowitsch MatCheckProduct(C, 1); 31928b400f6SJacob Faibussowitsch PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty"); 320fa80d070SPierre Jolivet A = C->product->A; 321fa80d070SPierre Jolivet B = C->product->B; 32205694d27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 32327d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 3249566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &m, &n)); 3259566063dSJacob Faibussowitsch PetscCall(MatGetSize(C, &M, &N)); 326fa80d070SPierre Jolivet if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) { 3279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, NULL, &n)); 3289566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, NULL, &N)); 3299566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 3309566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL)); 3319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m, n, M, N)); 332fa80d070SPierre Jolivet } 3339566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)B)->type_name)); 3349566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 3359566063dSJacob Faibussowitsch PetscCall(PetscNew(&contents)); 336fa80d070SPierre Jolivet C->product->data = contents; 337fa80d070SPierre Jolivet C->product->destroy = MatNormal_DenseDestroy; 338b9c875b8SPierre Jolivet if (right) PetscCall(MatProductCreate(a->A, C, NULL, contents->work)); 339b9c875b8SPierre Jolivet else PetscCall(MatProductCreate(a->A, B, NULL, contents->work)); 3409566063dSJacob Faibussowitsch PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB)); 3419566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(contents->work[0])); 3429566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(contents->work[0])); 3439566063dSJacob Faibussowitsch PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1)); 3449566063dSJacob Faibussowitsch PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB)); 3459566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(contents->work[1])); 3469566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(contents->work[1])); 3479566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 3489566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array)); 3499566063dSJacob Faibussowitsch PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array)); 3509566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 351fa80d070SPierre Jolivet C->ops->productnumeric = MatProductNumeric_Normal_Dense; 3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 353fa80d070SPierre Jolivet } 354fa80d070SPierre Jolivet 35566976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C) 356d71ae5a4SJacob Faibussowitsch { 357fa80d070SPierre Jolivet Mat_Product *product = C->product; 358fa80d070SPierre Jolivet 359fa80d070SPierre Jolivet PetscFunctionBegin; 360794904c8SPierre Jolivet if (product->type == MATPRODUCT_AB) C->ops->productsymbolic = MatProductSymbolic_Normal_Dense; 3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 362fa80d070SPierre Jolivet } 363fa80d070SPierre Jolivet 3642ef1f0ffSBarry Smith /*MC 3652ef1f0ffSBarry Smith MATNORMAL - a matrix that behaves like A'*A for `MatMult()` while only containing A 3662ef1f0ffSBarry Smith 3672ef1f0ffSBarry Smith Level: intermediate 3682ef1f0ffSBarry Smith 36927d8b95cSPierre Jolivet Developer Notes: 37027d8b95cSPierre Jolivet This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code 37127d8b95cSPierre Jolivet 37227d8b95cSPierre Jolivet Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage 37327d8b95cSPierre Jolivet 3741cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 3752ef1f0ffSBarry Smith M*/ 3762ef1f0ffSBarry Smith 377c8a8475eSBarry Smith /*@ 37811a5261eSBarry Smith MatCreateNormal - Creates a new `MATNORMAL` matrix object that behaves like A'*A. 379c8a8475eSBarry Smith 380c3339decSBarry Smith Collective 381c8a8475eSBarry Smith 382c8a8475eSBarry Smith Input Parameter: 383c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 384c8a8475eSBarry Smith 385c8a8475eSBarry Smith Output Parameter: 386c8a8475eSBarry Smith . N - the matrix that represents A'*A 387c8a8475eSBarry Smith 388c30e7cdbSSatish Balay Level: intermediate 389c30e7cdbSSatish Balay 39095452b02SPatrick Sanan Notes: 39195452b02SPatrick Sanan The product A'*A is NOT actually formed! Rather the new matrix 39211a5261eSBarry Smith object performs the matrix-vector product, `MatMult()`, by first multiplying by 393c8a8475eSBarry Smith A and then A' 39411a5261eSBarry Smith 3951cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 396c8a8475eSBarry Smith @*/ 397d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormal(Mat A, Mat *N) 398d71ae5a4SJacob Faibussowitsch { 399c8a8475eSBarry Smith Mat_Normal *Na; 4009ca0e1b6SStefano Zampini VecType vtype; 401c8a8475eSBarry Smith 402c8a8475eSBarry Smith PetscFunctionBegin; 4039566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N)); 4049566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap)); 4059566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap)); 40627d8b95cSPierre Jolivet PetscCall(MatSetType(*N, MATSHELL)); 4074dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na)); 40827d8b95cSPierre Jolivet PetscCall(MatShellSetContext(*N, Na)); 4099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 410c3122656SLisandro Dalcin Na->A = A; 4119566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &Na->w)); 4122205254eSKarl Rupp 41327d8b95cSPierre Jolivet PetscCall(MatSetBlockSizes(*N, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs))); 41427d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (void (*)(void))MatDestroy_Normal)); 41527d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_MULT, (void (*)(void))MatMult_Normal)); 41627d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Normal)); 41727d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_Normal)); 41827d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Normal)); 419a29b93afSAudic XU PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Normal)); 42027d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_COPY, (void (*)(void))MatCopy_Normal)); 421fa80d070SPierre Jolivet (*N)->ops->increaseoverlap = MatIncreaseOverlap_Normal; 422fa80d070SPierre Jolivet (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal; 423fa80d070SPierre Jolivet (*N)->ops->permute = MatPermute_Normal; 4249ca0e1b6SStefano Zampini 42527d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalGetMat_Normal)); 42627d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ)); 42727d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ)); 4282f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE) 42927d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_hypre_C", MatConvert_Normal_HYPRE)); 4302f0a5c5aSJose E. Roman #endif 43127d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense)); 43227d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense)); 43327d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable)); 43427d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 43527d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 4369566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE)); 4379566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A, &vtype)); 4389566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N, vtype)); 4399ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 4409566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N, A->boundtocpu)); 4419ca0e1b6SStefano Zampini #endif 44227d8b95cSPierre Jolivet PetscCall(MatSetUp(*N)); 44327d8b95cSPierre Jolivet PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL)); 4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 445c8a8475eSBarry Smith } 446