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 PetscErrorCode ierr; 23ebda224cSfranklin5 24ebda224cSfranklin5 PetscFunctionBegin; 25ebda224cSfranklin5 if (left) { 26ebda224cSfranklin5 if (!a->left) { 27ebda224cSfranklin5 ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 28ebda224cSfranklin5 ierr = VecCopy(left,a->left);CHKERRQ(ierr); 29ebda224cSfranklin5 } else { 30ebda224cSfranklin5 ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 31ebda224cSfranklin5 } 32ebda224cSfranklin5 } 33ebda224cSfranklin5 if (right) { 34ebda224cSfranklin5 if (!a->right) { 35ebda224cSfranklin5 ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 36ebda224cSfranklin5 ierr = VecCopy(right,a->right);CHKERRQ(ierr); 37ebda224cSfranklin5 } else { 38ebda224cSfranklin5 ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 39ebda224cSfranklin5 } 40ebda224cSfranklin5 } 41ebda224cSfranklin5 PetscFunctionReturn(0); 42ebda224cSfranklin5 } 43ebda224cSfranklin5 44ebda224cSfranklin5 PetscErrorCode MatMultHermitian_Normal(Mat N,Vec x,Vec y) 45ebda224cSfranklin5 { 46ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 47ebda224cSfranklin5 PetscErrorCode ierr; 48ebda224cSfranklin5 Vec in; 49ebda224cSfranklin5 50ebda224cSfranklin5 PetscFunctionBegin; 51ebda224cSfranklin5 in = x; 52ebda224cSfranklin5 if (Na->right) { 53ebda224cSfranklin5 if (!Na->rightwork) { 54ebda224cSfranklin5 ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 55ebda224cSfranklin5 } 56ebda224cSfranklin5 ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 57ebda224cSfranklin5 in = Na->rightwork; 58ebda224cSfranklin5 } 59ebda224cSfranklin5 ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 60ebda224cSfranklin5 ierr = MatMultHermitianTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 61ebda224cSfranklin5 if (Na->left) { 62ebda224cSfranklin5 ierr = VecPointwiseMult(y,Na->left,y);CHKERRQ(ierr); 63ebda224cSfranklin5 } 64ebda224cSfranklin5 ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 65ebda224cSfranklin5 PetscFunctionReturn(0); 66ebda224cSfranklin5 } 67ebda224cSfranklin5 68ebda224cSfranklin5 PetscErrorCode MatMultHermitianAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 69ebda224cSfranklin5 { 70ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 71ebda224cSfranklin5 PetscErrorCode ierr; 72ebda224cSfranklin5 Vec in; 73ebda224cSfranklin5 74ebda224cSfranklin5 PetscFunctionBegin; 75ebda224cSfranklin5 in = v1; 76ebda224cSfranklin5 if (Na->right) { 77ebda224cSfranklin5 if (!Na->rightwork) { 78ebda224cSfranklin5 ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 79ebda224cSfranklin5 } 80ebda224cSfranklin5 ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 81ebda224cSfranklin5 in = Na->rightwork; 82ebda224cSfranklin5 } 83ebda224cSfranklin5 ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 84ebda224cSfranklin5 ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 85ebda224cSfranklin5 if (Na->left) { 86ebda224cSfranklin5 ierr = MatMultHermitianTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 87ebda224cSfranklin5 ierr = VecPointwiseMult(v3,Na->left,v3);CHKERRQ(ierr); 88ebda224cSfranklin5 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 89ebda224cSfranklin5 } else { 90ebda224cSfranklin5 ierr = MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 91ebda224cSfranklin5 } 92ebda224cSfranklin5 PetscFunctionReturn(0); 93ebda224cSfranklin5 } 94ebda224cSfranklin5 95ebda224cSfranklin5 PetscErrorCode MatMultHermitianTranspose_Normal(Mat N,Vec x,Vec y) 96ebda224cSfranklin5 { 97ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 98ebda224cSfranklin5 PetscErrorCode ierr; 99ebda224cSfranklin5 Vec in; 100ebda224cSfranklin5 101ebda224cSfranklin5 PetscFunctionBegin; 102ebda224cSfranklin5 in = x; 103ebda224cSfranklin5 if (Na->left) { 104ebda224cSfranklin5 if (!Na->leftwork) { 105ebda224cSfranklin5 ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 106ebda224cSfranklin5 } 107ebda224cSfranklin5 ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 108ebda224cSfranklin5 in = Na->leftwork; 109ebda224cSfranklin5 } 110ebda224cSfranklin5 ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 111ebda224cSfranklin5 ierr = MatMultHermitianTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 112ebda224cSfranklin5 if (Na->right) { 113ebda224cSfranklin5 ierr = VecPointwiseMult(y,Na->right,y);CHKERRQ(ierr); 114ebda224cSfranklin5 } 115ebda224cSfranklin5 ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 116ebda224cSfranklin5 PetscFunctionReturn(0); 117ebda224cSfranklin5 } 118ebda224cSfranklin5 119ebda224cSfranklin5 PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 120ebda224cSfranklin5 { 121ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 122ebda224cSfranklin5 PetscErrorCode ierr; 123ebda224cSfranklin5 Vec in; 124ebda224cSfranklin5 125ebda224cSfranklin5 PetscFunctionBegin; 126ebda224cSfranklin5 in = v1; 127ebda224cSfranklin5 if (Na->left) { 128ebda224cSfranklin5 if (!Na->leftwork) { 129ebda224cSfranklin5 ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 130ebda224cSfranklin5 } 131ebda224cSfranklin5 ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 132ebda224cSfranklin5 in = Na->leftwork; 133ebda224cSfranklin5 } 134ebda224cSfranklin5 ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 135ebda224cSfranklin5 ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 136ebda224cSfranklin5 if (Na->right) { 137ebda224cSfranklin5 ierr = MatMultHermitianTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 138ebda224cSfranklin5 ierr = VecPointwiseMult(v3,Na->right,v3);CHKERRQ(ierr); 139ebda224cSfranklin5 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 140ebda224cSfranklin5 } else { 141ebda224cSfranklin5 ierr = MatMultHermitianTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 142ebda224cSfranklin5 } 143ebda224cSfranklin5 PetscFunctionReturn(0); 144ebda224cSfranklin5 } 145ebda224cSfranklin5 146ebda224cSfranklin5 PetscErrorCode MatDestroyHermitian_Normal(Mat N) 147ebda224cSfranklin5 { 148ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 149ebda224cSfranklin5 PetscErrorCode ierr; 150ebda224cSfranklin5 151ebda224cSfranklin5 PetscFunctionBegin; 152ebda224cSfranklin5 ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 153ebda224cSfranklin5 ierr = VecDestroy(&Na->w);CHKERRQ(ierr); 154ebda224cSfranklin5 ierr = VecDestroy(&Na->left);CHKERRQ(ierr); 155ebda224cSfranklin5 ierr = VecDestroy(&Na->right);CHKERRQ(ierr); 156ebda224cSfranklin5 ierr = VecDestroy(&Na->leftwork);CHKERRQ(ierr); 157ebda224cSfranklin5 ierr = VecDestroy(&Na->rightwork);CHKERRQ(ierr); 158ebda224cSfranklin5 ierr = PetscFree(N->data);CHKERRQ(ierr); 159ebda224cSfranklin5 PetscFunctionReturn(0); 160ebda224cSfranklin5 } 161ebda224cSfranklin5 162ebda224cSfranklin5 /* 163ebda224cSfranklin5 Slow, nonscalable version 164ebda224cSfranklin5 */ 165ebda224cSfranklin5 PetscErrorCode MatGetDiagonalHermitian_Normal(Mat N,Vec v) 166ebda224cSfranklin5 { 167ebda224cSfranklin5 Mat_Normal *Na = (Mat_Normal*)N->data; 168ebda224cSfranklin5 Mat A = Na->A; 169ebda224cSfranklin5 PetscErrorCode ierr; 170ebda224cSfranklin5 PetscInt i,j,rstart,rend,nnz; 171ebda224cSfranklin5 const PetscInt *cols; 172ebda224cSfranklin5 PetscScalar *diag,*work,*values; 173ebda224cSfranklin5 const PetscScalar *mvalues; 174ebda224cSfranklin5 175ebda224cSfranklin5 PetscFunctionBegin; 176ebda224cSfranklin5 ierr = PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work);CHKERRQ(ierr); 177ebda224cSfranklin5 ierr = PetscMemzero(work,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 178ebda224cSfranklin5 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 179ebda224cSfranklin5 for (i=rstart; i<rend; i++) { 180ebda224cSfranklin5 ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 181ebda224cSfranklin5 for (j=0; j<nnz; j++) { 182a8f6171dSBarry Smith work[cols[j]] += mvalues[j]*PetscConj(mvalues[j]); 183ebda224cSfranklin5 } 184ebda224cSfranklin5 ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 185ebda224cSfranklin5 } 186b2566f29SBarry Smith ierr = MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr); 187ebda224cSfranklin5 rstart = N->cmap->rstart; 188ebda224cSfranklin5 rend = N->cmap->rend; 189ebda224cSfranklin5 ierr = VecGetArray(v,&values);CHKERRQ(ierr); 190ebda224cSfranklin5 ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr); 191ebda224cSfranklin5 ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 192ebda224cSfranklin5 ierr = PetscFree2(diag,work);CHKERRQ(ierr); 193ebda224cSfranklin5 ierr = VecScale(v,Na->scale);CHKERRQ(ierr); 194ebda224cSfranklin5 PetscFunctionReturn(0); 195ebda224cSfranklin5 } 196ebda224cSfranklin5 197ebda224cSfranklin5 /*@ 198ebda224cSfranklin5 MatCreateNormalHermitian - Creates a new matrix object that behaves like (A*)'*A. 199ebda224cSfranklin5 200ebda224cSfranklin5 Collective on Mat 201ebda224cSfranklin5 202ebda224cSfranklin5 Input Parameter: 203ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix 204ebda224cSfranklin5 205ebda224cSfranklin5 Output Parameter: 206ebda224cSfranklin5 . N - the matrix that represents (A*)'*A 207ebda224cSfranklin5 208ebda224cSfranklin5 Level: intermediate 209ebda224cSfranklin5 210*95452b02SPatrick Sanan Notes: 211*95452b02SPatrick Sanan The product (A*)'*A is NOT actually formed! Rather the new matrix 212ebda224cSfranklin5 object performs the matrix-vector product by first multiplying by 213ebda224cSfranklin5 A and then (A*)' 214ebda224cSfranklin5 @*/ 215ebda224cSfranklin5 PetscErrorCode MatCreateNormalHermitian(Mat A,Mat *N) 216ebda224cSfranklin5 { 217ebda224cSfranklin5 PetscErrorCode ierr; 218ebda224cSfranklin5 PetscInt m,n; 219ebda224cSfranklin5 Mat_Normal *Na; 220ebda224cSfranklin5 221ebda224cSfranklin5 PetscFunctionBegin; 222ebda224cSfranklin5 ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 223ebda224cSfranklin5 ierr = MatCreate(PetscObjectComm((PetscObject)A),N);CHKERRQ(ierr); 224ebda224cSfranklin5 ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 225ebda224cSfranklin5 ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMALHERMITIAN);CHKERRQ(ierr); 226ebda224cSfranklin5 227ebda224cSfranklin5 ierr = PetscNewLog(*N,&Na);CHKERRQ(ierr); 228ebda224cSfranklin5 (*N)->data = (void*) Na; 229ebda224cSfranklin5 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 230ebda224cSfranklin5 Na->A = A; 231ebda224cSfranklin5 Na->scale = 1.0; 232ebda224cSfranklin5 233ebda224cSfranklin5 ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr); 234ebda224cSfranklin5 235ebda224cSfranklin5 (*N)->ops->destroy = MatDestroyHermitian_Normal; 236ebda224cSfranklin5 (*N)->ops->mult = MatMultHermitian_Normal; 237ebda224cSfranklin5 (*N)->ops->multtranspose = MatMultHermitianTranspose_Normal; 238ebda224cSfranklin5 (*N)->ops->multtransposeadd = MatMultHermitianTransposeAdd_Normal; 239ebda224cSfranklin5 (*N)->ops->multadd = MatMultHermitianAdd_Normal; 240ebda224cSfranklin5 (*N)->ops->getdiagonal = MatGetDiagonalHermitian_Normal; 241ebda224cSfranklin5 (*N)->ops->scale = MatScaleHermitian_Normal; 242ebda224cSfranklin5 (*N)->ops->diagonalscale = MatDiagonalScaleHermitian_Normal; 243ebda224cSfranklin5 (*N)->assembled = PETSC_TRUE; 244ebda224cSfranklin5 (*N)->cmap->N = A->cmap->N; 245ebda224cSfranklin5 (*N)->rmap->N = A->cmap->N; 246ebda224cSfranklin5 (*N)->cmap->n = A->cmap->n; 247ebda224cSfranklin5 (*N)->rmap->n = A->cmap->n; 248ebda224cSfranklin5 PetscFunctionReturn(0); 249ebda224cSfranklin5 } 250ebda224cSfranklin5 251