1ceb9aaf7SBarry Smith 2b45d2f2cSJed Brown #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/ 3c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> 4ceb9aaf7SBarry Smith 5ceb9aaf7SBarry Smith typedef struct { 6ab92ecdeSBarry Smith Mat A,U,V; 7ab92ecdeSBarry Smith Vec work1,work2; /* Sequential (big) vectors that hold partial products */ 8ab92ecdeSBarry Smith PetscMPIInt nwork; /* length of work vectors */ 9ab92ecdeSBarry Smith } Mat_LRC; 10ab92ecdeSBarry Smith 11ab92ecdeSBarry Smith 12ceb9aaf7SBarry Smith 13ceb9aaf7SBarry Smith #undef __FUNCT__ 14ab92ecdeSBarry Smith #define __FUNCT__ "MatMult_LRC" 15ab92ecdeSBarry Smith PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y) 16ceb9aaf7SBarry Smith { 17ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 18ceb9aaf7SBarry Smith PetscErrorCode ierr; 19ab92ecdeSBarry Smith PetscScalar *w1,*w2; 20ceb9aaf7SBarry Smith 21ceb9aaf7SBarry Smith PetscFunctionBegin; 22ab92ecdeSBarry Smith ierr = MatMult(Na->A,x,y);CHKERRQ(ierr); 23ab92ecdeSBarry Smith 24ab92ecdeSBarry Smith /* multiply the local part of V with the local part of x */ 25ab92ecdeSBarry Smith /* note in this call x is treated as a sequential vector */ 26ab92ecdeSBarry Smith ierr = MatMultTranspose_SeqDense(Na->V,x,Na->work1);CHKERRQ(ierr); 27ab92ecdeSBarry Smith 28ab92ecdeSBarry Smith /* Form the sum of all the local multiplies : this is work2 = V'*x = 29ab92ecdeSBarry Smith sum_{all processors} work1 */ 30ab92ecdeSBarry Smith 31ab92ecdeSBarry Smith ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr); 32ab92ecdeSBarry Smith ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr); 33*ce94432eSBarry Smith ierr = MPI_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr); 34ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr); 35ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr); 36ab92ecdeSBarry Smith 37ab92ecdeSBarry Smith /* multiply-sub y = y + U*work2 */ 38ab92ecdeSBarry Smith /* note in this call y is treated as a sequential vector */ 39ab92ecdeSBarry Smith ierr = MatMultAdd_SeqDense(Na->U,Na->work2,y,y);CHKERRQ(ierr); 40ceb9aaf7SBarry Smith PetscFunctionReturn(0); 41ceb9aaf7SBarry Smith } 42ceb9aaf7SBarry Smith 43ceb9aaf7SBarry Smith #undef __FUNCT__ 44ab92ecdeSBarry Smith #define __FUNCT__ "MatDestroy_LRC" 45ab92ecdeSBarry Smith PetscErrorCode MatDestroy_LRC(Mat N) 46ceb9aaf7SBarry Smith { 47ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 48ceb9aaf7SBarry Smith PetscErrorCode ierr; 49ceb9aaf7SBarry Smith 50ceb9aaf7SBarry Smith PetscFunctionBegin; 51bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 52bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->U);CHKERRQ(ierr); 53bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->V);CHKERRQ(ierr); 546bf464f9SBarry Smith ierr = VecDestroy(&Na->work1);CHKERRQ(ierr); 556bf464f9SBarry Smith ierr = VecDestroy(&Na->work2);CHKERRQ(ierr); 56bf0cc555SLisandro Dalcin ierr = PetscFree(N->data);CHKERRQ(ierr); 57ceb9aaf7SBarry Smith PetscFunctionReturn(0); 58ceb9aaf7SBarry Smith } 59ceb9aaf7SBarry Smith 60ceb9aaf7SBarry Smith 61ceb9aaf7SBarry Smith #undef __FUNCT__ 62ab92ecdeSBarry Smith #define __FUNCT__ "MatCreateLRC" 63ceb9aaf7SBarry Smith /*@ 64ab92ecdeSBarry Smith MatCreateLRC - Creates a new matrix object that behaves like A + U*V' 65ceb9aaf7SBarry Smith 66ceb9aaf7SBarry Smith Collective on Mat 67ceb9aaf7SBarry Smith 68ceb9aaf7SBarry Smith Input Parameter: 69ab92ecdeSBarry Smith + A - the (sparse) matrix 70ab92ecdeSBarry Smith - U. V - two dense rectangular (tall and skinny) matrices 71ceb9aaf7SBarry Smith 72ceb9aaf7SBarry Smith Output Parameter: 73ab92ecdeSBarry Smith . N - the matrix that represents A + U*V' 74ceb9aaf7SBarry Smith 75ceb9aaf7SBarry Smith Level: intermediate 76ceb9aaf7SBarry Smith 773f32ffd6SMatthew Knepley Notes: The matrix A + U*V' is not formed! Rather the new matrix 78ceb9aaf7SBarry Smith object performs the matrix-vector product by first multiplying by 79ab92ecdeSBarry Smith A and then adding the other term 80ceb9aaf7SBarry Smith @*/ 817087cfbeSBarry Smith PetscErrorCode MatCreateLRC(Mat A,Mat U, Mat V,Mat *N) 82ceb9aaf7SBarry Smith { 83ceb9aaf7SBarry Smith PetscErrorCode ierr; 84ceb9aaf7SBarry Smith PetscInt m,n; 85ab92ecdeSBarry Smith Mat_LRC *Na; 86ceb9aaf7SBarry Smith 87ceb9aaf7SBarry Smith PetscFunctionBegin; 88ceb9aaf7SBarry Smith ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 89*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),N);CHKERRQ(ierr); 90ceb9aaf7SBarry Smith ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 91ab92ecdeSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr); 92ceb9aaf7SBarry Smith 9338f2d2fdSLisandro Dalcin ierr = PetscNewLog(*N,Mat_LRC,&Na);CHKERRQ(ierr); 9438f2d2fdSLisandro Dalcin (*N)->data = (void*) Na; 95ceb9aaf7SBarry Smith Na->A = A; 96ab92ecdeSBarry Smith 97ab92ecdeSBarry Smith ierr = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr); 98ab92ecdeSBarry Smith ierr = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr); 99ceb9aaf7SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 100ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr); 101ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr); 102ceb9aaf7SBarry Smith 103d0f46423SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr); 104ab92ecdeSBarry Smith ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr); 105d0f46423SBarry Smith Na->nwork = U->cmap->N; 106ab92ecdeSBarry Smith 107ab92ecdeSBarry Smith (*N)->ops->destroy = MatDestroy_LRC; 108ab92ecdeSBarry Smith (*N)->ops->mult = MatMult_LRC; 109ceb9aaf7SBarry Smith (*N)->assembled = PETSC_TRUE; 110d0f46423SBarry Smith (*N)->cmap->N = A->cmap->N; 111d0f46423SBarry Smith (*N)->rmap->N = A->cmap->N; 112d0f46423SBarry Smith (*N)->cmap->n = A->cmap->n; 113d0f46423SBarry Smith (*N)->rmap->n = A->cmap->n; 114ceb9aaf7SBarry Smith PetscFunctionReturn(0); 115ceb9aaf7SBarry Smith } 116ceb9aaf7SBarry Smith 117