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 442e5cc420SJose E. Roman PetscErrorCode MatMult_NormalHermitian(Mat N,Vec x,Vec y) 45ebda224cSfranklin5 { 46ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 47ebda224cSfranklin5 Vec in; 48ebda224cSfranklin5 49ebda224cSfranklin5 PetscFunctionBegin; 50ebda224cSfranklin5 in = x; 51ebda224cSfranklin5 if (Na->right) { 52ebda224cSfranklin5 if (!Na->rightwork) { 539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 54ebda224cSfranklin5 } 559566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 56ebda224cSfranklin5 in = Na->rightwork; 57ebda224cSfranklin5 } 589566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 599566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y)); 60*1baa6e33SBarry Smith if (Na->left) PetscCall(VecPointwiseMult(y,Na->left,y)); 619566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 62ebda224cSfranklin5 PetscFunctionReturn(0); 63ebda224cSfranklin5 } 64ebda224cSfranklin5 65ebda224cSfranklin5 PetscErrorCode MatMultHermitianAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 66ebda224cSfranklin5 { 67ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 68ebda224cSfranklin5 Vec in; 69ebda224cSfranklin5 70ebda224cSfranklin5 PetscFunctionBegin; 71ebda224cSfranklin5 in = v1; 72ebda224cSfranklin5 if (Na->right) { 73ebda224cSfranklin5 if (!Na->rightwork) { 749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 75ebda224cSfranklin5 } 769566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 77ebda224cSfranklin5 in = Na->rightwork; 78ebda224cSfranklin5 } 799566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 809566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 81ebda224cSfranklin5 if (Na->left) { 829566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3)); 839566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->left,v3)); 849566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 85ebda224cSfranklin5 } else { 869566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3)); 87ebda224cSfranklin5 } 88ebda224cSfranklin5 PetscFunctionReturn(0); 89ebda224cSfranklin5 } 90ebda224cSfranklin5 91ebda224cSfranklin5 PetscErrorCode MatMultHermitianTranspose_Normal(Mat N,Vec x,Vec y) 92ebda224cSfranklin5 { 93ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 94ebda224cSfranklin5 Vec in; 95ebda224cSfranklin5 96ebda224cSfranklin5 PetscFunctionBegin; 97ebda224cSfranklin5 in = x; 98ebda224cSfranklin5 if (Na->left) { 99ebda224cSfranklin5 if (!Na->leftwork) { 1009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 101ebda224cSfranklin5 } 1029566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 103ebda224cSfranklin5 in = Na->leftwork; 104ebda224cSfranklin5 } 1059566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1069566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y)); 107*1baa6e33SBarry Smith if (Na->right) PetscCall(VecPointwiseMult(y,Na->right,y)); 1089566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 109ebda224cSfranklin5 PetscFunctionReturn(0); 110ebda224cSfranklin5 } 111ebda224cSfranklin5 112ebda224cSfranklin5 PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 113ebda224cSfranklin5 { 114ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 115ebda224cSfranklin5 Vec in; 116ebda224cSfranklin5 117ebda224cSfranklin5 PetscFunctionBegin; 118ebda224cSfranklin5 in = v1; 119ebda224cSfranklin5 if (Na->left) { 120ebda224cSfranklin5 if (!Na->leftwork) { 1219566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 122ebda224cSfranklin5 } 1239566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 124ebda224cSfranklin5 in = Na->leftwork; 125ebda224cSfranklin5 } 1269566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1279566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 128ebda224cSfranklin5 if (Na->right) { 1299566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3)); 1309566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->right,v3)); 1319566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 132ebda224cSfranklin5 } else { 1339566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3)); 134ebda224cSfranklin5 } 135ebda224cSfranklin5 PetscFunctionReturn(0); 136ebda224cSfranklin5 } 137ebda224cSfranklin5 1382e5cc420SJose E. Roman PetscErrorCode MatDestroy_NormalHermitian(Mat N) 139ebda224cSfranklin5 { 140ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 141ebda224cSfranklin5 142ebda224cSfranklin5 PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 144a48507d3SJose E. Roman PetscCall(MatDestroy(&Na->D)); 1459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->w)); 1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->left)); 1479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->right)); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->leftwork)); 1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rightwork)); 1509566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 1519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMatHermitian_C",NULL)); 152ebda224cSfranklin5 PetscFunctionReturn(0); 153ebda224cSfranklin5 } 154ebda224cSfranklin5 155ebda224cSfranklin5 /* 156ebda224cSfranklin5 Slow, nonscalable version 157ebda224cSfranklin5 */ 1582e5cc420SJose E. Roman PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N,Vec v) 159ebda224cSfranklin5 { 160ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 161ebda224cSfranklin5 Mat A = Na->A; 162ebda224cSfranklin5 PetscInt i,j,rstart,rend,nnz; 163ebda224cSfranklin5 const PetscInt *cols; 164ebda224cSfranklin5 PetscScalar *diag,*work,*values; 165ebda224cSfranklin5 const PetscScalar *mvalues; 166ebda224cSfranklin5 167ebda224cSfranklin5 PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work)); 1699566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(work,A->cmap->N)); 1709566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 171ebda224cSfranklin5 for (i=rstart; i<rend; i++) { 1729566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nnz,&cols,&mvalues)); 173ebda224cSfranklin5 for (j=0; j<nnz; j++) { 174a8f6171dSBarry Smith work[cols[j]] += mvalues[j]*PetscConj(mvalues[j]); 175ebda224cSfranklin5 } 1769566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nnz,&cols,&mvalues)); 177ebda224cSfranklin5 } 1781c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N))); 179ebda224cSfranklin5 rstart = N->cmap->rstart; 180ebda224cSfranklin5 rend = N->cmap->rend; 1819566063dSJacob Faibussowitsch PetscCall(VecGetArray(v,&values)); 1829566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(values,diag+rstart,rend-rstart)); 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&values)); 1849566063dSJacob Faibussowitsch PetscCall(PetscFree2(diag,work)); 1859566063dSJacob Faibussowitsch PetscCall(VecScale(v,Na->scale)); 186ebda224cSfranklin5 PetscFunctionReturn(0); 187ebda224cSfranklin5 } 188ebda224cSfranklin5 1892e5cc420SJose E. Roman PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N,Mat *D) 190a48507d3SJose E. Roman { 191a48507d3SJose E. Roman Mat_Normal *Na = (Mat_Normal*)N->data; 192a48507d3SJose E. Roman Mat M,A = Na->A; 193a48507d3SJose E. Roman 194a48507d3SJose E. Roman PetscFunctionBegin; 195a48507d3SJose E. Roman PetscCall(MatGetDiagonalBlock(A,&M)); 196a48507d3SJose E. Roman PetscCall(MatCreateNormalHermitian(M,&Na->D)); 197a48507d3SJose E. Roman *D = Na->D; 198a48507d3SJose E. Roman PetscFunctionReturn(0); 199a48507d3SJose E. Roman } 200a48507d3SJose E. Roman 2012e5cc420SJose E. Roman PetscErrorCode MatNormalGetMat_NormalHermitian(Mat A,Mat *M) 20265f45395SPierre Jolivet { 20365f45395SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 20465f45395SPierre Jolivet 20565f45395SPierre Jolivet PetscFunctionBegin; 20665f45395SPierre Jolivet *M = Aa->A; 20765f45395SPierre Jolivet PetscFunctionReturn(0); 20865f45395SPierre Jolivet } 20965f45395SPierre Jolivet 21065f45395SPierre Jolivet /*@ 21165f45395SPierre Jolivet MatNormalHermitianGetMat - Gets the Mat object stored inside a MATNORMALHERMITIAN 21265f45395SPierre Jolivet 21365f45395SPierre Jolivet Logically collective on Mat 21465f45395SPierre Jolivet 21565f45395SPierre Jolivet Input Parameter: 21665f45395SPierre Jolivet . A - the MATNORMALHERMITIAN matrix 21765f45395SPierre Jolivet 21865f45395SPierre Jolivet Output Parameter: 21965f45395SPierre Jolivet . M - the matrix object stored inside A 22065f45395SPierre Jolivet 22165f45395SPierre Jolivet Level: intermediate 22265f45395SPierre Jolivet 223db781477SPatrick Sanan .seealso: `MatCreateNormalHermitian()` 22465f45395SPierre Jolivet 22565f45395SPierre Jolivet @*/ 22665f45395SPierre Jolivet PetscErrorCode MatNormalHermitianGetMat(Mat A,Mat *M) 22765f45395SPierre Jolivet { 22865f45395SPierre Jolivet PetscFunctionBegin; 22965f45395SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 23065f45395SPierre Jolivet PetscValidType(A,1); 23165f45395SPierre Jolivet PetscValidPointer(M,2); 232cac4c232SBarry Smith PetscUseMethod(A,"MatNormalGetMatHermitian_C",(Mat,Mat*),(A,M)); 23365f45395SPierre Jolivet PetscFunctionReturn(0); 23465f45395SPierre Jolivet } 23565f45395SPierre Jolivet 236ebda224cSfranklin5 /*@ 237ebda224cSfranklin5 MatCreateNormalHermitian - Creates a new matrix object that behaves like (A*)'*A. 238ebda224cSfranklin5 239ebda224cSfranklin5 Collective on Mat 240ebda224cSfranklin5 241ebda224cSfranklin5 Input Parameter: 242ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix 243ebda224cSfranklin5 244ebda224cSfranklin5 Output Parameter: 245ebda224cSfranklin5 . N - the matrix that represents (A*)'*A 246ebda224cSfranklin5 247ebda224cSfranklin5 Level: intermediate 248ebda224cSfranklin5 24995452b02SPatrick Sanan Notes: 25095452b02SPatrick Sanan The product (A*)'*A is NOT actually formed! Rather the new matrix 251ebda224cSfranklin5 object performs the matrix-vector product by first multiplying by 252ebda224cSfranklin5 A and then (A*)' 253ebda224cSfranklin5 @*/ 254ebda224cSfranklin5 PetscErrorCode MatCreateNormalHermitian(Mat A,Mat *N) 255ebda224cSfranklin5 { 256ebda224cSfranklin5 PetscInt m,n; 257ebda224cSfranklin5 Mat_Normal *Na; 2585cb049f5SJose E. Roman VecType vtype; 259ebda224cSfranklin5 260ebda224cSfranklin5 PetscFunctionBegin; 2619566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,&n)); 2629566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),N)); 2639566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N,MATNORMALHERMITIAN)); 2659566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->rmap)); 2669566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->cmap)); 267ebda224cSfranklin5 2689566063dSJacob Faibussowitsch PetscCall(PetscNewLog(*N,&Na)); 269ebda224cSfranklin5 (*N)->data = (void*) Na; 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 271ebda224cSfranklin5 Na->A = A; 272ebda224cSfranklin5 Na->scale = 1.0; 273ebda224cSfranklin5 2749566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&Na->w)); 275ebda224cSfranklin5 2762e5cc420SJose E. Roman (*N)->ops->destroy = MatDestroy_NormalHermitian; 2772e5cc420SJose E. Roman (*N)->ops->mult = MatMult_NormalHermitian; 278ebda224cSfranklin5 (*N)->ops->multtranspose = MatMultHermitianTranspose_Normal; 279ebda224cSfranklin5 (*N)->ops->multtransposeadd = MatMultHermitianTransposeAdd_Normal; 280ebda224cSfranklin5 (*N)->ops->multadd = MatMultHermitianAdd_Normal; 2812e5cc420SJose E. Roman (*N)->ops->getdiagonal = MatGetDiagonal_NormalHermitian; 2822e5cc420SJose E. Roman (*N)->ops->getdiagonalblock = MatGetDiagonalBlock_NormalHermitian; 2832e5cc420SJose E. Roman (*N)->ops->scale = MatScale_NormalHermitian; 2842e5cc420SJose E. Roman (*N)->ops->diagonalscale = MatDiagonalScale_NormalHermitian; 285ebda224cSfranklin5 (*N)->assembled = PETSC_TRUE; 2865cb049f5SJose E. Roman (*N)->preallocated = PETSC_TRUE; 28712ab1b64SPierre Jolivet 2882e5cc420SJose E. Roman PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMatHermitian_C",MatNormalGetMat_NormalHermitian)); 2899566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N,MAT_HERMITIAN,PETSC_TRUE)); 2909566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A,&vtype)); 2919566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N,vtype)); 2925cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE) 2939566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N,A->boundtocpu)); 2945cb049f5SJose E. Roman #endif 295ebda224cSfranklin5 PetscFunctionReturn(0); 296ebda224cSfranklin5 } 297