1be1d678aSKris Buschelman #define PETSCMAT_DLL 2c8a8475eSBarry Smith 3b9147fbbSdalcinl #include "include/private/matimpl.h" /*I "petscmat.h" I*/ 4c8a8475eSBarry Smith 5c8a8475eSBarry Smith typedef struct { 6c8a8475eSBarry Smith Mat A; 7c8a8475eSBarry Smith Vec w; 8c8a8475eSBarry Smith } Mat_Normal; 9c8a8475eSBarry Smith 10c8a8475eSBarry Smith #undef __FUNCT__ 11c8a8475eSBarry Smith #define __FUNCT__ "MatMult_Normal" 12dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y) 13c8a8475eSBarry Smith { 14c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 15dfbe8321SBarry Smith PetscErrorCode ierr; 16c8a8475eSBarry Smith 17c8a8475eSBarry Smith PetscFunctionBegin; 18c8a8475eSBarry Smith ierr = MatMult(Na->A,x,Na->w);CHKERRQ(ierr); 19c8a8475eSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 20c8a8475eSBarry Smith PetscFunctionReturn(0); 21c8a8475eSBarry Smith } 22c8a8475eSBarry Smith 23c8a8475eSBarry Smith #undef __FUNCT__ 24db090513SMatthew Knepley #define __FUNCT__ "MatMultAdd_Normal" 25db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 26db090513SMatthew Knepley { 27db090513SMatthew Knepley Mat_Normal *Na = (Mat_Normal*)N->data; 28db090513SMatthew Knepley PetscErrorCode ierr; 29db090513SMatthew Knepley 30db090513SMatthew Knepley PetscFunctionBegin; 31db090513SMatthew Knepley ierr = MatMult(Na->A,v1,Na->w);CHKERRQ(ierr); 32db090513SMatthew Knepley ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 33db090513SMatthew Knepley PetscFunctionReturn(0); 34db090513SMatthew Knepley } 35db090513SMatthew Knepley 36db090513SMatthew Knepley #undef __FUNCT__ 37c8a8475eSBarry Smith #define __FUNCT__ "MatDestroy_Normal" 38dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N) 39c8a8475eSBarry Smith { 40c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 41dfbe8321SBarry Smith PetscErrorCode ierr; 42c8a8475eSBarry Smith 43c8a8475eSBarry Smith PetscFunctionBegin; 44c3122656SLisandro Dalcin if (Na->A) { ierr = MatDestroy(Na->A);CHKERRQ(ierr); } 45c3122656SLisandro Dalcin if (Na->w) { ierr = VecDestroy(Na->w);CHKERRQ(ierr); } 46c8a8475eSBarry Smith ierr = PetscFree(Na);CHKERRQ(ierr); 47c8a8475eSBarry Smith PetscFunctionReturn(0); 48c8a8475eSBarry Smith } 49c8a8475eSBarry Smith 5017a6fd94SBarry Smith /* 5117a6fd94SBarry Smith Slow, nonscalable version 5217a6fd94SBarry Smith */ 5317a6fd94SBarry Smith #undef __FUNCT__ 5417a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal" 55dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 5617a6fd94SBarry Smith { 5717a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 5817a6fd94SBarry Smith Mat A = Na->A; 59dfbe8321SBarry Smith PetscErrorCode ierr; 60b24ad042SBarry Smith PetscInt i,j,rstart,rend,nnz; 61b24ad042SBarry Smith const PetscInt *cols; 6217a6fd94SBarry Smith PetscScalar *diag,*work,*values; 63b3cc6726SBarry Smith const PetscScalar *mvalues; 6417a6fd94SBarry Smith 6517a6fd94SBarry Smith PetscFunctionBegin; 66899cda47SBarry Smith ierr = PetscMalloc(2*A->cmap.N*sizeof(PetscScalar),&diag);CHKERRQ(ierr); 67899cda47SBarry Smith work = diag + A->cmap.N; 68899cda47SBarry Smith ierr = PetscMemzero(work,A->cmap.N*sizeof(PetscScalar));CHKERRQ(ierr); 6917a6fd94SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 7017a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 71b3cc6726SBarry Smith ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 7217a6fd94SBarry Smith for (j=0; j<nnz; j++) { 73b3cc6726SBarry Smith work[cols[j]] += mvalues[j]*mvalues[j]; 7417a6fd94SBarry Smith } 75b3cc6726SBarry Smith ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 7617a6fd94SBarry Smith } 77899cda47SBarry Smith ierr = MPI_Allreduce(work,diag,A->cmap.N,MPIU_SCALAR,MPI_SUM,N->comm);CHKERRQ(ierr); 78357abbc8SBarry Smith rstart = N->cmap.rstart; 79357abbc8SBarry Smith rend = N->cmap.rend; 8017a6fd94SBarry Smith ierr = VecGetArray(v,&values);CHKERRQ(ierr); 8117a6fd94SBarry Smith ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr); 8217a6fd94SBarry Smith ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 8317a6fd94SBarry Smith ierr = PetscFree(diag);CHKERRQ(ierr); 8417a6fd94SBarry Smith PetscFunctionReturn(0); 8517a6fd94SBarry Smith } 86c8a8475eSBarry Smith 87c8a8475eSBarry Smith #undef __FUNCT__ 88c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal" 89c8a8475eSBarry Smith /*@ 90c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 91c8a8475eSBarry Smith 92c8a8475eSBarry Smith Collective on Mat 93c8a8475eSBarry Smith 94c8a8475eSBarry Smith Input Parameter: 95c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 96c8a8475eSBarry Smith 97c8a8475eSBarry Smith Output Parameter: 98c8a8475eSBarry Smith . N - the matrix that represents A'*A 99c8a8475eSBarry Smith 100c30e7cdbSSatish Balay Level: intermediate 101c30e7cdbSSatish Balay 102c8a8475eSBarry Smith Notes: The product A'*A is NOT actually formed! Rather the new matrix 103c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 104c8a8475eSBarry Smith A and then A' 105c8a8475eSBarry Smith @*/ 106be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat A,Mat *N) 107c8a8475eSBarry Smith { 108dfbe8321SBarry Smith PetscErrorCode ierr; 109b24ad042SBarry Smith PetscInt m,n; 110c8a8475eSBarry Smith Mat_Normal *Na; 111c8a8475eSBarry Smith 112c8a8475eSBarry Smith PetscFunctionBegin; 113c8a8475eSBarry Smith ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 114f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,N);CHKERRQ(ierr); 115f69a0ea3SMatthew Knepley ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 116c8a8475eSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr); 117c8a8475eSBarry Smith 118*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(*N,Mat_Normal,&Na);CHKERRQ(ierr); 119*38f2d2fdSLisandro Dalcin (*N)->data = (void*) Na; 120c8a8475eSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 121c3122656SLisandro Dalcin Na->A = A; 122c8a8475eSBarry Smith 123c8a8475eSBarry Smith ierr = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr); 124c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 125c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 126db090513SMatthew Knepley (*N)->ops->multadd = MatMultAdd_Normal; 12717a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 128c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 129899cda47SBarry Smith (*N)->cmap.N = A->cmap.N; 130899cda47SBarry Smith (*N)->rmap.N = A->cmap.N; 131899cda47SBarry Smith (*N)->cmap.n = A->cmap.n; 132899cda47SBarry Smith (*N)->rmap.n = A->cmap.n; 133c8a8475eSBarry Smith PetscFunctionReturn(0); 134c8a8475eSBarry Smith } 135c8a8475eSBarry Smith 136