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 112e5cc420SJose E. Roman PetscErrorCode MatScale_NormalHermitian(Mat inA,PetscScalar scale) 12ebda224cSfranklin5 { 13ebda224cSfranklin5 Mat_Normal *a = (Mat_Normal*)inA->data; 14ebda224cSfranklin5 15ebda224cSfranklin5 PetscFunctionBegin; 16ebda224cSfranklin5 a->scale *= scale; 17ebda224cSfranklin5 PetscFunctionReturn(0); 18ebda224cSfranklin5 } 19ebda224cSfranklin5 202e5cc420SJose E. Roman PetscErrorCode MatDiagonalScale_NormalHermitian(Mat inA,Vec left,Vec right) 21ebda224cSfranklin5 { 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 44*ab704866SPierre Jolivet PetscErrorCode MatCreateSubMatrices_NormalHermitian(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 45*ab704866SPierre Jolivet { 46*ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal*)mat->data; 47*ab704866SPierre Jolivet Mat B = a->A, *suba; 48*ab704866SPierre Jolivet IS *row; 49*ab704866SPierre Jolivet PetscInt M; 50*ab704866SPierre Jolivet 51*ab704866SPierre Jolivet PetscFunctionBegin; 52*ab704866SPierre Jolivet PetscCheck(!a->left && !a->right && irow == icol,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not implemented"); 53*ab704866SPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 54*ab704866SPierre Jolivet PetscCall(PetscCalloc1(n,submat)); 55*ab704866SPierre Jolivet } 56*ab704866SPierre Jolivet PetscCall(MatGetSize(B,&M,NULL)); 57*ab704866SPierre Jolivet PetscCall(PetscMalloc1(n,&row)); 58*ab704866SPierre Jolivet PetscCall(ISCreateStride(PETSC_COMM_SELF,M,0,1,&row[0])); 59*ab704866SPierre Jolivet PetscCall(ISSetIdentity(row[0])); 60*ab704866SPierre Jolivet for (M = 1; M < n; ++M) row[M] = row[0]; 61*ab704866SPierre Jolivet PetscCall(MatCreateSubMatrices(B,n,row,icol,MAT_INITIAL_MATRIX,&suba)); 62*ab704866SPierre Jolivet for (M = 0; M < n; ++M) { 63*ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(suba[M],*submat+M)); 64*ab704866SPierre Jolivet ((Mat_Normal*)(*submat)[M]->data)->scale = a->scale; 65*ab704866SPierre Jolivet } 66*ab704866SPierre Jolivet PetscCall(ISDestroy(&row[0])); 67*ab704866SPierre Jolivet PetscCall(PetscFree(row)); 68*ab704866SPierre Jolivet PetscCall(MatDestroySubMatrices(n,&suba)); 69*ab704866SPierre Jolivet PetscFunctionReturn(0); 70*ab704866SPierre Jolivet } 71*ab704866SPierre Jolivet 72*ab704866SPierre Jolivet PetscErrorCode MatPermute_NormalHermitian(Mat A,IS rowp,IS colp,Mat *B) 73*ab704866SPierre Jolivet { 74*ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 75*ab704866SPierre Jolivet Mat C,Aa = a->A; 76*ab704866SPierre Jolivet IS row; 77*ab704866SPierre Jolivet 78*ab704866SPierre Jolivet PetscFunctionBegin; 79*ab704866SPierre Jolivet PetscCheck(rowp == colp,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Row permutation and column permutation must be the same"); 80*ab704866SPierre Jolivet PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa),Aa->rmap->n,Aa->rmap->rstart,1,&row)); 81*ab704866SPierre Jolivet PetscCall(ISSetIdentity(row)); 82*ab704866SPierre Jolivet PetscCall(MatPermute(Aa,row,colp,&C)); 83*ab704866SPierre Jolivet PetscCall(ISDestroy(&row)); 84*ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(C,B)); 85*ab704866SPierre Jolivet PetscCall(MatDestroy(&C)); 86*ab704866SPierre Jolivet PetscFunctionReturn(0); 87*ab704866SPierre Jolivet } 88*ab704866SPierre Jolivet 89*ab704866SPierre Jolivet PetscErrorCode MatDuplicate_NormalHermitian(Mat A, MatDuplicateOption op, Mat *B) 90*ab704866SPierre Jolivet { 91*ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 92*ab704866SPierre Jolivet Mat C; 93*ab704866SPierre Jolivet 94*ab704866SPierre Jolivet PetscFunctionBegin; 95*ab704866SPierre Jolivet PetscCheck(!a->left && !a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 96*ab704866SPierre Jolivet PetscCall(MatDuplicate(a->A,op,&C)); 97*ab704866SPierre Jolivet PetscCall(MatCreateNormalHermitian(C,B)); 98*ab704866SPierre Jolivet PetscCall(MatDestroy(&C)); 99*ab704866SPierre Jolivet if (op == MAT_COPY_VALUES) ((Mat_Normal*)(*B)->data)->scale = a->scale; 100*ab704866SPierre Jolivet PetscFunctionReturn(0); 101*ab704866SPierre Jolivet } 102*ab704866SPierre Jolivet 103*ab704866SPierre Jolivet PetscErrorCode MatCopy_NormalHermitian(Mat A,Mat B,MatStructure str) 104*ab704866SPierre Jolivet { 105*ab704866SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data,*b = (Mat_Normal*)B->data; 106*ab704866SPierre Jolivet 107*ab704866SPierre Jolivet PetscFunctionBegin; 108*ab704866SPierre Jolivet PetscCheck(!a->left && !a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 109*ab704866SPierre Jolivet PetscCall(MatCopy(a->A,b->A,str)); 110*ab704866SPierre Jolivet b->scale = a->scale; 111*ab704866SPierre Jolivet PetscCall(VecDestroy(&b->left)); 112*ab704866SPierre Jolivet PetscCall(VecDestroy(&b->right)); 113*ab704866SPierre Jolivet PetscCall(VecDestroy(&b->leftwork)); 114*ab704866SPierre Jolivet PetscCall(VecDestroy(&b->rightwork)); 115*ab704866SPierre Jolivet PetscFunctionReturn(0); 116*ab704866SPierre Jolivet } 117*ab704866SPierre Jolivet 1182e5cc420SJose E. Roman PetscErrorCode MatMult_NormalHermitian(Mat N,Vec x,Vec y) 119ebda224cSfranklin5 { 120ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 121ebda224cSfranklin5 Vec in; 122ebda224cSfranklin5 123ebda224cSfranklin5 PetscFunctionBegin; 124ebda224cSfranklin5 in = x; 125ebda224cSfranklin5 if (Na->right) { 126ebda224cSfranklin5 if (!Na->rightwork) { 1279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 128ebda224cSfranklin5 } 1299566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 130ebda224cSfranklin5 in = Na->rightwork; 131ebda224cSfranklin5 } 1329566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1339566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y)); 1341baa6e33SBarry Smith if (Na->left) PetscCall(VecPointwiseMult(y,Na->left,y)); 1359566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 136ebda224cSfranklin5 PetscFunctionReturn(0); 137ebda224cSfranklin5 } 138ebda224cSfranklin5 139ebda224cSfranklin5 PetscErrorCode MatMultHermitianAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 140ebda224cSfranklin5 { 141ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 142ebda224cSfranklin5 Vec in; 143ebda224cSfranklin5 144ebda224cSfranklin5 PetscFunctionBegin; 145ebda224cSfranklin5 in = v1; 146ebda224cSfranklin5 if (Na->right) { 147ebda224cSfranklin5 if (!Na->rightwork) { 1489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 149ebda224cSfranklin5 } 1509566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 151ebda224cSfranklin5 in = Na->rightwork; 152ebda224cSfranklin5 } 1539566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1549566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 155ebda224cSfranklin5 if (Na->left) { 1569566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3)); 1579566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->left,v3)); 1589566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 159ebda224cSfranklin5 } else { 1609566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3)); 161ebda224cSfranklin5 } 162ebda224cSfranklin5 PetscFunctionReturn(0); 163ebda224cSfranklin5 } 164ebda224cSfranklin5 165ebda224cSfranklin5 PetscErrorCode MatMultHermitianTranspose_Normal(Mat N,Vec x,Vec y) 166ebda224cSfranklin5 { 167ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 168ebda224cSfranklin5 Vec in; 169ebda224cSfranklin5 170ebda224cSfranklin5 PetscFunctionBegin; 171ebda224cSfranklin5 in = x; 172ebda224cSfranklin5 if (Na->left) { 173ebda224cSfranklin5 if (!Na->leftwork) { 1749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 175ebda224cSfranklin5 } 1769566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 177ebda224cSfranklin5 in = Na->leftwork; 178ebda224cSfranklin5 } 1799566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1809566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y)); 1811baa6e33SBarry Smith if (Na->right) PetscCall(VecPointwiseMult(y,Na->right,y)); 1829566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 183ebda224cSfranklin5 PetscFunctionReturn(0); 184ebda224cSfranklin5 } 185ebda224cSfranklin5 186ebda224cSfranklin5 PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 187ebda224cSfranklin5 { 188ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 189ebda224cSfranklin5 Vec in; 190ebda224cSfranklin5 191ebda224cSfranklin5 PetscFunctionBegin; 192ebda224cSfranklin5 in = v1; 193ebda224cSfranklin5 if (Na->left) { 194ebda224cSfranklin5 if (!Na->leftwork) { 1959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 196ebda224cSfranklin5 } 1979566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 198ebda224cSfranklin5 in = Na->leftwork; 199ebda224cSfranklin5 } 2009566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 2019566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 202ebda224cSfranklin5 if (Na->right) { 2039566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3)); 2049566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->right,v3)); 2059566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 206ebda224cSfranklin5 } else { 2079566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3)); 208ebda224cSfranklin5 } 209ebda224cSfranklin5 PetscFunctionReturn(0); 210ebda224cSfranklin5 } 211ebda224cSfranklin5 2122e5cc420SJose E. Roman PetscErrorCode MatDestroy_NormalHermitian(Mat N) 213ebda224cSfranklin5 { 214ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 215ebda224cSfranklin5 216ebda224cSfranklin5 PetscFunctionBegin; 2179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 218a48507d3SJose E. Roman PetscCall(MatDestroy(&Na->D)); 2199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->w)); 2209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->left)); 2219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->right)); 2229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->leftwork)); 2239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rightwork)); 2249566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMatHermitian_C",NULL)); 226*ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normalh_seqaij_C",NULL)); 227*ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normalh_mpiaij_C",NULL)); 228ebda224cSfranklin5 PetscFunctionReturn(0); 229ebda224cSfranklin5 } 230ebda224cSfranklin5 231ebda224cSfranklin5 /* 232ebda224cSfranklin5 Slow, nonscalable version 233ebda224cSfranklin5 */ 2342e5cc420SJose E. Roman PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N,Vec v) 235ebda224cSfranklin5 { 236ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 237ebda224cSfranklin5 Mat A = Na->A; 238ebda224cSfranklin5 PetscInt i,j,rstart,rend,nnz; 239ebda224cSfranklin5 const PetscInt *cols; 240ebda224cSfranklin5 PetscScalar *diag,*work,*values; 241ebda224cSfranklin5 const PetscScalar *mvalues; 242ebda224cSfranklin5 243ebda224cSfranklin5 PetscFunctionBegin; 2449566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work)); 2459566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(work,A->cmap->N)); 2469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 247ebda224cSfranklin5 for (i=rstart; i<rend; i++) { 2489566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nnz,&cols,&mvalues)); 249ebda224cSfranklin5 for (j=0; j<nnz; j++) { 250a8f6171dSBarry Smith work[cols[j]] += mvalues[j]*PetscConj(mvalues[j]); 251ebda224cSfranklin5 } 2529566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nnz,&cols,&mvalues)); 253ebda224cSfranklin5 } 2541c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N))); 255ebda224cSfranklin5 rstart = N->cmap->rstart; 256ebda224cSfranklin5 rend = N->cmap->rend; 2579566063dSJacob Faibussowitsch PetscCall(VecGetArray(v,&values)); 2589566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(values,diag+rstart,rend-rstart)); 2599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&values)); 2609566063dSJacob Faibussowitsch PetscCall(PetscFree2(diag,work)); 2619566063dSJacob Faibussowitsch PetscCall(VecScale(v,Na->scale)); 262ebda224cSfranklin5 PetscFunctionReturn(0); 263ebda224cSfranklin5 } 264ebda224cSfranklin5 2652e5cc420SJose E. Roman PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N,Mat *D) 266a48507d3SJose E. Roman { 267a48507d3SJose E. Roman Mat_Normal *Na = (Mat_Normal*)N->data; 268a48507d3SJose E. Roman Mat M,A = Na->A; 269a48507d3SJose E. Roman 270a48507d3SJose E. Roman PetscFunctionBegin; 271a48507d3SJose E. Roman PetscCall(MatGetDiagonalBlock(A,&M)); 272a48507d3SJose E. Roman PetscCall(MatCreateNormalHermitian(M,&Na->D)); 273a48507d3SJose E. Roman *D = Na->D; 274a48507d3SJose E. Roman PetscFunctionReturn(0); 275a48507d3SJose E. Roman } 276a48507d3SJose E. Roman 2772e5cc420SJose E. Roman PetscErrorCode MatNormalGetMat_NormalHermitian(Mat A,Mat *M) 27865f45395SPierre Jolivet { 27965f45395SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 28065f45395SPierre Jolivet 28165f45395SPierre Jolivet PetscFunctionBegin; 28265f45395SPierre Jolivet *M = Aa->A; 28365f45395SPierre Jolivet PetscFunctionReturn(0); 28465f45395SPierre Jolivet } 28565f45395SPierre Jolivet 28665f45395SPierre Jolivet /*@ 28765f45395SPierre Jolivet MatNormalHermitianGetMat - Gets the Mat object stored inside a MATNORMALHERMITIAN 28865f45395SPierre Jolivet 28965f45395SPierre Jolivet Logically collective on Mat 29065f45395SPierre Jolivet 29165f45395SPierre Jolivet Input Parameter: 29265f45395SPierre Jolivet . A - the MATNORMALHERMITIAN matrix 29365f45395SPierre Jolivet 29465f45395SPierre Jolivet Output Parameter: 29565f45395SPierre Jolivet . M - the matrix object stored inside A 29665f45395SPierre Jolivet 29765f45395SPierre Jolivet Level: intermediate 29865f45395SPierre Jolivet 299db781477SPatrick Sanan .seealso: `MatCreateNormalHermitian()` 30065f45395SPierre Jolivet 30165f45395SPierre Jolivet @*/ 30265f45395SPierre Jolivet PetscErrorCode MatNormalHermitianGetMat(Mat A,Mat *M) 30365f45395SPierre Jolivet { 30465f45395SPierre Jolivet PetscFunctionBegin; 30565f45395SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 30665f45395SPierre Jolivet PetscValidType(A,1); 30765f45395SPierre Jolivet PetscValidPointer(M,2); 308cac4c232SBarry Smith PetscUseMethod(A,"MatNormalGetMatHermitian_C",(Mat,Mat*),(A,M)); 30965f45395SPierre Jolivet PetscFunctionReturn(0); 31065f45395SPierre Jolivet } 31165f45395SPierre Jolivet 312*ab704866SPierre Jolivet PetscErrorCode MatConvert_NormalHermitian_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 313*ab704866SPierre Jolivet { 314*ab704866SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 315*ab704866SPierre Jolivet Mat B,conjugate; 316*ab704866SPierre Jolivet PetscInt m,n,M,N; 317*ab704866SPierre Jolivet 318*ab704866SPierre Jolivet PetscFunctionBegin; 319*ab704866SPierre Jolivet PetscCall(MatGetSize(A,&M,&N)); 320*ab704866SPierre Jolivet PetscCall(MatGetLocalSize(A,&m,&n)); 321*ab704866SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 322*ab704866SPierre Jolivet B = *newmat; 323*ab704866SPierre Jolivet PetscCall(MatProductReplaceMats(Aa->A,Aa->A,NULL,B)); 324*ab704866SPierre Jolivet } else { 325*ab704866SPierre Jolivet PetscCall(MatProductCreate(Aa->A,Aa->A,NULL,&B)); 326*ab704866SPierre Jolivet PetscCall(MatProductSetType(B,MATPRODUCT_AtB)); 327*ab704866SPierre Jolivet PetscCall(MatProductSetFromOptions(B)); 328*ab704866SPierre Jolivet PetscCall(MatProductSymbolic(B)); 329*ab704866SPierre Jolivet PetscCall(MatSetOption(B,!PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN,PETSC_TRUE)); 330*ab704866SPierre Jolivet } 331*ab704866SPierre Jolivet if (PetscDefined(USE_COMPLEX)) { 332*ab704866SPierre Jolivet PetscCall(MatDuplicate(Aa->A,MAT_COPY_VALUES,&conjugate)); 333*ab704866SPierre Jolivet PetscCall(MatConjugate(conjugate)); 334*ab704866SPierre Jolivet PetscCall(MatProductReplaceMats(conjugate,Aa->A,NULL,B)); 335*ab704866SPierre Jolivet } 336*ab704866SPierre Jolivet PetscCall(MatProductNumeric(B)); 337*ab704866SPierre Jolivet if (PetscDefined(USE_COMPLEX)) PetscCall(MatDestroy(&conjugate)); 338*ab704866SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 339*ab704866SPierre Jolivet PetscCall(MatHeaderReplace(A,&B)); 340*ab704866SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 341*ab704866SPierre Jolivet PetscCall(MatConvert(*newmat,MATAIJ,MAT_INPLACE_MATRIX,newmat)); 342*ab704866SPierre Jolivet PetscFunctionReturn(0); 343*ab704866SPierre Jolivet } 344*ab704866SPierre Jolivet 345ebda224cSfranklin5 /*@ 346ebda224cSfranklin5 MatCreateNormalHermitian - Creates a new matrix object that behaves like (A*)'*A. 347ebda224cSfranklin5 348ebda224cSfranklin5 Collective on Mat 349ebda224cSfranklin5 350ebda224cSfranklin5 Input Parameter: 351ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix 352ebda224cSfranklin5 353ebda224cSfranklin5 Output Parameter: 354ebda224cSfranklin5 . N - the matrix that represents (A*)'*A 355ebda224cSfranklin5 356ebda224cSfranklin5 Level: intermediate 357ebda224cSfranklin5 35895452b02SPatrick Sanan Notes: 35995452b02SPatrick Sanan The product (A*)'*A is NOT actually formed! Rather the new matrix 360ebda224cSfranklin5 object performs the matrix-vector product by first multiplying by 361ebda224cSfranklin5 A and then (A*)' 362ebda224cSfranklin5 @*/ 363ebda224cSfranklin5 PetscErrorCode MatCreateNormalHermitian(Mat A,Mat *N) 364ebda224cSfranklin5 { 365ebda224cSfranklin5 PetscInt m,n; 366ebda224cSfranklin5 Mat_Normal *Na; 3675cb049f5SJose E. Roman VecType vtype; 368ebda224cSfranklin5 369ebda224cSfranklin5 PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,&n)); 3719566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),N)); 3729566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE)); 3739566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N,MATNORMALHERMITIAN)); 3749566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->rmap)); 3759566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->cmap)); 376ebda224cSfranklin5 3779566063dSJacob Faibussowitsch PetscCall(PetscNewLog(*N,&Na)); 378ebda224cSfranklin5 (*N)->data = (void*) Na; 3799566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 380ebda224cSfranklin5 Na->A = A; 381ebda224cSfranklin5 Na->scale = 1.0; 382ebda224cSfranklin5 3839566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&Na->w)); 384ebda224cSfranklin5 3852e5cc420SJose E. Roman (*N)->ops->destroy = MatDestroy_NormalHermitian; 3862e5cc420SJose E. Roman (*N)->ops->mult = MatMult_NormalHermitian; 387ebda224cSfranklin5 (*N)->ops->multtranspose = MatMultHermitianTranspose_Normal; 388ebda224cSfranklin5 (*N)->ops->multtransposeadd = MatMultHermitianTransposeAdd_Normal; 389ebda224cSfranklin5 (*N)->ops->multadd = MatMultHermitianAdd_Normal; 3902e5cc420SJose E. Roman (*N)->ops->getdiagonal = MatGetDiagonal_NormalHermitian; 3912e5cc420SJose E. Roman (*N)->ops->getdiagonalblock = MatGetDiagonalBlock_NormalHermitian; 3922e5cc420SJose E. Roman (*N)->ops->scale = MatScale_NormalHermitian; 3932e5cc420SJose E. Roman (*N)->ops->diagonalscale = MatDiagonalScale_NormalHermitian; 394*ab704866SPierre Jolivet (*N)->ops->createsubmatrices= MatCreateSubMatrices_NormalHermitian; 395*ab704866SPierre Jolivet (*N)->ops->permute = MatPermute_NormalHermitian; 396*ab704866SPierre Jolivet (*N)->ops->duplicate = MatDuplicate_NormalHermitian; 397*ab704866SPierre Jolivet (*N)->ops->copy = MatCopy_NormalHermitian; 398ebda224cSfranklin5 (*N)->assembled = PETSC_TRUE; 3995cb049f5SJose E. Roman (*N)->preallocated = PETSC_TRUE; 40012ab1b64SPierre Jolivet 4012e5cc420SJose E. Roman PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMatHermitian_C",MatNormalGetMat_NormalHermitian)); 402*ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normalh_seqaij_C",MatConvert_NormalHermitian_AIJ)); 403*ab704866SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normalh_mpiaij_C",MatConvert_NormalHermitian_AIJ)); 4049566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N,MAT_HERMITIAN,PETSC_TRUE)); 4059566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A,&vtype)); 4069566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N,vtype)); 4075cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE) 4089566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N,A->boundtocpu)); 4095cb049f5SJose E. Roman #endif 410ebda224cSfranklin5 PetscFunctionReturn(0); 411ebda224cSfranklin5 } 412