1ebda224cSfranklin5 2ebda224cSfranklin5 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3ebda224cSfranklin5 4ebda224cSfranklin5 typedef struct { 5ebda224cSfranklin5 Mat A; 6ebda224cSfranklin5 Vec w,left,right,leftwork,rightwork; 7ebda224cSfranklin5 PetscScalar scale; 8ebda224cSfranklin5 } Mat_Normal; 9ebda224cSfranklin5 10ebda224cSfranklin5 PetscErrorCode MatScaleHermitian_Normal(Mat inA,PetscScalar scale) 11ebda224cSfranklin5 { 12ebda224cSfranklin5 Mat_Normal *a = (Mat_Normal*)inA->data; 13ebda224cSfranklin5 14ebda224cSfranklin5 PetscFunctionBegin; 15ebda224cSfranklin5 a->scale *= scale; 16ebda224cSfranklin5 PetscFunctionReturn(0); 17ebda224cSfranklin5 } 18ebda224cSfranklin5 19ebda224cSfranklin5 PetscErrorCode MatDiagonalScaleHermitian_Normal(Mat inA,Vec left,Vec right) 20ebda224cSfranklin5 { 21ebda224cSfranklin5 Mat_Normal *a = (Mat_Normal*)inA->data; 22ebda224cSfranklin5 23ebda224cSfranklin5 PetscFunctionBegin; 24ebda224cSfranklin5 if (left) { 25ebda224cSfranklin5 if (!a->left) { 269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&a->left)); 279566063dSJacob Faibussowitsch PetscCall(VecCopy(left,a->left)); 28ebda224cSfranklin5 } else { 299566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->left,left,a->left)); 30ebda224cSfranklin5 } 31ebda224cSfranklin5 } 32ebda224cSfranklin5 if (right) { 33ebda224cSfranklin5 if (!a->right) { 349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&a->right)); 359566063dSJacob Faibussowitsch PetscCall(VecCopy(right,a->right)); 36ebda224cSfranklin5 } else { 379566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->right,right,a->right)); 38ebda224cSfranklin5 } 39ebda224cSfranklin5 } 40ebda224cSfranklin5 PetscFunctionReturn(0); 41ebda224cSfranklin5 } 42ebda224cSfranklin5 43ebda224cSfranklin5 PetscErrorCode MatMultHermitian_Normal(Mat N,Vec x,Vec y) 44ebda224cSfranklin5 { 45ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 46ebda224cSfranklin5 Vec in; 47ebda224cSfranklin5 48ebda224cSfranklin5 PetscFunctionBegin; 49ebda224cSfranklin5 in = x; 50ebda224cSfranklin5 if (Na->right) { 51ebda224cSfranklin5 if (!Na->rightwork) { 529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 53ebda224cSfranklin5 } 549566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 55ebda224cSfranklin5 in = Na->rightwork; 56ebda224cSfranklin5 } 579566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 589566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y)); 59ebda224cSfranklin5 if (Na->left) { 609566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,Na->left,y)); 61ebda224cSfranklin5 } 629566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 63ebda224cSfranklin5 PetscFunctionReturn(0); 64ebda224cSfranklin5 } 65ebda224cSfranklin5 66ebda224cSfranklin5 PetscErrorCode MatMultHermitianAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 67ebda224cSfranklin5 { 68ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 69ebda224cSfranklin5 Vec in; 70ebda224cSfranklin5 71ebda224cSfranklin5 PetscFunctionBegin; 72ebda224cSfranklin5 in = v1; 73ebda224cSfranklin5 if (Na->right) { 74ebda224cSfranklin5 if (!Na->rightwork) { 759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 76ebda224cSfranklin5 } 779566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 78ebda224cSfranklin5 in = Na->rightwork; 79ebda224cSfranklin5 } 809566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 819566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 82ebda224cSfranklin5 if (Na->left) { 839566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3)); 849566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->left,v3)); 859566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 86ebda224cSfranklin5 } else { 879566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3)); 88ebda224cSfranklin5 } 89ebda224cSfranklin5 PetscFunctionReturn(0); 90ebda224cSfranklin5 } 91ebda224cSfranklin5 92ebda224cSfranklin5 PetscErrorCode MatMultHermitianTranspose_Normal(Mat N,Vec x,Vec y) 93ebda224cSfranklin5 { 94ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 95ebda224cSfranklin5 Vec in; 96ebda224cSfranklin5 97ebda224cSfranklin5 PetscFunctionBegin; 98ebda224cSfranklin5 in = x; 99ebda224cSfranklin5 if (Na->left) { 100ebda224cSfranklin5 if (!Na->leftwork) { 1019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 102ebda224cSfranklin5 } 1039566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 104ebda224cSfranklin5 in = Na->leftwork; 105ebda224cSfranklin5 } 1069566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1079566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,y)); 108ebda224cSfranklin5 if (Na->right) { 1099566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,Na->right,y)); 110ebda224cSfranklin5 } 1119566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 112ebda224cSfranklin5 PetscFunctionReturn(0); 113ebda224cSfranklin5 } 114ebda224cSfranklin5 115ebda224cSfranklin5 PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 116ebda224cSfranklin5 { 117ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 118ebda224cSfranklin5 Vec in; 119ebda224cSfranklin5 120ebda224cSfranklin5 PetscFunctionBegin; 121ebda224cSfranklin5 in = v1; 122ebda224cSfranklin5 if (Na->left) { 123ebda224cSfranklin5 if (!Na->leftwork) { 1249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 125ebda224cSfranklin5 } 1269566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 127ebda224cSfranklin5 in = Na->leftwork; 128ebda224cSfranklin5 } 1299566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1309566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 131ebda224cSfranklin5 if (Na->right) { 1329566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Na->A,Na->w,v3)); 1339566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->right,v3)); 1349566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 135ebda224cSfranklin5 } else { 1369566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3)); 137ebda224cSfranklin5 } 138ebda224cSfranklin5 PetscFunctionReturn(0); 139ebda224cSfranklin5 } 140ebda224cSfranklin5 141ebda224cSfranklin5 PetscErrorCode MatDestroyHermitian_Normal(Mat N) 142ebda224cSfranklin5 { 143ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 144ebda224cSfranklin5 145ebda224cSfranklin5 PetscFunctionBegin; 1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 1479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->w)); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->left)); 1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->right)); 1509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->leftwork)); 1519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rightwork)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMatHermitian_C",NULL)); 154ebda224cSfranklin5 PetscFunctionReturn(0); 155ebda224cSfranklin5 } 156ebda224cSfranklin5 157ebda224cSfranklin5 /* 158ebda224cSfranklin5 Slow, nonscalable version 159ebda224cSfranklin5 */ 160ebda224cSfranklin5 PetscErrorCode MatGetDiagonalHermitian_Normal(Mat N,Vec v) 161ebda224cSfranklin5 { 162ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 163ebda224cSfranklin5 Mat A = Na->A; 164ebda224cSfranklin5 PetscInt i,j,rstart,rend,nnz; 165ebda224cSfranklin5 const PetscInt *cols; 166ebda224cSfranklin5 PetscScalar *diag,*work,*values; 167ebda224cSfranklin5 const PetscScalar *mvalues; 168ebda224cSfranklin5 169ebda224cSfranklin5 PetscFunctionBegin; 1709566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work)); 1719566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(work,A->cmap->N)); 1729566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 173ebda224cSfranklin5 for (i=rstart; i<rend; i++) { 1749566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nnz,&cols,&mvalues)); 175ebda224cSfranklin5 for (j=0; j<nnz; j++) { 176a8f6171dSBarry Smith work[cols[j]] += mvalues[j]*PetscConj(mvalues[j]); 177ebda224cSfranklin5 } 1789566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nnz,&cols,&mvalues)); 179ebda224cSfranklin5 } 1809566063dSJacob Faibussowitsch PetscCallMPI(MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N))); 181ebda224cSfranklin5 rstart = N->cmap->rstart; 182ebda224cSfranklin5 rend = N->cmap->rend; 1839566063dSJacob Faibussowitsch PetscCall(VecGetArray(v,&values)); 1849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(values,diag+rstart,rend-rstart)); 1859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&values)); 1869566063dSJacob Faibussowitsch PetscCall(PetscFree2(diag,work)); 1879566063dSJacob Faibussowitsch PetscCall(VecScale(v,Na->scale)); 188ebda224cSfranklin5 PetscFunctionReturn(0); 189ebda224cSfranklin5 } 190ebda224cSfranklin5 19165f45395SPierre Jolivet PetscErrorCode MatNormalGetMatHermitian_Normal(Mat A,Mat *M) 19265f45395SPierre Jolivet { 19365f45395SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 19465f45395SPierre Jolivet 19565f45395SPierre Jolivet PetscFunctionBegin; 19665f45395SPierre Jolivet *M = Aa->A; 19765f45395SPierre Jolivet PetscFunctionReturn(0); 19865f45395SPierre Jolivet } 19965f45395SPierre Jolivet 20065f45395SPierre Jolivet /*@ 20165f45395SPierre Jolivet MatNormalHermitianGetMat - Gets the Mat object stored inside a MATNORMALHERMITIAN 20265f45395SPierre Jolivet 20365f45395SPierre Jolivet Logically collective on Mat 20465f45395SPierre Jolivet 20565f45395SPierre Jolivet Input Parameter: 20665f45395SPierre Jolivet . A - the MATNORMALHERMITIAN matrix 20765f45395SPierre Jolivet 20865f45395SPierre Jolivet Output Parameter: 20965f45395SPierre Jolivet . M - the matrix object stored inside A 21065f45395SPierre Jolivet 21165f45395SPierre Jolivet Level: intermediate 21265f45395SPierre Jolivet 21365f45395SPierre Jolivet .seealso: MatCreateNormalHermitian() 21465f45395SPierre Jolivet 21565f45395SPierre Jolivet @*/ 21665f45395SPierre Jolivet PetscErrorCode MatNormalHermitianGetMat(Mat A,Mat *M) 21765f45395SPierre Jolivet { 21865f45395SPierre Jolivet PetscFunctionBegin; 21965f45395SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 22065f45395SPierre Jolivet PetscValidType(A,1); 22165f45395SPierre Jolivet PetscValidPointer(M,2); 222*cac4c232SBarry Smith PetscUseMethod(A,"MatNormalGetMatHermitian_C",(Mat,Mat*),(A,M)); 22365f45395SPierre Jolivet PetscFunctionReturn(0); 22465f45395SPierre Jolivet } 22565f45395SPierre Jolivet 226ebda224cSfranklin5 /*@ 227ebda224cSfranklin5 MatCreateNormalHermitian - Creates a new matrix object that behaves like (A*)'*A. 228ebda224cSfranklin5 229ebda224cSfranklin5 Collective on Mat 230ebda224cSfranklin5 231ebda224cSfranklin5 Input Parameter: 232ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix 233ebda224cSfranklin5 234ebda224cSfranklin5 Output Parameter: 235ebda224cSfranklin5 . N - the matrix that represents (A*)'*A 236ebda224cSfranklin5 237ebda224cSfranklin5 Level: intermediate 238ebda224cSfranklin5 23995452b02SPatrick Sanan Notes: 24095452b02SPatrick Sanan The product (A*)'*A is NOT actually formed! Rather the new matrix 241ebda224cSfranklin5 object performs the matrix-vector product by first multiplying by 242ebda224cSfranklin5 A and then (A*)' 243ebda224cSfranklin5 @*/ 244ebda224cSfranklin5 PetscErrorCode MatCreateNormalHermitian(Mat A,Mat *N) 245ebda224cSfranklin5 { 246ebda224cSfranklin5 PetscInt m,n; 247ebda224cSfranklin5 Mat_Normal *Na; 2485cb049f5SJose E. Roman VecType vtype; 249ebda224cSfranklin5 250ebda224cSfranklin5 PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,&n)); 2529566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),N)); 2539566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N,MATNORMALHERMITIAN)); 2559566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->rmap)); 2569566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->cmap)); 257ebda224cSfranklin5 2589566063dSJacob Faibussowitsch PetscCall(PetscNewLog(*N,&Na)); 259ebda224cSfranklin5 (*N)->data = (void*) Na; 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 261ebda224cSfranklin5 Na->A = A; 262ebda224cSfranklin5 Na->scale = 1.0; 263ebda224cSfranklin5 2649566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&Na->w)); 265ebda224cSfranklin5 266ebda224cSfranklin5 (*N)->ops->destroy = MatDestroyHermitian_Normal; 267ebda224cSfranklin5 (*N)->ops->mult = MatMultHermitian_Normal; 268ebda224cSfranklin5 (*N)->ops->multtranspose = MatMultHermitianTranspose_Normal; 269ebda224cSfranklin5 (*N)->ops->multtransposeadd = MatMultHermitianTransposeAdd_Normal; 270ebda224cSfranklin5 (*N)->ops->multadd = MatMultHermitianAdd_Normal; 271ebda224cSfranklin5 (*N)->ops->getdiagonal = MatGetDiagonalHermitian_Normal; 272ebda224cSfranklin5 (*N)->ops->scale = MatScaleHermitian_Normal; 273ebda224cSfranklin5 (*N)->ops->diagonalscale = MatDiagonalScaleHermitian_Normal; 274ebda224cSfranklin5 (*N)->assembled = PETSC_TRUE; 2755cb049f5SJose E. Roman (*N)->preallocated = PETSC_TRUE; 27612ab1b64SPierre Jolivet 2779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMatHermitian_C",MatNormalGetMatHermitian_Normal)); 2789566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N,MAT_HERMITIAN,PETSC_TRUE)); 2799566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A,&vtype)); 2809566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N,vtype)); 2815cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE) 2829566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N,A->boundtocpu)); 2835cb049f5SJose E. Roman #endif 284ebda224cSfranklin5 PetscFunctionReturn(0); 285ebda224cSfranklin5 } 286