1*be1d678aSKris Buschelman #define PETSCMAT_DLL 2c8a8475eSBarry Smith 3c8a8475eSBarry Smith #include "src/mat/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__ 24c8a8475eSBarry Smith #define __FUNCT__ "MatDestroy_Normal" 25dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N) 26c8a8475eSBarry Smith { 27c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 28dfbe8321SBarry Smith PetscErrorCode ierr; 29c8a8475eSBarry Smith 30c8a8475eSBarry Smith PetscFunctionBegin; 31c8a8475eSBarry Smith ierr = PetscObjectDereference((PetscObject)Na->A);CHKERRQ(ierr); 32c8a8475eSBarry Smith ierr = VecDestroy(Na->w);CHKERRQ(ierr); 33c8a8475eSBarry Smith ierr = PetscFree(Na);CHKERRQ(ierr); 34c8a8475eSBarry Smith PetscFunctionReturn(0); 35c8a8475eSBarry Smith } 36c8a8475eSBarry Smith 3717a6fd94SBarry Smith /* 3817a6fd94SBarry Smith Slow, nonscalable version 3917a6fd94SBarry Smith */ 4017a6fd94SBarry Smith #undef __FUNCT__ 4117a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal" 42dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 4317a6fd94SBarry Smith { 4417a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 4517a6fd94SBarry Smith Mat A = Na->A; 46dfbe8321SBarry Smith PetscErrorCode ierr; 47b24ad042SBarry Smith PetscInt i,j,rstart,rend,nnz; 48b24ad042SBarry Smith const PetscInt *cols; 4917a6fd94SBarry Smith PetscScalar *diag,*work,*values; 50b3cc6726SBarry Smith const PetscScalar *mvalues; 5117a6fd94SBarry Smith PetscMap cmap; 5217a6fd94SBarry Smith 5317a6fd94SBarry Smith PetscFunctionBegin; 5417a6fd94SBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr); 5517a6fd94SBarry Smith work = diag + A->N; 5617a6fd94SBarry Smith ierr = PetscMemzero(work,A->N*sizeof(PetscScalar));CHKERRQ(ierr); 5717a6fd94SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 5817a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 59b3cc6726SBarry Smith ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 6017a6fd94SBarry Smith for (j=0; j<nnz; j++) { 61b3cc6726SBarry Smith work[cols[j]] += mvalues[j]*mvalues[j]; 6217a6fd94SBarry Smith } 63b3cc6726SBarry Smith ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 6417a6fd94SBarry Smith } 6517a6fd94SBarry Smith ierr = MPI_Allreduce(work,diag,A->N,MPIU_SCALAR,MPI_SUM,N->comm);CHKERRQ(ierr); 6617a6fd94SBarry Smith ierr = MatGetPetscMaps(A,PETSC_NULL,&cmap);CHKERRQ(ierr); 6717a6fd94SBarry Smith ierr = PetscMapGetLocalRange(cmap,&rstart,&rend);CHKERRQ(ierr); 6817a6fd94SBarry Smith ierr = VecGetArray(v,&values);CHKERRQ(ierr); 6917a6fd94SBarry Smith ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr); 7017a6fd94SBarry Smith ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 7117a6fd94SBarry Smith ierr = PetscFree(diag);CHKERRQ(ierr); 7217a6fd94SBarry Smith PetscFunctionReturn(0); 7317a6fd94SBarry Smith } 74c8a8475eSBarry Smith 75c8a8475eSBarry Smith #undef __FUNCT__ 76c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal" 77c8a8475eSBarry Smith /*@ 78c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 79c8a8475eSBarry Smith 80c8a8475eSBarry Smith Collective on Mat 81c8a8475eSBarry Smith 82c8a8475eSBarry Smith Input Parameter: 83c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 84c8a8475eSBarry Smith 85c8a8475eSBarry Smith Output Parameter: 86c8a8475eSBarry Smith . N - the matrix that represents A'*A 87c8a8475eSBarry Smith 88c30e7cdbSSatish Balay Level: intermediate 89c30e7cdbSSatish Balay 90c8a8475eSBarry Smith Notes: The product A'*A is NOT actually formed! Rather the new matrix 91c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 92c8a8475eSBarry Smith A and then A' 93c8a8475eSBarry Smith @*/ 94*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat A,Mat *N) 95c8a8475eSBarry Smith { 96dfbe8321SBarry Smith PetscErrorCode ierr; 97b24ad042SBarry Smith PetscInt m,n; 98c8a8475eSBarry Smith Mat_Normal *Na; 99c8a8475eSBarry Smith 100c8a8475eSBarry Smith PetscFunctionBegin; 101c8a8475eSBarry Smith ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 102c8a8475eSBarry Smith ierr = MatCreate(A->comm,n,n,PETSC_DECIDE,PETSC_DECIDE,N);CHKERRQ(ierr); 103c8a8475eSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr); 104c8a8475eSBarry Smith 105c8a8475eSBarry Smith ierr = PetscNew(Mat_Normal,&Na);CHKERRQ(ierr); 106c8a8475eSBarry Smith Na->A = A; 107c8a8475eSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 108c8a8475eSBarry Smith (*N)->data = (void*) Na; 109c8a8475eSBarry Smith 110c8a8475eSBarry Smith ierr = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr); 111c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 112c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 11317a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 114c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 115c8a8475eSBarry Smith (*N)->N = A->N; 116c8a8475eSBarry Smith (*N)->M = A->N; 117c8a8475eSBarry Smith (*N)->n = A->n; 118c8a8475eSBarry Smith (*N)->m = A->n; 119c8a8475eSBarry Smith PetscFunctionReturn(0); 120c8a8475eSBarry Smith } 121c8a8475eSBarry Smith 122