1c8a8475eSBarry Smith /*$Id: bvec2.c,v 1.202 2001/09/12 03:26:24 bsmith Exp $*/ 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" 12c8a8475eSBarry Smith int MatMult_Normal(Mat N,Vec x,Vec y) 13c8a8475eSBarry Smith { 14c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 15c8a8475eSBarry Smith int 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" 25c8a8475eSBarry Smith int MatDestroy_Normal(Mat N) 26c8a8475eSBarry Smith { 27c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 28c8a8475eSBarry Smith int 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 37*17a6fd94SBarry Smith /* 38*17a6fd94SBarry Smith Slow, nonscalable version 39*17a6fd94SBarry Smith */ 40*17a6fd94SBarry Smith #undef __FUNCT__ 41*17a6fd94SBarry Smith #define __FUNCT__ "MatGetDiagonal_Normal" 42*17a6fd94SBarry Smith int MatGetDiagonal_Normal(Mat N,Vec v) 43*17a6fd94SBarry Smith { 44*17a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 45*17a6fd94SBarry Smith Mat A = Na->A; 46*17a6fd94SBarry Smith int ierr,i,j,rstart,rend,nnz,*cols; 47*17a6fd94SBarry Smith PetscScalar *diag,*work,*values; 48*17a6fd94SBarry Smith PetscMap cmap; 49*17a6fd94SBarry Smith 50*17a6fd94SBarry Smith PetscFunctionBegin; 51*17a6fd94SBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr); 52*17a6fd94SBarry Smith work = diag + A->N; 53*17a6fd94SBarry Smith ierr = PetscMemzero(work,A->N*sizeof(PetscScalar));CHKERRQ(ierr); 54*17a6fd94SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 55*17a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 56*17a6fd94SBarry Smith ierr = MatGetRow(A,i,&nnz,&cols,&values);CHKERRQ(ierr); 57*17a6fd94SBarry Smith for (j=0; j<nnz; j++) { 58*17a6fd94SBarry Smith work[cols[j]] += values[j]*values[j]; 59*17a6fd94SBarry Smith } 60*17a6fd94SBarry Smith ierr = MatRestoreRow(A,i,&nnz,&cols,&values);CHKERRQ(ierr); 61*17a6fd94SBarry Smith } 62*17a6fd94SBarry Smith ierr = MPI_Allreduce(work,diag,A->N,MPIU_SCALAR,MPI_SUM,N->comm);CHKERRQ(ierr); 63*17a6fd94SBarry Smith ierr = MatGetPetscMaps(A,PETSC_NULL,&cmap);CHKERRQ(ierr); 64*17a6fd94SBarry Smith ierr = PetscMapGetLocalRange(cmap,&rstart,&rend);CHKERRQ(ierr); 65*17a6fd94SBarry Smith ierr = VecGetArray(v,&values);CHKERRQ(ierr); 66*17a6fd94SBarry Smith ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr); 67*17a6fd94SBarry Smith ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 68*17a6fd94SBarry Smith ierr = PetscFree(diag);CHKERRQ(ierr); 69*17a6fd94SBarry Smith PetscFunctionReturn(0); 70*17a6fd94SBarry Smith } 71c8a8475eSBarry Smith 72c8a8475eSBarry Smith #undef __FUNCT__ 73c8a8475eSBarry Smith #define __FUNCT__ "MatCreateNormal" 74c8a8475eSBarry Smith /*@ 75c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 76c8a8475eSBarry Smith 77c8a8475eSBarry Smith Collective on Mat 78c8a8475eSBarry Smith 79c8a8475eSBarry Smith Input Parameter: 80c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 81c8a8475eSBarry Smith 82c8a8475eSBarry Smith Output Parameter: 83c8a8475eSBarry Smith . N - the matrix that represents A'*A 84c8a8475eSBarry Smith 85c8a8475eSBarry Smith Notes: The product A'*A is NOT actually formed! Rather the new matrix 86c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 87c8a8475eSBarry Smith A and then A' 88c8a8475eSBarry Smith @*/ 89c8a8475eSBarry Smith int MatCreateNormal(Mat A,Mat *N) 90c8a8475eSBarry Smith { 91c8a8475eSBarry Smith int ierr,m,n; 92c8a8475eSBarry Smith Mat_Normal *Na; 93c8a8475eSBarry Smith 94c8a8475eSBarry Smith PetscFunctionBegin; 95c8a8475eSBarry Smith ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 96c8a8475eSBarry Smith ierr = MatCreate(A->comm,n,n,PETSC_DECIDE,PETSC_DECIDE,N);CHKERRQ(ierr); 97c8a8475eSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr); 98c8a8475eSBarry Smith 99c8a8475eSBarry Smith ierr = PetscNew(Mat_Normal,&Na);CHKERRQ(ierr); 100c8a8475eSBarry Smith Na->A = A; 101c8a8475eSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 102c8a8475eSBarry Smith (*N)->data = (void*) Na; 103c8a8475eSBarry Smith 104c8a8475eSBarry Smith ierr = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr); 105c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 106c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 107*17a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 108c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 109c8a8475eSBarry Smith (*N)->N = A->N; 110c8a8475eSBarry Smith (*N)->M = A->N; 111c8a8475eSBarry Smith (*N)->n = A->n; 112c8a8475eSBarry Smith (*N)->m = A->n; 113c8a8475eSBarry Smith PetscFunctionReturn(0); 114c8a8475eSBarry Smith } 115c8a8475eSBarry Smith 116