xref: /petsc/src/mat/impls/lrc/lrc.c (revision 38f2d2fdb3b6f522a3102c6eb796cebecf3224c0)
1ceb9aaf7SBarry Smith #define PETSCMAT_DLL
2ceb9aaf7SBarry Smith 
3b9147fbbSdalcinl #include "include/private/matimpl.h"          /*I "petscmat.h" I*/
4ab92ecdeSBarry Smith #include "src/mat/impls/dense/seq/dense.h"
5ceb9aaf7SBarry Smith 
6ceb9aaf7SBarry Smith typedef struct {
7ab92ecdeSBarry Smith   Mat         A,U,V;
8ab92ecdeSBarry Smith   Vec         work1,work2;/* Sequential (big) vectors that hold partial products */
9ab92ecdeSBarry Smith   PetscMPIInt nwork;      /* length of work vectors */
10ab92ecdeSBarry Smith } Mat_LRC;
11ab92ecdeSBarry Smith 
12ab92ecdeSBarry Smith 
13ceb9aaf7SBarry Smith 
14ceb9aaf7SBarry Smith #undef __FUNCT__
15ab92ecdeSBarry Smith #define __FUNCT__ "MatMult_LRC"
16ab92ecdeSBarry Smith PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y)
17ceb9aaf7SBarry Smith {
18ab92ecdeSBarry Smith   Mat_LRC        *Na = (Mat_LRC*)N->data;
19ceb9aaf7SBarry Smith   PetscErrorCode ierr;
20ab92ecdeSBarry Smith   PetscScalar    *w1,*w2;
21ceb9aaf7SBarry Smith 
22ceb9aaf7SBarry Smith   PetscFunctionBegin;
23ab92ecdeSBarry Smith   ierr = MatMult(Na->A,x,y);CHKERRQ(ierr);
24ab92ecdeSBarry Smith 
25ab92ecdeSBarry Smith   /* multiply the local part of V with the local part of x */
26ab92ecdeSBarry Smith   /* note in this call x is treated as a sequential vector  */
27ab92ecdeSBarry Smith   ierr = MatMultTranspose_SeqDense(Na->V,x,Na->work1);CHKERRQ(ierr);
28ab92ecdeSBarry Smith 
29ab92ecdeSBarry Smith   /* Form the sum of all the local multiplies : this is work2 = V'*x =
30ab92ecdeSBarry Smith      sum_{all processors} work1 */
31ab92ecdeSBarry Smith 
32ab92ecdeSBarry Smith   ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr);
33ab92ecdeSBarry Smith   ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr);
34ab92ecdeSBarry Smith   ierr = MPI_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,PetscSum_Op,N->comm);CHKERRQ(ierr);
35ab92ecdeSBarry Smith   ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr);
36ab92ecdeSBarry Smith   ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr);
37ab92ecdeSBarry Smith 
38ab92ecdeSBarry Smith   /* multiply-sub y = y  + U*work2 */
39ab92ecdeSBarry Smith   /* note in this call y is treated as a sequential vector  */
40ab92ecdeSBarry Smith   ierr = MatMultAdd_SeqDense(Na->U,Na->work2,y,y);CHKERRQ(ierr);
41ceb9aaf7SBarry Smith   PetscFunctionReturn(0);
42ceb9aaf7SBarry Smith }
43ceb9aaf7SBarry Smith 
44ceb9aaf7SBarry Smith #undef __FUNCT__
45ab92ecdeSBarry Smith #define __FUNCT__ "MatDestroy_LRC"
46ab92ecdeSBarry Smith PetscErrorCode MatDestroy_LRC(Mat N)
47ceb9aaf7SBarry Smith {
48ab92ecdeSBarry 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);
53ab92ecdeSBarry Smith   ierr = PetscObjectDereference((PetscObject)Na->U);CHKERRQ(ierr);
54ab92ecdeSBarry Smith   ierr = PetscObjectDereference((PetscObject)Na->V);CHKERRQ(ierr);
55ab92ecdeSBarry Smith   ierr = VecDestroy(Na->work1);CHKERRQ(ierr);
56ab92ecdeSBarry 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__
63ab92ecdeSBarry Smith #define __FUNCT__ "MatCreateLRC"
64ceb9aaf7SBarry Smith /*@
65ab92ecdeSBarry 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:
70ab92ecdeSBarry Smith +   A  - the (sparse) matrix
71ab92ecdeSBarry Smith -   U. V - two dense rectangular (tall and skinny) matrices
72ceb9aaf7SBarry Smith 
73ceb9aaf7SBarry Smith    Output Parameter:
74ab92ecdeSBarry Smith .   N - the matrix that represents A + U*V'
75ceb9aaf7SBarry Smith 
76ceb9aaf7SBarry Smith    Level: intermediate
77ceb9aaf7SBarry Smith 
78ab92ecdeSBarry Smith    Notes: The matrix A + U*V' formed! Rather the new matrix
79ceb9aaf7SBarry Smith           object performs the matrix-vector product by first multiplying by
80ab92ecdeSBarry Smith           A and then adding the other term
81ceb9aaf7SBarry Smith @*/
82ab92ecdeSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat A,Mat U, Mat V,Mat *N)
83ceb9aaf7SBarry Smith {
84ceb9aaf7SBarry Smith   PetscErrorCode ierr;
85ceb9aaf7SBarry Smith   PetscInt       m,n;
86ab92ecdeSBarry 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);
92ab92ecdeSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr);
93ceb9aaf7SBarry Smith 
94*38f2d2fdSLisandro Dalcin   ierr       = PetscNewLog(*N,Mat_LRC,&Na);CHKERRQ(ierr);
95*38f2d2fdSLisandro Dalcin   (*N)->data = (void*) Na;
96ceb9aaf7SBarry Smith   Na->A      = A;
97ab92ecdeSBarry Smith 
98ab92ecdeSBarry Smith   ierr      = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr);
99ab92ecdeSBarry Smith   ierr      = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr);
100ceb9aaf7SBarry Smith   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
101ab92ecdeSBarry Smith   ierr      = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr);
102ab92ecdeSBarry Smith   ierr      = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr);
103ceb9aaf7SBarry Smith 
104899cda47SBarry Smith   ierr                   = VecCreateSeq(PETSC_COMM_SELF,U->cmap.N,&Na->work1);CHKERRQ(ierr);
105ab92ecdeSBarry Smith   ierr                   = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr);
106899cda47SBarry Smith   Na->nwork              = U->cmap.N;
107ab92ecdeSBarry Smith 
108ab92ecdeSBarry Smith   (*N)->ops->destroy     = MatDestroy_LRC;
109ab92ecdeSBarry Smith   (*N)->ops->mult        = MatMult_LRC;
110ceb9aaf7SBarry Smith   (*N)->assembled        = PETSC_TRUE;
111899cda47SBarry Smith   (*N)->cmap.N                = A->cmap.N;
112899cda47SBarry Smith   (*N)->rmap.N                = A->cmap.N;
113899cda47SBarry Smith   (*N)->cmap.n                = A->cmap.n;
114899cda47SBarry Smith   (*N)->rmap.n                = A->cmap.n;
115ceb9aaf7SBarry Smith   PetscFunctionReturn(0);
116ceb9aaf7SBarry Smith }
117ceb9aaf7SBarry Smith 
118