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