1c8a8475eSBarry Smith 2af0996ceSBarry Smith #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 1069ef1043SBarry Smith PetscErrorCode MatScale_Normal(Mat inA,PetscScalar scale) 1169ef1043SBarry Smith { 1269ef1043SBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 1369ef1043SBarry Smith 1469ef1043SBarry Smith PetscFunctionBegin; 1569ef1043SBarry Smith a->scale *= scale; 1669ef1043SBarry Smith PetscFunctionReturn(0); 1769ef1043SBarry Smith } 1869ef1043SBarry Smith 197cfd0b8cSBarry Smith PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right) 207cfd0b8cSBarry Smith { 217cfd0b8cSBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 227cfd0b8cSBarry Smith PetscErrorCode ierr; 237cfd0b8cSBarry Smith 247cfd0b8cSBarry Smith PetscFunctionBegin; 257cfd0b8cSBarry Smith if (left) { 267cfd0b8cSBarry Smith if (!a->left) { 277cfd0b8cSBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 287cfd0b8cSBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 297cfd0b8cSBarry Smith } else { 307cfd0b8cSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 317cfd0b8cSBarry Smith } 327cfd0b8cSBarry Smith } 337cfd0b8cSBarry Smith if (right) { 347cfd0b8cSBarry Smith if (!a->right) { 357cfd0b8cSBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 367cfd0b8cSBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 377cfd0b8cSBarry Smith } else { 387cfd0b8cSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 397cfd0b8cSBarry Smith } 407cfd0b8cSBarry Smith } 417cfd0b8cSBarry Smith PetscFunctionReturn(0); 427cfd0b8cSBarry Smith } 437cfd0b8cSBarry Smith 44dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y) 45c8a8475eSBarry Smith { 46c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 47dfbe8321SBarry Smith PetscErrorCode ierr; 487cfd0b8cSBarry Smith Vec in; 49c8a8475eSBarry Smith 50c8a8475eSBarry Smith PetscFunctionBegin; 517cfd0b8cSBarry Smith in = x; 527cfd0b8cSBarry Smith if (Na->right) { 537cfd0b8cSBarry Smith if (!Na->rightwork) { 547cfd0b8cSBarry Smith ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 557cfd0b8cSBarry Smith } 567cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 577cfd0b8cSBarry Smith in = Na->rightwork; 587cfd0b8cSBarry Smith } 597cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 60c8a8475eSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 617cfd0b8cSBarry Smith if (Na->left) { 627cfd0b8cSBarry Smith ierr = VecPointwiseMult(y,Na->left,y);CHKERRQ(ierr); 637cfd0b8cSBarry Smith } 64b20f02adSBarry Smith ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 65c8a8475eSBarry Smith PetscFunctionReturn(0); 66c8a8475eSBarry Smith } 67c8a8475eSBarry Smith 68db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 69db090513SMatthew Knepley { 70db090513SMatthew Knepley Mat_Normal *Na = (Mat_Normal*)N->data; 71db090513SMatthew Knepley PetscErrorCode ierr; 727cfd0b8cSBarry Smith Vec in; 73db090513SMatthew Knepley 74db090513SMatthew Knepley PetscFunctionBegin; 757cfd0b8cSBarry Smith in = v1; 767cfd0b8cSBarry Smith if (Na->right) { 777cfd0b8cSBarry Smith if (!Na->rightwork) { 787cfd0b8cSBarry Smith ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 797cfd0b8cSBarry Smith } 807cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 817cfd0b8cSBarry Smith in = Na->rightwork; 827cfd0b8cSBarry Smith } 837cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 84b20f02adSBarry Smith ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 857cfd0b8cSBarry Smith if (Na->left) { 867cfd0b8cSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 877cfd0b8cSBarry Smith ierr = VecPointwiseMult(v3,Na->left,v3);CHKERRQ(ierr); 887cfd0b8cSBarry Smith ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 897cfd0b8cSBarry Smith } else { 90b20f02adSBarry Smith ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 917cfd0b8cSBarry Smith } 92b20f02adSBarry Smith PetscFunctionReturn(0); 93b20f02adSBarry Smith } 94b20f02adSBarry Smith 95b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y) 96b20f02adSBarry Smith { 97b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 98b20f02adSBarry Smith PetscErrorCode ierr; 997cfd0b8cSBarry Smith Vec in; 100b20f02adSBarry Smith 101b20f02adSBarry Smith PetscFunctionBegin; 1027cfd0b8cSBarry Smith in = x; 1037cfd0b8cSBarry Smith if (Na->left) { 1047cfd0b8cSBarry Smith if (!Na->leftwork) { 1057cfd0b8cSBarry Smith ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 1067cfd0b8cSBarry Smith } 1077cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 1087cfd0b8cSBarry Smith in = Na->leftwork; 1097cfd0b8cSBarry Smith } 1107cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 111b20f02adSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 1127cfd0b8cSBarry Smith if (Na->right) { 1137cfd0b8cSBarry Smith ierr = VecPointwiseMult(y,Na->right,y);CHKERRQ(ierr); 1147cfd0b8cSBarry Smith } 115b20f02adSBarry Smith ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 116b20f02adSBarry Smith PetscFunctionReturn(0); 117b20f02adSBarry Smith } 118b20f02adSBarry Smith 119b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 120b20f02adSBarry Smith { 121b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 122b20f02adSBarry Smith PetscErrorCode ierr; 1237cfd0b8cSBarry Smith Vec in; 124b20f02adSBarry Smith 125b20f02adSBarry Smith PetscFunctionBegin; 1267cfd0b8cSBarry Smith in = v1; 1277cfd0b8cSBarry Smith if (Na->left) { 1287cfd0b8cSBarry Smith if (!Na->leftwork) { 1297cfd0b8cSBarry Smith ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 1307cfd0b8cSBarry Smith } 1317cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 1327cfd0b8cSBarry Smith in = Na->leftwork; 1337cfd0b8cSBarry Smith } 1347cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 135b20f02adSBarry Smith ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 1367cfd0b8cSBarry Smith if (Na->right) { 1377cfd0b8cSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 1387cfd0b8cSBarry Smith ierr = VecPointwiseMult(v3,Na->right,v3);CHKERRQ(ierr); 1397cfd0b8cSBarry Smith ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1407cfd0b8cSBarry Smith } else { 141db090513SMatthew Knepley ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 1427cfd0b8cSBarry Smith } 143db090513SMatthew Knepley PetscFunctionReturn(0); 144db090513SMatthew Knepley } 145db090513SMatthew Knepley 146dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N) 147c8a8475eSBarry Smith { 148c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 149dfbe8321SBarry Smith PetscErrorCode ierr; 150c8a8475eSBarry Smith 151c8a8475eSBarry Smith PetscFunctionBegin; 1526bf464f9SBarry Smith ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 1536bf464f9SBarry Smith ierr = VecDestroy(&Na->w);CHKERRQ(ierr); 1546bf464f9SBarry Smith ierr = VecDestroy(&Na->left);CHKERRQ(ierr); 1556bf464f9SBarry Smith ierr = VecDestroy(&Na->right);CHKERRQ(ierr); 1566bf464f9SBarry Smith ierr = VecDestroy(&Na->leftwork);CHKERRQ(ierr); 1576bf464f9SBarry Smith ierr = VecDestroy(&Na->rightwork);CHKERRQ(ierr); 158bf0cc555SLisandro Dalcin ierr = PetscFree(N->data);CHKERRQ(ierr); 159c8a8475eSBarry Smith PetscFunctionReturn(0); 160c8a8475eSBarry Smith } 161c8a8475eSBarry Smith 16217a6fd94SBarry Smith /* 16317a6fd94SBarry Smith Slow, nonscalable version 16417a6fd94SBarry Smith */ 165dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 16617a6fd94SBarry Smith { 16717a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 16817a6fd94SBarry Smith Mat A = Na->A; 169dfbe8321SBarry Smith PetscErrorCode ierr; 170b24ad042SBarry Smith PetscInt i,j,rstart,rend,nnz; 171b24ad042SBarry Smith const PetscInt *cols; 17217a6fd94SBarry Smith PetscScalar *diag,*work,*values; 173b3cc6726SBarry Smith const PetscScalar *mvalues; 17417a6fd94SBarry Smith 17517a6fd94SBarry Smith PetscFunctionBegin; 176dcca6d9dSJed Brown ierr = PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work);CHKERRQ(ierr); 177580bdb30SBarry Smith ierr = PetscArrayzero(work,A->cmap->N);CHKERRQ(ierr); 17817a6fd94SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 17917a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 180b3cc6726SBarry Smith ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 18117a6fd94SBarry Smith for (j=0; j<nnz; j++) { 182b3cc6726SBarry Smith work[cols[j]] += mvalues[j]*mvalues[j]; 18317a6fd94SBarry Smith } 184b3cc6726SBarry Smith ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 18517a6fd94SBarry Smith } 186b2566f29SBarry Smith ierr = MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr); 187d0f46423SBarry Smith rstart = N->cmap->rstart; 188d0f46423SBarry Smith rend = N->cmap->rend; 18917a6fd94SBarry Smith ierr = VecGetArray(v,&values);CHKERRQ(ierr); 190580bdb30SBarry Smith ierr = PetscArraycpy(values,diag+rstart,rend-rstart);CHKERRQ(ierr); 19117a6fd94SBarry Smith ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 19274ed9c26SBarry Smith ierr = PetscFree2(diag,work);CHKERRQ(ierr); 193b20f02adSBarry Smith ierr = VecScale(v,Na->scale);CHKERRQ(ierr); 19417a6fd94SBarry Smith PetscFunctionReturn(0); 19517a6fd94SBarry Smith } 196c8a8475eSBarry Smith 197c8a8475eSBarry Smith /*@ 198c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 199c8a8475eSBarry Smith 200c8a8475eSBarry Smith Collective on Mat 201c8a8475eSBarry Smith 202c8a8475eSBarry Smith Input Parameter: 203c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 204c8a8475eSBarry Smith 205c8a8475eSBarry Smith Output Parameter: 206c8a8475eSBarry Smith . N - the matrix that represents A'*A 207c8a8475eSBarry Smith 208c30e7cdbSSatish Balay Level: intermediate 209c30e7cdbSSatish Balay 21095452b02SPatrick Sanan Notes: 21195452b02SPatrick Sanan The product A'*A is NOT actually formed! Rather the new matrix 212c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 213c8a8475eSBarry Smith A and then A' 214c8a8475eSBarry Smith @*/ 2157087cfbeSBarry Smith PetscErrorCode MatCreateNormal(Mat A,Mat *N) 216c8a8475eSBarry Smith { 217dfbe8321SBarry Smith PetscErrorCode ierr; 218*9ca0e1b6SStefano Zampini PetscInt n,nn; 219c8a8475eSBarry Smith Mat_Normal *Na; 220*9ca0e1b6SStefano Zampini VecType vtype; 221c8a8475eSBarry Smith 222c8a8475eSBarry Smith PetscFunctionBegin; 223*9ca0e1b6SStefano Zampini ierr = MatGetSize(A,NULL,&nn);CHKERRQ(ierr); 224*9ca0e1b6SStefano Zampini ierr = MatGetLocalSize(A,NULL,&n);CHKERRQ(ierr); 225ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),N);CHKERRQ(ierr); 226*9ca0e1b6SStefano Zampini ierr = MatSetSizes(*N,n,n,nn,nn);CHKERRQ(ierr); 227c8a8475eSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr); 228c8a8475eSBarry Smith 229b00a9115SJed Brown ierr = PetscNewLog(*N,&Na);CHKERRQ(ierr); 23038f2d2fdSLisandro Dalcin (*N)->data = (void*) Na; 231c8a8475eSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 232c3122656SLisandro Dalcin Na->A = A; 233b20f02adSBarry Smith Na->scale = 1.0; 234c8a8475eSBarry Smith 235*9ca0e1b6SStefano Zampini ierr = MatCreateVecs(A,NULL,&Na->w);CHKERRQ(ierr); 2362205254eSKarl Rupp 237c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 238c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 239b20f02adSBarry Smith (*N)->ops->multtranspose = MatMultTranspose_Normal; 240b20f02adSBarry Smith (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal; 241db090513SMatthew Knepley (*N)->ops->multadd = MatMultAdd_Normal; 24217a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 24369ef1043SBarry Smith (*N)->ops->scale = MatScale_Normal; 24469ef1043SBarry Smith (*N)->ops->diagonalscale = MatDiagonalScale_Normal; 245c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 24648331cefSBarry Smith (*N)->preallocated = PETSC_TRUE; 247*9ca0e1b6SStefano Zampini 248*9ca0e1b6SStefano Zampini ierr = MatGetVecType(A,&vtype);CHKERRQ(ierr); 249*9ca0e1b6SStefano Zampini ierr = MatSetVecType(*N,vtype);CHKERRQ(ierr); 250*9ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 251*9ca0e1b6SStefano Zampini ierr = MatBindToCPU(*N,A->boundtocpu);CHKERRQ(ierr); 252*9ca0e1b6SStefano Zampini #endif 253c8a8475eSBarry Smith PetscFunctionReturn(0); 254c8a8475eSBarry Smith } 255