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