1c8a8475eSBarry Smith 2b45d2f2cSJed Brown #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/ 3c8a8475eSBarry Smith 4c8a8475eSBarry Smith typedef struct { 57cfd0b8cSBarry Smith Mat A; 67cfd0b8cSBarry Smith Vec w,left,right,leftwork,rightwork; 7b20f02adSBarry Smith PetscScalar scale; 8c8a8475eSBarry Smith } Mat_Normal; 9c8a8475eSBarry Smith 10c8a8475eSBarry Smith #undef __FUNCT__ 1169ef1043SBarry Smith #define __FUNCT__ "MatScale_Normal" 1269ef1043SBarry Smith PetscErrorCode MatScale_Normal(Mat inA,PetscScalar scale) 1369ef1043SBarry Smith { 1469ef1043SBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 1569ef1043SBarry Smith 1669ef1043SBarry Smith PetscFunctionBegin; 1769ef1043SBarry Smith a->scale *= scale; 1869ef1043SBarry Smith PetscFunctionReturn(0); 1969ef1043SBarry Smith } 2069ef1043SBarry Smith 2169ef1043SBarry Smith #undef __FUNCT__ 227cfd0b8cSBarry Smith #define __FUNCT__ "MatDiagonalScale_Normal" 237cfd0b8cSBarry Smith PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right) 247cfd0b8cSBarry Smith { 257cfd0b8cSBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 267cfd0b8cSBarry Smith PetscErrorCode ierr; 277cfd0b8cSBarry Smith 287cfd0b8cSBarry Smith PetscFunctionBegin; 297cfd0b8cSBarry Smith if (left) { 307cfd0b8cSBarry Smith if (!a->left) { 317cfd0b8cSBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 327cfd0b8cSBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 337cfd0b8cSBarry Smith } else { 347cfd0b8cSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 357cfd0b8cSBarry Smith } 367cfd0b8cSBarry Smith } 377cfd0b8cSBarry Smith if (right) { 387cfd0b8cSBarry Smith if (!a->right) { 397cfd0b8cSBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 407cfd0b8cSBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 417cfd0b8cSBarry Smith } else { 427cfd0b8cSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 437cfd0b8cSBarry Smith } 447cfd0b8cSBarry Smith } 457cfd0b8cSBarry Smith PetscFunctionReturn(0); 467cfd0b8cSBarry Smith } 477cfd0b8cSBarry Smith 487cfd0b8cSBarry Smith #undef __FUNCT__ 49c8a8475eSBarry Smith #define __FUNCT__ "MatMult_Normal" 50dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y) 51c8a8475eSBarry Smith { 52c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 53dfbe8321SBarry Smith PetscErrorCode ierr; 547cfd0b8cSBarry Smith Vec in; 55c8a8475eSBarry Smith 56c8a8475eSBarry Smith PetscFunctionBegin; 577cfd0b8cSBarry Smith in = x; 587cfd0b8cSBarry Smith if (Na->right) { 597cfd0b8cSBarry Smith if (!Na->rightwork) { 607cfd0b8cSBarry Smith ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 617cfd0b8cSBarry Smith } 627cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 637cfd0b8cSBarry Smith in = Na->rightwork; 647cfd0b8cSBarry Smith } 657cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 66c8a8475eSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 677cfd0b8cSBarry Smith if (Na->left) { 687cfd0b8cSBarry Smith ierr = VecPointwiseMult(y,Na->left,y);CHKERRQ(ierr); 697cfd0b8cSBarry Smith } 70b20f02adSBarry Smith ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 71c8a8475eSBarry Smith PetscFunctionReturn(0); 72c8a8475eSBarry Smith } 73c8a8475eSBarry Smith 74c8a8475eSBarry Smith #undef __FUNCT__ 75db090513SMatthew Knepley #define __FUNCT__ "MatMultAdd_Normal" 76db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 77db090513SMatthew Knepley { 78db090513SMatthew Knepley Mat_Normal *Na = (Mat_Normal*)N->data; 79db090513SMatthew Knepley PetscErrorCode ierr; 807cfd0b8cSBarry Smith Vec in; 81db090513SMatthew Knepley 82db090513SMatthew Knepley PetscFunctionBegin; 837cfd0b8cSBarry Smith in = v1; 847cfd0b8cSBarry Smith if (Na->right) { 857cfd0b8cSBarry Smith if (!Na->rightwork) { 867cfd0b8cSBarry Smith ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 877cfd0b8cSBarry Smith } 887cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 897cfd0b8cSBarry Smith in = Na->rightwork; 907cfd0b8cSBarry Smith } 917cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 92b20f02adSBarry Smith ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 937cfd0b8cSBarry Smith if (Na->left) { 947cfd0b8cSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 957cfd0b8cSBarry Smith ierr = VecPointwiseMult(v3,Na->left,v3);CHKERRQ(ierr); 967cfd0b8cSBarry Smith ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 977cfd0b8cSBarry Smith } else { 98b20f02adSBarry Smith ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 997cfd0b8cSBarry Smith } 100b20f02adSBarry Smith PetscFunctionReturn(0); 101b20f02adSBarry Smith } 102b20f02adSBarry Smith 103b20f02adSBarry Smith #undef __FUNCT__ 104b20f02adSBarry Smith #define __FUNCT__ "MatMultTranspose_Normal" 105b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y) 106b20f02adSBarry Smith { 107b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 108b20f02adSBarry Smith PetscErrorCode ierr; 1097cfd0b8cSBarry Smith Vec in; 110b20f02adSBarry Smith 111b20f02adSBarry Smith PetscFunctionBegin; 1127cfd0b8cSBarry Smith in = x; 1137cfd0b8cSBarry Smith if (Na->left) { 1147cfd0b8cSBarry Smith if (!Na->leftwork) { 1157cfd0b8cSBarry Smith ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 1167cfd0b8cSBarry Smith } 1177cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 1187cfd0b8cSBarry Smith in = Na->leftwork; 1197cfd0b8cSBarry Smith } 1207cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 121b20f02adSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 1227cfd0b8cSBarry Smith if (Na->right) { 1237cfd0b8cSBarry Smith ierr = VecPointwiseMult(y,Na->right,y);CHKERRQ(ierr); 1247cfd0b8cSBarry Smith } 125b20f02adSBarry Smith ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 126b20f02adSBarry Smith PetscFunctionReturn(0); 127b20f02adSBarry Smith } 128b20f02adSBarry Smith 129b20f02adSBarry Smith #undef __FUNCT__ 130b20f02adSBarry Smith #define __FUNCT__ "MatMultTransposeAdd_Normal" 131b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 132b20f02adSBarry Smith { 133b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 134b20f02adSBarry Smith PetscErrorCode ierr; 1357cfd0b8cSBarry Smith Vec in; 136b20f02adSBarry Smith 137b20f02adSBarry Smith PetscFunctionBegin; 1387cfd0b8cSBarry Smith in = v1; 1397cfd0b8cSBarry Smith if (Na->left) { 1407cfd0b8cSBarry Smith if (!Na->leftwork) { 1417cfd0b8cSBarry Smith ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 1427cfd0b8cSBarry Smith } 1437cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 1447cfd0b8cSBarry Smith in = Na->leftwork; 1457cfd0b8cSBarry Smith } 1467cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 147b20f02adSBarry Smith ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 1487cfd0b8cSBarry Smith if (Na->right) { 1497cfd0b8cSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 1507cfd0b8cSBarry Smith ierr = VecPointwiseMult(v3,Na->right,v3);CHKERRQ(ierr); 1517cfd0b8cSBarry Smith ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1527cfd0b8cSBarry Smith } else { 153db090513SMatthew Knepley ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 1547cfd0b8cSBarry Smith } 155db090513SMatthew Knepley PetscFunctionReturn(0); 156db090513SMatthew Knepley } 157db090513SMatthew Knepley 158db090513SMatthew Knepley #undef __FUNCT__ 159c8a8475eSBarry Smith #define __FUNCT__ "MatDestroy_Normal" 160dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N) 161c8a8475eSBarry Smith { 162c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 163dfbe8321SBarry Smith PetscErrorCode ierr; 164c8a8475eSBarry Smith 165c8a8475eSBarry Smith PetscFunctionBegin; 1666bf464f9SBarry Smith ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 1676bf464f9SBarry Smith ierr = VecDestroy(&Na->w);CHKERRQ(ierr); 1686bf464f9SBarry Smith ierr = VecDestroy(&Na->left);CHKERRQ(ierr); 1696bf464f9SBarry Smith ierr = VecDestroy(&Na->right);CHKERRQ(ierr); 1706bf464f9SBarry Smith ierr = VecDestroy(&Na->leftwork);CHKERRQ(ierr); 1716bf464f9SBarry Smith ierr = VecDestroy(&Na->rightwork);CHKERRQ(ierr); 172bf0cc555SLisandro Dalcin ierr = PetscFree(N->data);CHKERRQ(ierr); 173c8a8475eSBarry Smith PetscFunctionReturn(0); 174c8a8475eSBarry Smith } 175c8a8475eSBarry Smith 17617a6fd94SBarry Smith /* 17717a6fd94SBarry Smith Slow, nonscalable version 17817a6fd94SBarry Smith */ 17917a6fd94SBarry Smith #undef __FUNCT__ 18017a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal" 181dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 18217a6fd94SBarry Smith { 18317a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 18417a6fd94SBarry Smith Mat A = Na->A; 185dfbe8321SBarry Smith PetscErrorCode ierr; 186b24ad042SBarry Smith PetscInt i,j,rstart,rend,nnz; 187b24ad042SBarry Smith const PetscInt *cols; 18817a6fd94SBarry Smith PetscScalar *diag,*work,*values; 189b3cc6726SBarry Smith const PetscScalar *mvalues; 19017a6fd94SBarry Smith 19117a6fd94SBarry Smith PetscFunctionBegin; 19274ed9c26SBarry Smith ierr = PetscMalloc2(A->cmap->N,PetscScalar,&diag,A->cmap->N,PetscScalar,&work);CHKERRQ(ierr); 193d0f46423SBarry Smith ierr = PetscMemzero(work,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 19417a6fd94SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 19517a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 196b3cc6726SBarry Smith ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 19717a6fd94SBarry Smith for (j=0; j<nnz; j++) { 198b3cc6726SBarry Smith work[cols[j]] += mvalues[j]*mvalues[j]; 19917a6fd94SBarry Smith } 200b3cc6726SBarry Smith ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 20117a6fd94SBarry Smith } 202a790c004SSatish Balay ierr = MPI_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,((PetscObject)N)->comm);CHKERRQ(ierr); 203d0f46423SBarry Smith rstart = N->cmap->rstart; 204d0f46423SBarry Smith rend = N->cmap->rend; 20517a6fd94SBarry Smith ierr = VecGetArray(v,&values);CHKERRQ(ierr); 20617a6fd94SBarry Smith ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr); 20717a6fd94SBarry Smith ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 20874ed9c26SBarry Smith ierr = PetscFree2(diag,work);CHKERRQ(ierr); 209b20f02adSBarry Smith ierr = VecScale(v,Na->scale);CHKERRQ(ierr); 21017a6fd94SBarry Smith PetscFunctionReturn(0); 21117a6fd94SBarry Smith } 212c8a8475eSBarry Smith 213c8a8475eSBarry Smith #undef __FUNCT__ 214c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal" 215c8a8475eSBarry Smith /*@ 216c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 217c8a8475eSBarry Smith 218c8a8475eSBarry Smith Collective on Mat 219c8a8475eSBarry Smith 220c8a8475eSBarry Smith Input Parameter: 221c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 222c8a8475eSBarry Smith 223c8a8475eSBarry Smith Output Parameter: 224c8a8475eSBarry Smith . N - the matrix that represents A'*A 225c8a8475eSBarry Smith 226c30e7cdbSSatish Balay Level: intermediate 227c30e7cdbSSatish Balay 228c8a8475eSBarry Smith Notes: The product A'*A is NOT actually formed! Rather the new matrix 229c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 230c8a8475eSBarry Smith A and then A' 231c8a8475eSBarry Smith @*/ 2327087cfbeSBarry Smith PetscErrorCode MatCreateNormal(Mat A,Mat *N) 233c8a8475eSBarry Smith { 234dfbe8321SBarry Smith PetscErrorCode ierr; 235b24ad042SBarry Smith PetscInt m,n; 236c8a8475eSBarry Smith Mat_Normal *Na; 237c8a8475eSBarry Smith 238c8a8475eSBarry Smith PetscFunctionBegin; 239c8a8475eSBarry Smith ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 2407adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,N);CHKERRQ(ierr); 241f69a0ea3SMatthew Knepley ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 242c8a8475eSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr); 243c8a8475eSBarry Smith 24438f2d2fdSLisandro Dalcin ierr = PetscNewLog(*N,Mat_Normal,&Na);CHKERRQ(ierr); 24538f2d2fdSLisandro Dalcin (*N)->data = (void*) Na; 246c8a8475eSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 247c3122656SLisandro Dalcin Na->A = A; 248b20f02adSBarry Smith Na->scale = 1.0; 249c8a8475eSBarry Smith 2507adad957SLisandro Dalcin ierr = VecCreateMPI(((PetscObject)A)->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr); 251*2205254eSKarl Rupp 252c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 253c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 254b20f02adSBarry Smith (*N)->ops->multtranspose = MatMultTranspose_Normal; 255b20f02adSBarry Smith (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal; 256db090513SMatthew Knepley (*N)->ops->multadd = MatMultAdd_Normal; 25717a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 25869ef1043SBarry Smith (*N)->ops->scale = MatScale_Normal; 25969ef1043SBarry Smith (*N)->ops->diagonalscale = MatDiagonalScale_Normal; 260c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 261d0f46423SBarry Smith (*N)->cmap->N = A->cmap->N; 262d0f46423SBarry Smith (*N)->rmap->N = A->cmap->N; 263d0f46423SBarry Smith (*N)->cmap->n = A->cmap->n; 264d0f46423SBarry Smith (*N)->rmap->n = A->cmap->n; 265c8a8475eSBarry Smith PetscFunctionReturn(0); 266c8a8475eSBarry Smith } 267c8a8475eSBarry Smith 268