1*27d8b95cSPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/ 2ebda224cSfranklin5 3ebda224cSfranklin5 typedef struct { 4ebda224cSfranklin5 Mat A; 5a48507d3SJose E. Roman Mat D; /* local submatrix for diagonal part */ 6*27d8b95cSPierre Jolivet Vec w; 713b343c1SBarry Smith } Mat_NormalHermitian; 8ebda224cSfranklin5 966976f2fSJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_NormalHermitian(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 10d71ae5a4SJacob Faibussowitsch { 11*27d8b95cSPierre Jolivet Mat_NormalHermitian *a; 12*27d8b95cSPierre Jolivet Mat B, *suba; 13ab704866SPierre Jolivet IS *row; 14ab704866SPierre Jolivet PetscInt M; 15ab704866SPierre Jolivet 16ab704866SPierre Jolivet PetscFunctionBegin; 17*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)mat->data)->zrows && !((Mat_Shell *)mat->data)->zcols, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatZeroRows()/MatZeroRowsColumns() after the SubMatrices creation 18*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)mat->data)->axpy, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatAXPY() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatAXPY() after the SubMatrices creation 19*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)mat->data)->left && !((Mat_Shell *)mat->data)->right, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalScale() after the SubMatrices creation with a SubVector 20*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)mat->data)->dshift, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalSet() after the SubMatrices creation with a SubVector 21*27d8b95cSPierre Jolivet PetscCheck(irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented"); 22*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(mat, &a)); 23*27d8b95cSPierre Jolivet B = a->A; 2448a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat)); 25ab704866SPierre Jolivet PetscCall(MatGetSize(B, &M, NULL)); 26ab704866SPierre Jolivet PetscCall(PetscMalloc1(n, &row)); 27ab704866SPierre Jolivet PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0])); 28ab704866SPierre Jolivet PetscCall(ISSetIdentity(row[0])); 29ab704866SPierre Jolivet for (M = 1; M < n; ++M) row[M] = row[0]; 30ab704866SPierre Jolivet PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba)); 31ab704866SPierre Jolivet for (M = 0; M < n; ++M) { 32ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(suba[M], *submat + M)); 33*27d8b95cSPierre Jolivet ((Mat_Shell *)(*submat)[M]->data)->vscale = ((Mat_Shell *)mat->data)->vscale; 34*27d8b95cSPierre Jolivet ((Mat_Shell *)(*submat)[M]->data)->vshift = ((Mat_Shell *)mat->data)->vshift; 35ab704866SPierre Jolivet } 36ab704866SPierre Jolivet PetscCall(ISDestroy(&row[0])); 37ab704866SPierre Jolivet PetscCall(PetscFree(row)); 38ab704866SPierre Jolivet PetscCall(MatDestroySubMatrices(n, &suba)); 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40ab704866SPierre Jolivet } 41ab704866SPierre Jolivet 4266976f2fSJacob Faibussowitsch static PetscErrorCode MatPermute_NormalHermitian(Mat A, IS rowp, IS colp, Mat *B) 43d71ae5a4SJacob Faibussowitsch { 44*27d8b95cSPierre Jolivet Mat_NormalHermitian *a; 45*27d8b95cSPierre Jolivet Mat C, Aa; 46ab704866SPierre Jolivet IS row; 47ab704866SPierre Jolivet 48ab704866SPierre Jolivet PetscFunctionBegin; 49*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatZeroRows()/MatZeroRowsColumns() after the permutation with a permuted zrows and zcols 50*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatAXPY() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatAXPY() after the permutation with a permuted axpy 51*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)A->data)->left && !((Mat_Shell *)A->data)->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalScale() after the permutation with a permuted left and right 52*27d8b95cSPierre Jolivet PetscCheck(!((Mat_Shell *)A->data)->dshift, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalSet() after the permutation with a permuted dshift 53ab704866SPierre Jolivet PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same"); 54*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 55*27d8b95cSPierre Jolivet Aa = a->A; 56ab704866SPierre Jolivet PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row)); 57ab704866SPierre Jolivet PetscCall(ISSetIdentity(row)); 58ab704866SPierre Jolivet PetscCall(MatPermute(Aa, row, colp, &C)); 59ab704866SPierre Jolivet PetscCall(ISDestroy(&row)); 60ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(C, B)); 61ab704866SPierre Jolivet PetscCall(MatDestroy(&C)); 62*27d8b95cSPierre Jolivet ((Mat_Shell *)(*B)->data)->vscale = ((Mat_Shell *)A->data)->vscale; 63*27d8b95cSPierre Jolivet ((Mat_Shell *)(*B)->data)->vshift = ((Mat_Shell *)A->data)->vshift; 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65ab704866SPierre Jolivet } 66ab704866SPierre Jolivet 6766976f2fSJacob Faibussowitsch static PetscErrorCode MatDuplicate_NormalHermitian(Mat A, MatDuplicateOption op, Mat *B) 68d71ae5a4SJacob Faibussowitsch { 69*27d8b95cSPierre Jolivet Mat_NormalHermitian *a; 70ab704866SPierre Jolivet Mat C; 71ab704866SPierre Jolivet 72ab704866SPierre Jolivet PetscFunctionBegin; 73*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 74ab704866SPierre Jolivet PetscCall(MatDuplicate(a->A, op, &C)); 75ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(C, B)); 76ab704866SPierre Jolivet PetscCall(MatDestroy(&C)); 77*27d8b95cSPierre Jolivet if (op == MAT_COPY_VALUES) PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN)); 783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79ab704866SPierre Jolivet } 80ab704866SPierre Jolivet 8166976f2fSJacob Faibussowitsch static PetscErrorCode MatCopy_NormalHermitian(Mat A, Mat B, MatStructure str) 82d71ae5a4SJacob Faibussowitsch { 83*27d8b95cSPierre Jolivet Mat_NormalHermitian *a, *b; 84ab704866SPierre Jolivet 85ab704866SPierre Jolivet PetscFunctionBegin; 86*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 87*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(B, &b)); 88ab704866SPierre Jolivet PetscCall(MatCopy(a->A, b->A, str)); 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90ab704866SPierre Jolivet } 91ab704866SPierre Jolivet 9266976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_NormalHermitian(Mat N, Vec x, Vec y) 93d71ae5a4SJacob Faibussowitsch { 94*27d8b95cSPierre Jolivet Mat_NormalHermitian *Na; 95ebda224cSfranklin5 96ebda224cSfranklin5 PetscFunctionBegin; 97*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 98*27d8b95cSPierre Jolivet PetscCall(MatMult(Na->A, x, Na->w)); 999566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y)); 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 101ebda224cSfranklin5 } 102ebda224cSfranklin5 10366976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_NormalHermitian(Mat N) 104d71ae5a4SJacob Faibussowitsch { 105*27d8b95cSPierre Jolivet Mat_NormalHermitian *Na; 106ebda224cSfranklin5 107ebda224cSfranklin5 PetscFunctionBegin; 108*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 110a48507d3SJose E. Roman PetscCall(MatDestroy(&Na->D)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->w)); 112*27d8b95cSPierre Jolivet PetscCall(PetscFree(Na)); 11314e4dea2SJose E. Roman PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalHermitianGetMat_C", NULL)); 114*27d8b95cSPierre Jolivet #if !defined(PETSC_USE_COMPLEX) 115*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL)); 116*27d8b95cSPierre Jolivet #endif 117ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_seqaij_C", NULL)); 118ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_mpiaij_C", NULL)); 1192f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE) 1202f0a5c5aSJose E. Roman PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_hypre_C", NULL)); 1212f0a5c5aSJose E. Roman #endif 122*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatShellSetContext_C", NULL)); 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124ebda224cSfranklin5 } 125ebda224cSfranklin5 126ebda224cSfranklin5 /* 127ebda224cSfranklin5 Slow, nonscalable version 128ebda224cSfranklin5 */ 12966976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N, Vec v) 130d71ae5a4SJacob Faibussowitsch { 131*27d8b95cSPierre Jolivet Mat_NormalHermitian *Na; 132*27d8b95cSPierre Jolivet Mat A; 133ebda224cSfranklin5 PetscInt i, j, rstart, rend, nnz; 134ebda224cSfranklin5 const PetscInt *cols; 135ebda224cSfranklin5 PetscScalar *diag, *work, *values; 136ebda224cSfranklin5 const PetscScalar *mvalues; 137ebda224cSfranklin5 138ebda224cSfranklin5 PetscFunctionBegin; 139*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 140*27d8b95cSPierre Jolivet A = Na->A; 1419566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work)); 1429566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(work, A->cmap->N)); 1439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 144ebda224cSfranklin5 for (i = rstart; i < rend; i++) { 1459566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues)); 146ad540459SPierre Jolivet for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * PetscConj(mvalues[j]); 1479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues)); 148ebda224cSfranklin5 } 1491c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N))); 150ebda224cSfranklin5 rstart = N->cmap->rstart; 151ebda224cSfranklin5 rend = N->cmap->rend; 1529566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &values)); 1539566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart)); 1549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &values)); 1559566063dSJacob Faibussowitsch PetscCall(PetscFree2(diag, work)); 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 157ebda224cSfranklin5 } 158ebda224cSfranklin5 15966976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N, Mat *D) 160d71ae5a4SJacob Faibussowitsch { 161*27d8b95cSPierre Jolivet Mat_NormalHermitian *Na; 162*27d8b95cSPierre Jolivet Mat M, A; 163a48507d3SJose E. Roman 164a48507d3SJose E. Roman PetscFunctionBegin; 165*27d8b95cSPierre 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 166*27d8b95cSPierre 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 167*27d8b95cSPierre 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 168*27d8b95cSPierre 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 169*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(N, &Na)); 170*27d8b95cSPierre Jolivet A = Na->A; 171a48507d3SJose E. Roman PetscCall(MatGetDiagonalBlock(A, &M)); 172a48507d3SJose E. Roman PetscCall(MatCreateNormalHermitian(M, &Na->D)); 173a48507d3SJose E. Roman *D = Na->D; 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 175a48507d3SJose E. Roman } 176a48507d3SJose E. Roman 17714e4dea2SJose E. Roman static PetscErrorCode MatNormalHermitianGetMat_NormalHermitian(Mat A, Mat *M) 178d71ae5a4SJacob Faibussowitsch { 179*27d8b95cSPierre Jolivet Mat_NormalHermitian *Aa; 18065f45395SPierre Jolivet 18165f45395SPierre Jolivet PetscFunctionBegin; 182*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &Aa)); 18365f45395SPierre Jolivet *M = Aa->A; 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18565f45395SPierre Jolivet } 18665f45395SPierre Jolivet 18765f45395SPierre Jolivet /*@ 18811a5261eSBarry Smith MatNormalHermitianGetMat - Gets the `Mat` object stored inside a `MATNORMALHERMITIAN` 18965f45395SPierre Jolivet 1902ef1f0ffSBarry Smith Logically Collective 19165f45395SPierre Jolivet 19265f45395SPierre Jolivet Input Parameter: 19311a5261eSBarry Smith . A - the `MATNORMALHERMITIAN` matrix 19465f45395SPierre Jolivet 19565f45395SPierre Jolivet Output Parameter: 196*27d8b95cSPierre Jolivet . M - the matrix object stored inside `A` 19765f45395SPierre Jolivet 19865f45395SPierre Jolivet Level: intermediate 19965f45395SPierre Jolivet 2001cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 20165f45395SPierre Jolivet @*/ 202d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalHermitianGetMat(Mat A, Mat *M) 203d71ae5a4SJacob Faibussowitsch { 20465f45395SPierre Jolivet PetscFunctionBegin; 20565f45395SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 20665f45395SPierre Jolivet PetscValidType(A, 1); 2074f572ea9SToby Isaac PetscAssertPointer(M, 2); 20814e4dea2SJose E. Roman PetscUseMethod(A, "MatNormalHermitianGetMat_C", (Mat, Mat *), (A, M)); 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21065f45395SPierre Jolivet } 21165f45395SPierre Jolivet 21266976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_NormalHermitian_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 213d71ae5a4SJacob Faibussowitsch { 214*27d8b95cSPierre Jolivet Mat_NormalHermitian *Aa; 215ab704866SPierre Jolivet Mat B, conjugate; 216ab704866SPierre Jolivet PetscInt m, n, M, N; 217ab704866SPierre Jolivet 218ab704866SPierre Jolivet PetscFunctionBegin; 219*27d8b95cSPierre Jolivet PetscCall(MatShellGetContext(A, &Aa)); 220ab704866SPierre Jolivet PetscCall(MatGetSize(A, &M, &N)); 221ab704866SPierre Jolivet PetscCall(MatGetLocalSize(A, &m, &n)); 222ab704866SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 223ab704866SPierre Jolivet B = *newmat; 224ab704866SPierre Jolivet PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B)); 225ab704866SPierre Jolivet } else { 226ab704866SPierre Jolivet PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B)); 227ab704866SPierre Jolivet PetscCall(MatProductSetType(B, MATPRODUCT_AtB)); 228ab704866SPierre Jolivet PetscCall(MatProductSetFromOptions(B)); 229ab704866SPierre Jolivet PetscCall(MatProductSymbolic(B)); 230ab704866SPierre Jolivet PetscCall(MatSetOption(B, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE)); 231ab704866SPierre Jolivet } 232ab704866SPierre Jolivet if (PetscDefined(USE_COMPLEX)) { 233ab704866SPierre Jolivet PetscCall(MatDuplicate(Aa->A, MAT_COPY_VALUES, &conjugate)); 234ab704866SPierre Jolivet PetscCall(MatConjugate(conjugate)); 235ab704866SPierre Jolivet PetscCall(MatProductReplaceMats(conjugate, Aa->A, NULL, B)); 236ab704866SPierre Jolivet } 237ab704866SPierre Jolivet PetscCall(MatProductNumeric(B)); 238ab704866SPierre Jolivet if (PetscDefined(USE_COMPLEX)) PetscCall(MatDestroy(&conjugate)); 239ab704866SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 240ab704866SPierre Jolivet PetscCall(MatHeaderReplace(A, &B)); 241ab704866SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 242ab704866SPierre Jolivet PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat)); 2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 244ab704866SPierre Jolivet } 245ab704866SPierre Jolivet 2462f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE) 24766976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_NormalHermitian_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 2482f0a5c5aSJose E. Roman { 2492f0a5c5aSJose E. Roman PetscFunctionBegin; 2502f0a5c5aSJose E. Roman if (reuse == MAT_INITIAL_MATRIX) { 2512f0a5c5aSJose E. Roman PetscCall(MatConvert(A, MATAIJ, reuse, B)); 2522f0a5c5aSJose E. Roman PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B)); 2532f0a5c5aSJose E. Roman } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */ 2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2552f0a5c5aSJose E. Roman } 2562f0a5c5aSJose E. Roman #endif 2572f0a5c5aSJose E. Roman 2582ef1f0ffSBarry Smith /*MC 2592ef1f0ffSBarry Smith MATNORMALHERMITIAN - a matrix that behaves like (A*)'*A for `MatMult()` while only containing A 2602ef1f0ffSBarry Smith 2612ef1f0ffSBarry Smith Level: intermediate 2622ef1f0ffSBarry Smith 263*27d8b95cSPierre Jolivet Developer Notes: 264*27d8b95cSPierre Jolivet This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code 265*27d8b95cSPierre Jolivet 266*27d8b95cSPierre Jolivet Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage 267*27d8b95cSPierre Jolivet 2681cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormalHermitian()`, `MatMult()`, `MatNormalHermitianGetMat()`, `MATNORMAL`, `MatCreateNormal()` 2692ef1f0ffSBarry Smith M*/ 2702ef1f0ffSBarry Smith 271ebda224cSfranklin5 /*@ 27211a5261eSBarry Smith MatCreateNormalHermitian - Creates a new matrix object `MATNORMALHERMITIAN` that behaves like (A*)'*A. 273ebda224cSfranklin5 274c3339decSBarry Smith Collective 275ebda224cSfranklin5 276ebda224cSfranklin5 Input Parameter: 277ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix 278ebda224cSfranklin5 279ebda224cSfranklin5 Output Parameter: 280ebda224cSfranklin5 . N - the matrix that represents (A*)'*A 281ebda224cSfranklin5 282ebda224cSfranklin5 Level: intermediate 283ebda224cSfranklin5 28411a5261eSBarry Smith Note: 28595452b02SPatrick Sanan The product (A*)'*A is NOT actually formed! Rather the new matrix 28611a5261eSBarry Smith object performs the matrix-vector product, `MatMult()`, by first multiplying by 287ebda224cSfranklin5 A and then (A*)' 28811a5261eSBarry Smith 2891cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatNormalHermitianGetMat()` 290ebda224cSfranklin5 @*/ 291d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormalHermitian(Mat A, Mat *N) 292d71ae5a4SJacob Faibussowitsch { 29313b343c1SBarry Smith Mat_NormalHermitian *Na; 2945cb049f5SJose E. Roman VecType vtype; 295ebda224cSfranklin5 296ebda224cSfranklin5 PetscFunctionBegin; 2979566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N)); 2989566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap)); 2999566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap)); 300*27d8b95cSPierre Jolivet PetscCall(MatSetType(*N, MATSHELL)); 3014dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na)); 302*27d8b95cSPierre Jolivet PetscCall(MatShellSetContext(*N, Na)); 3039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 304ebda224cSfranklin5 Na->A = A; 3059566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &Na->w)); 306ebda224cSfranklin5 307*27d8b95cSPierre Jolivet PetscCall(MatSetBlockSizes(*N, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs))); 308*27d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (void (*)(void))MatDestroy_NormalHermitian)); 309*27d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_MULT, (void (*)(void))MatMult_NormalHermitian)); 310*27d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_MULT_HERMITIAN_TRANSPOSE, (void (*)(void))MatMult_NormalHermitian)); 311*27d8b95cSPierre Jolivet #if !defined(PETSC_USE_COMPLEX) 312*27d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_NormalHermitian)); 313*27d8b95cSPierre Jolivet #endif 314*27d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_NormalHermitian)); 315*27d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NormalHermitian)); 316*27d8b95cSPierre Jolivet PetscCall(MatShellSetOperation(*N, MATOP_COPY, (void (*)(void))MatCopy_NormalHermitian)); 3172e5cc420SJose E. Roman (*N)->ops->getdiagonalblock = MatGetDiagonalBlock_NormalHermitian; 318ab704866SPierre Jolivet (*N)->ops->createsubmatrices = MatCreateSubMatrices_NormalHermitian; 319ab704866SPierre Jolivet (*N)->ops->permute = MatPermute_NormalHermitian; 32012ab1b64SPierre Jolivet 321*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalHermitianGetMat_C", MatNormalHermitianGetMat_NormalHermitian)); 322*27d8b95cSPierre Jolivet #if !defined(PETSC_USE_COMPLEX) 323*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalHermitianGetMat_NormalHermitian)); 3242f0a5c5aSJose E. Roman #endif 325*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normalh_seqaij_C", MatConvert_NormalHermitian_AIJ)); 326*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normalh_mpiaij_C", MatConvert_NormalHermitian_AIJ)); 327*27d8b95cSPierre Jolivet #if defined(PETSC_HAVE_HYPRE) 328*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normalh_hypre_C", MatConvert_NormalHermitian_HYPRE)); 329*27d8b95cSPierre Jolivet #endif 330*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable)); 331*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 332*27d8b95cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 333*27d8b95cSPierre Jolivet PetscCall(MatSetOption(*N, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE)); 3349566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A, &vtype)); 3359566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N, vtype)); 3365cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE) 3379566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N, A->boundtocpu)); 3385cb049f5SJose E. Roman #endif 339*27d8b95cSPierre Jolivet PetscCall(MatSetUp(*N)); 340*27d8b95cSPierre Jolivet PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMALHERMITIAN)); 3413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 342ebda224cSfranklin5 } 343