1ebda224cSfranklin5 2ebda224cSfranklin5 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3ebda224cSfranklin5 4ebda224cSfranklin5 typedef struct { 5ebda224cSfranklin5 Mat A; 6a48507d3SJose E. Roman Mat D; /* local submatrix for diagonal part */ 7ebda224cSfranklin5 Vec w, left, right, leftwork, rightwork; 8ebda224cSfranklin5 PetscScalar scale; 9ebda224cSfranklin5 } Mat_Normal; 10ebda224cSfranklin5 11d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_NormalHermitian(Mat inA, PetscScalar scale) 12d71ae5a4SJacob Faibussowitsch { 13ebda224cSfranklin5 Mat_Normal *a = (Mat_Normal *)inA->data; 14ebda224cSfranklin5 15ebda224cSfranklin5 PetscFunctionBegin; 16ebda224cSfranklin5 a->scale *= scale; 17ebda224cSfranklin5 PetscFunctionReturn(0); 18ebda224cSfranklin5 } 19ebda224cSfranklin5 20d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_NormalHermitian(Mat inA, Vec left, Vec right) 21d71ae5a4SJacob Faibussowitsch { 22ebda224cSfranklin5 Mat_Normal *a = (Mat_Normal *)inA->data; 23ebda224cSfranklin5 24ebda224cSfranklin5 PetscFunctionBegin; 25ebda224cSfranklin5 if (left) { 26ebda224cSfranklin5 if (!a->left) { 279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left, &a->left)); 289566063dSJacob Faibussowitsch PetscCall(VecCopy(left, a->left)); 29ebda224cSfranklin5 } else { 309566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->left, left, a->left)); 31ebda224cSfranklin5 } 32ebda224cSfranklin5 } 33ebda224cSfranklin5 if (right) { 34ebda224cSfranklin5 if (!a->right) { 359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right, &a->right)); 369566063dSJacob Faibussowitsch PetscCall(VecCopy(right, a->right)); 37ebda224cSfranklin5 } else { 389566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->right, right, a->right)); 39ebda224cSfranklin5 } 40ebda224cSfranklin5 } 41ebda224cSfranklin5 PetscFunctionReturn(0); 42ebda224cSfranklin5 } 43ebda224cSfranklin5 44d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_NormalHermitian(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 45d71ae5a4SJacob Faibussowitsch { 46ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal *)mat->data; 47ab704866SPierre Jolivet Mat B = a->A, *suba; 48ab704866SPierre Jolivet IS *row; 49ab704866SPierre Jolivet PetscInt M; 50ab704866SPierre Jolivet 51ab704866SPierre Jolivet PetscFunctionBegin; 52ab704866SPierre Jolivet PetscCheck(!a->left && !a->right && irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented"); 5348a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat)); 54ab704866SPierre Jolivet PetscCall(MatGetSize(B, &M, NULL)); 55ab704866SPierre Jolivet PetscCall(PetscMalloc1(n, &row)); 56ab704866SPierre Jolivet PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0])); 57ab704866SPierre Jolivet PetscCall(ISSetIdentity(row[0])); 58ab704866SPierre Jolivet for (M = 1; M < n; ++M) row[M] = row[0]; 59ab704866SPierre Jolivet PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba)); 60ab704866SPierre Jolivet for (M = 0; M < n; ++M) { 61ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(suba[M], *submat + M)); 62ab704866SPierre Jolivet ((Mat_Normal *)(*submat)[M]->data)->scale = a->scale; 63ab704866SPierre Jolivet } 64ab704866SPierre Jolivet PetscCall(ISDestroy(&row[0])); 65ab704866SPierre Jolivet PetscCall(PetscFree(row)); 66ab704866SPierre Jolivet PetscCall(MatDestroySubMatrices(n, &suba)); 67ab704866SPierre Jolivet PetscFunctionReturn(0); 68ab704866SPierre Jolivet } 69ab704866SPierre Jolivet 70d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPermute_NormalHermitian(Mat A, IS rowp, IS colp, Mat *B) 71d71ae5a4SJacob Faibussowitsch { 72ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal *)A->data; 73ab704866SPierre Jolivet Mat C, Aa = a->A; 74ab704866SPierre Jolivet IS row; 75ab704866SPierre Jolivet 76ab704866SPierre Jolivet PetscFunctionBegin; 77ab704866SPierre Jolivet PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same"); 78ab704866SPierre Jolivet PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row)); 79ab704866SPierre Jolivet PetscCall(ISSetIdentity(row)); 80ab704866SPierre Jolivet PetscCall(MatPermute(Aa, row, colp, &C)); 81ab704866SPierre Jolivet PetscCall(ISDestroy(&row)); 82ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(C, B)); 83ab704866SPierre Jolivet PetscCall(MatDestroy(&C)); 84ab704866SPierre Jolivet PetscFunctionReturn(0); 85ab704866SPierre Jolivet } 86ab704866SPierre Jolivet 87d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_NormalHermitian(Mat A, MatDuplicateOption op, Mat *B) 88d71ae5a4SJacob Faibussowitsch { 89ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal *)A->data; 90ab704866SPierre Jolivet Mat C; 91ab704866SPierre Jolivet 92ab704866SPierre Jolivet PetscFunctionBegin; 93ab704866SPierre Jolivet PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented"); 94ab704866SPierre Jolivet PetscCall(MatDuplicate(a->A, op, &C)); 95ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(C, B)); 96ab704866SPierre Jolivet PetscCall(MatDestroy(&C)); 97ab704866SPierre Jolivet if (op == MAT_COPY_VALUES) ((Mat_Normal *)(*B)->data)->scale = a->scale; 98ab704866SPierre Jolivet PetscFunctionReturn(0); 99ab704866SPierre Jolivet } 100ab704866SPierre Jolivet 101d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_NormalHermitian(Mat A, Mat B, MatStructure str) 102d71ae5a4SJacob Faibussowitsch { 103ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal *)A->data, *b = (Mat_Normal *)B->data; 104ab704866SPierre Jolivet 105ab704866SPierre Jolivet PetscFunctionBegin; 106ab704866SPierre Jolivet PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented"); 107ab704866SPierre Jolivet PetscCall(MatCopy(a->A, b->A, str)); 108ab704866SPierre Jolivet b->scale = a->scale; 109ab704866SPierre Jolivet PetscCall(VecDestroy(&b->left)); 110ab704866SPierre Jolivet PetscCall(VecDestroy(&b->right)); 111ab704866SPierre Jolivet PetscCall(VecDestroy(&b->leftwork)); 112ab704866SPierre Jolivet PetscCall(VecDestroy(&b->rightwork)); 113ab704866SPierre Jolivet PetscFunctionReturn(0); 114ab704866SPierre Jolivet } 115ab704866SPierre Jolivet 116d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_NormalHermitian(Mat N, Vec x, Vec y) 117d71ae5a4SJacob Faibussowitsch { 118ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal *)N->data; 119ebda224cSfranklin5 Vec in; 120ebda224cSfranklin5 121ebda224cSfranklin5 PetscFunctionBegin; 122ebda224cSfranklin5 in = x; 123ebda224cSfranklin5 if (Na->right) { 12448a46eb9SPierre Jolivet if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork)); 1259566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in)); 126ebda224cSfranklin5 in = Na->rightwork; 127ebda224cSfranklin5 } 1289566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, in, Na->w)); 1299566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y)); 1301baa6e33SBarry Smith if (Na->left) PetscCall(VecPointwiseMult(y, Na->left, y)); 1319566063dSJacob Faibussowitsch PetscCall(VecScale(y, Na->scale)); 132ebda224cSfranklin5 PetscFunctionReturn(0); 133ebda224cSfranklin5 } 134ebda224cSfranklin5 135d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3) 136d71ae5a4SJacob Faibussowitsch { 137ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal *)N->data; 138ebda224cSfranklin5 Vec in; 139ebda224cSfranklin5 140ebda224cSfranklin5 PetscFunctionBegin; 141ebda224cSfranklin5 in = v1; 142ebda224cSfranklin5 if (Na->right) { 14348a46eb9SPierre Jolivet if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork)); 1449566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in)); 145ebda224cSfranklin5 in = Na->rightwork; 146ebda224cSfranklin5 } 1479566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, in, Na->w)); 1489566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w, Na->scale)); 149ebda224cSfranklin5 if (Na->left) { 1509566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3)); 1519566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3, Na->left, v3)); 1529566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3, 1.0, v2)); 153ebda224cSfranklin5 } else { 1549566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3)); 155ebda224cSfranklin5 } 156ebda224cSfranklin5 PetscFunctionReturn(0); 157ebda224cSfranklin5 } 158ebda224cSfranklin5 159d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTranspose_Normal(Mat N, Vec x, Vec y) 160d71ae5a4SJacob Faibussowitsch { 161ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal *)N->data; 162ebda224cSfranklin5 Vec in; 163ebda224cSfranklin5 164ebda224cSfranklin5 PetscFunctionBegin; 165ebda224cSfranklin5 in = x; 166ebda224cSfranklin5 if (Na->left) { 16748a46eb9SPierre Jolivet if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork)); 1689566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in)); 169ebda224cSfranklin5 in = Na->leftwork; 170ebda224cSfranklin5 } 1719566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, in, Na->w)); 1729566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y)); 1731baa6e33SBarry Smith if (Na->right) PetscCall(VecPointwiseMult(y, Na->right, y)); 1749566063dSJacob Faibussowitsch PetscCall(VecScale(y, Na->scale)); 175ebda224cSfranklin5 PetscFunctionReturn(0); 176ebda224cSfranklin5 } 177ebda224cSfranklin5 178d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3) 179d71ae5a4SJacob Faibussowitsch { 180ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal *)N->data; 181ebda224cSfranklin5 Vec in; 182ebda224cSfranklin5 183ebda224cSfranklin5 PetscFunctionBegin; 184ebda224cSfranklin5 in = v1; 185ebda224cSfranklin5 if (Na->left) { 18648a46eb9SPierre Jolivet if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork)); 1879566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in)); 188ebda224cSfranklin5 in = Na->leftwork; 189ebda224cSfranklin5 } 1909566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, in, Na->w)); 1919566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w, Na->scale)); 192ebda224cSfranklin5 if (Na->right) { 1939566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3)); 1949566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3, Na->right, v3)); 1959566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3, 1.0, v2)); 196ebda224cSfranklin5 } else { 1979566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3)); 198ebda224cSfranklin5 } 199ebda224cSfranklin5 PetscFunctionReturn(0); 200ebda224cSfranklin5 } 201ebda224cSfranklin5 202d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_NormalHermitian(Mat N) 203d71ae5a4SJacob Faibussowitsch { 204ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal *)N->data; 205ebda224cSfranklin5 206ebda224cSfranklin5 PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 208a48507d3SJose E. Roman PetscCall(MatDestroy(&Na->D)); 2099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->w)); 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->left)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->right)); 2129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->leftwork)); 2139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rightwork)); 2149566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 2159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMatHermitian_C", NULL)); 216ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_seqaij_C", NULL)); 217ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_mpiaij_C", NULL)); 218ebda224cSfranklin5 PetscFunctionReturn(0); 219ebda224cSfranklin5 } 220ebda224cSfranklin5 221ebda224cSfranklin5 /* 222ebda224cSfranklin5 Slow, nonscalable version 223ebda224cSfranklin5 */ 224d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N, Vec v) 225d71ae5a4SJacob Faibussowitsch { 226ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal *)N->data; 227ebda224cSfranklin5 Mat A = Na->A; 228ebda224cSfranklin5 PetscInt i, j, rstart, rend, nnz; 229ebda224cSfranklin5 const PetscInt *cols; 230ebda224cSfranklin5 PetscScalar *diag, *work, *values; 231ebda224cSfranklin5 const PetscScalar *mvalues; 232ebda224cSfranklin5 233ebda224cSfranklin5 PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work)); 2359566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(work, A->cmap->N)); 2369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 237ebda224cSfranklin5 for (i = rstart; i < rend; i++) { 2389566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues)); 239ad540459SPierre Jolivet for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * PetscConj(mvalues[j]); 2409566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues)); 241ebda224cSfranklin5 } 2421c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N))); 243ebda224cSfranklin5 rstart = N->cmap->rstart; 244ebda224cSfranklin5 rend = N->cmap->rend; 2459566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &values)); 2469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart)); 2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &values)); 2489566063dSJacob Faibussowitsch PetscCall(PetscFree2(diag, work)); 2499566063dSJacob Faibussowitsch PetscCall(VecScale(v, Na->scale)); 250ebda224cSfranklin5 PetscFunctionReturn(0); 251ebda224cSfranklin5 } 252ebda224cSfranklin5 253d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N, Mat *D) 254d71ae5a4SJacob Faibussowitsch { 255a48507d3SJose E. Roman Mat_Normal *Na = (Mat_Normal *)N->data; 256a48507d3SJose E. Roman Mat M, A = Na->A; 257a48507d3SJose E. Roman 258a48507d3SJose E. Roman PetscFunctionBegin; 259a48507d3SJose E. Roman PetscCall(MatGetDiagonalBlock(A, &M)); 260a48507d3SJose E. Roman PetscCall(MatCreateNormalHermitian(M, &Na->D)); 261a48507d3SJose E. Roman *D = Na->D; 262a48507d3SJose E. Roman PetscFunctionReturn(0); 263a48507d3SJose E. Roman } 264a48507d3SJose E. Roman 265d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat_NormalHermitian(Mat A, Mat *M) 266d71ae5a4SJacob Faibussowitsch { 26765f45395SPierre Jolivet Mat_Normal *Aa = (Mat_Normal *)A->data; 26865f45395SPierre Jolivet 26965f45395SPierre Jolivet PetscFunctionBegin; 27065f45395SPierre Jolivet *M = Aa->A; 27165f45395SPierre Jolivet PetscFunctionReturn(0); 27265f45395SPierre Jolivet } 27365f45395SPierre Jolivet 27465f45395SPierre Jolivet /*@ 27511a5261eSBarry Smith MatNormalHermitianGetMat - Gets the `Mat` object stored inside a `MATNORMALHERMITIAN` 27665f45395SPierre Jolivet 277*c3339decSBarry Smith Logically collective 27865f45395SPierre Jolivet 27965f45395SPierre Jolivet Input Parameter: 28011a5261eSBarry Smith . A - the `MATNORMALHERMITIAN` matrix 28165f45395SPierre Jolivet 28265f45395SPierre Jolivet Output Parameter: 28365f45395SPierre Jolivet . M - the matrix object stored inside A 28465f45395SPierre Jolivet 28565f45395SPierre Jolivet Level: intermediate 28665f45395SPierre Jolivet 28711a5261eSBarry Smith .seealso: `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 28865f45395SPierre Jolivet @*/ 289d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalHermitianGetMat(Mat A, Mat *M) 290d71ae5a4SJacob Faibussowitsch { 29165f45395SPierre Jolivet PetscFunctionBegin; 29265f45395SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 29365f45395SPierre Jolivet PetscValidType(A, 1); 29465f45395SPierre Jolivet PetscValidPointer(M, 2); 295cac4c232SBarry Smith PetscUseMethod(A, "MatNormalGetMatHermitian_C", (Mat, Mat *), (A, M)); 29665f45395SPierre Jolivet PetscFunctionReturn(0); 29765f45395SPierre Jolivet } 29865f45395SPierre Jolivet 299d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_NormalHermitian_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 300d71ae5a4SJacob Faibussowitsch { 301ab704866SPierre Jolivet Mat_Normal *Aa = (Mat_Normal *)A->data; 302ab704866SPierre Jolivet Mat B, conjugate; 303ab704866SPierre Jolivet PetscInt m, n, M, N; 304ab704866SPierre Jolivet 305ab704866SPierre Jolivet PetscFunctionBegin; 306ab704866SPierre Jolivet PetscCall(MatGetSize(A, &M, &N)); 307ab704866SPierre Jolivet PetscCall(MatGetLocalSize(A, &m, &n)); 308ab704866SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 309ab704866SPierre Jolivet B = *newmat; 310ab704866SPierre Jolivet PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B)); 311ab704866SPierre Jolivet } else { 312ab704866SPierre Jolivet PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B)); 313ab704866SPierre Jolivet PetscCall(MatProductSetType(B, MATPRODUCT_AtB)); 314ab704866SPierre Jolivet PetscCall(MatProductSetFromOptions(B)); 315ab704866SPierre Jolivet PetscCall(MatProductSymbolic(B)); 316ab704866SPierre Jolivet PetscCall(MatSetOption(B, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE)); 317ab704866SPierre Jolivet } 318ab704866SPierre Jolivet if (PetscDefined(USE_COMPLEX)) { 319ab704866SPierre Jolivet PetscCall(MatDuplicate(Aa->A, MAT_COPY_VALUES, &conjugate)); 320ab704866SPierre Jolivet PetscCall(MatConjugate(conjugate)); 321ab704866SPierre Jolivet PetscCall(MatProductReplaceMats(conjugate, Aa->A, NULL, B)); 322ab704866SPierre Jolivet } 323ab704866SPierre Jolivet PetscCall(MatProductNumeric(B)); 324ab704866SPierre Jolivet if (PetscDefined(USE_COMPLEX)) PetscCall(MatDestroy(&conjugate)); 325ab704866SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 326ab704866SPierre Jolivet PetscCall(MatHeaderReplace(A, &B)); 327ab704866SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 328ab704866SPierre Jolivet PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat)); 329ab704866SPierre Jolivet PetscFunctionReturn(0); 330ab704866SPierre Jolivet } 331ab704866SPierre Jolivet 332ebda224cSfranklin5 /*@ 33311a5261eSBarry Smith MatCreateNormalHermitian - Creates a new matrix object `MATNORMALHERMITIAN` that behaves like (A*)'*A. 334ebda224cSfranklin5 335*c3339decSBarry Smith Collective 336ebda224cSfranklin5 337ebda224cSfranklin5 Input Parameter: 338ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix 339ebda224cSfranklin5 340ebda224cSfranklin5 Output Parameter: 341ebda224cSfranklin5 . N - the matrix that represents (A*)'*A 342ebda224cSfranklin5 343ebda224cSfranklin5 Level: intermediate 344ebda224cSfranklin5 34511a5261eSBarry Smith Note: 34695452b02SPatrick Sanan The product (A*)'*A is NOT actually formed! Rather the new matrix 34711a5261eSBarry Smith object performs the matrix-vector product, `MatMult()`, by first multiplying by 348ebda224cSfranklin5 A and then (A*)' 34911a5261eSBarry Smith 35011a5261eSBarry Smith .seealso: `MATNORMAL`, `MATNORMALHERMITIAN`, `MatNormalHermitianGetMat()` 351ebda224cSfranklin5 @*/ 352d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormalHermitian(Mat A, Mat *N) 353d71ae5a4SJacob Faibussowitsch { 354ebda224cSfranklin5 PetscInt m, n; 355ebda224cSfranklin5 Mat_Normal *Na; 3565cb049f5SJose E. Roman VecType vtype; 357ebda224cSfranklin5 358ebda224cSfranklin5 PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 3609566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N)); 3619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N, n, n, PETSC_DECIDE, PETSC_DECIDE)); 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMALHERMITIAN)); 3639566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap)); 3649566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap)); 365ebda224cSfranklin5 3664dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na)); 367ebda224cSfranklin5 (*N)->data = (void *)Na; 3689566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 369ebda224cSfranklin5 Na->A = A; 370ebda224cSfranklin5 Na->scale = 1.0; 371ebda224cSfranklin5 3729566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &Na->w)); 373ebda224cSfranklin5 3742e5cc420SJose E. Roman (*N)->ops->destroy = MatDestroy_NormalHermitian; 3752e5cc420SJose E. Roman (*N)->ops->mult = MatMult_NormalHermitian; 376ebda224cSfranklin5 (*N)->ops->multtranspose = MatMultHermitianTranspose_Normal; 377ebda224cSfranklin5 (*N)->ops->multtransposeadd = MatMultHermitianTransposeAdd_Normal; 378ebda224cSfranklin5 (*N)->ops->multadd = MatMultHermitianAdd_Normal; 3792e5cc420SJose E. Roman (*N)->ops->getdiagonal = MatGetDiagonal_NormalHermitian; 3802e5cc420SJose E. Roman (*N)->ops->getdiagonalblock = MatGetDiagonalBlock_NormalHermitian; 3812e5cc420SJose E. Roman (*N)->ops->scale = MatScale_NormalHermitian; 3822e5cc420SJose E. Roman (*N)->ops->diagonalscale = MatDiagonalScale_NormalHermitian; 383ab704866SPierre Jolivet (*N)->ops->createsubmatrices = MatCreateSubMatrices_NormalHermitian; 384ab704866SPierre Jolivet (*N)->ops->permute = MatPermute_NormalHermitian; 385ab704866SPierre Jolivet (*N)->ops->duplicate = MatDuplicate_NormalHermitian; 386ab704866SPierre Jolivet (*N)->ops->copy = MatCopy_NormalHermitian; 387ebda224cSfranklin5 (*N)->assembled = PETSC_TRUE; 3885cb049f5SJose E. Roman (*N)->preallocated = PETSC_TRUE; 38912ab1b64SPierre Jolivet 3902e5cc420SJose E. Roman PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatNormalGetMatHermitian_C", MatNormalGetMat_NormalHermitian)); 391ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_seqaij_C", MatConvert_NormalHermitian_AIJ)); 392ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_mpiaij_C", MatConvert_NormalHermitian_AIJ)); 3939566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N, MAT_HERMITIAN, PETSC_TRUE)); 3949566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A, &vtype)); 3959566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N, vtype)); 3965cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE) 3979566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N, A->boundtocpu)); 3985cb049f5SJose E. Roman #endif 399ebda224cSfranklin5 PetscFunctionReturn(0); 400ebda224cSfranklin5 } 401