xref: /petsc/src/mat/impls/lrc/lrc.c (revision d751fb92f56777f4a566f5e5fedf76f657538c08)
1ceb9aaf7SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
3ceb9aaf7SBarry Smith 
4*d751fb92SStefano Zampini PETSC_EXTERN PetscErrorCode VecGetRootType_Private(Vec,VecType*);
5*d751fb92SStefano Zampini 
6ceb9aaf7SBarry Smith typedef struct {
7986c4d72SJose E. Roman   Mat A;           /* sparse matrix */
8986c4d72SJose E. Roman   Mat U,V;         /* dense tall-skinny matrices */
9986c4d72SJose E. Roman   Vec c;           /* sequential vector containing the diagonal of C */
10986c4d72SJose E. Roman   Vec work1,work2; /* sequential vectors that hold partial products */
110575051bSJose E. Roman   Vec xl,yl;       /* auxiliary sequential vectors for matmult operation */
12ab92ecdeSBarry Smith } Mat_LRC;
13ab92ecdeSBarry Smith 
14*d751fb92SStefano Zampini static PetscErrorCode MatMult_LRC_kernel(Mat N,Vec x,Vec y,PetscBool transpose)
15ceb9aaf7SBarry Smith {
16ab92ecdeSBarry Smith   Mat_LRC        *Na = (Mat_LRC*)N->data;
17ceb9aaf7SBarry Smith   PetscErrorCode ierr;
18*d751fb92SStefano Zampini   PetscMPIInt    size;
19*d751fb92SStefano Zampini   Mat            U,V;
20ceb9aaf7SBarry Smith 
21ceb9aaf7SBarry Smith   PetscFunctionBegin;
22*d751fb92SStefano Zampini   U = transpose ? Na->V : Na->U;
23*d751fb92SStefano Zampini   V = transpose ? Na->U : Na->V;
24*d751fb92SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)N),&size);CHKERRMPI(ierr);
25*d751fb92SStefano Zampini   if (size == 1) {
26*d751fb92SStefano Zampini     ierr = MatMultHermitianTranspose(V,x,Na->work1);CHKERRQ(ierr);
27*d751fb92SStefano Zampini     if (Na->c) {
28*d751fb92SStefano Zampini       ierr = VecPointwiseMult(Na->work1,Na->c,Na->work1);CHKERRQ(ierr);
29*d751fb92SStefano Zampini     }
30*d751fb92SStefano Zampini     if (Na->A) {
31*d751fb92SStefano Zampini       if (transpose) {
32*d751fb92SStefano Zampini         ierr = MatMultTranspose(Na->A,x,y);CHKERRQ(ierr);
33*d751fb92SStefano Zampini       } else {
34*d751fb92SStefano Zampini         ierr = MatMult(Na->A,x,y);CHKERRQ(ierr);
35*d751fb92SStefano Zampini       }
36*d751fb92SStefano Zampini       ierr = MatMultAdd(U,Na->work1,y,y);CHKERRQ(ierr);
37*d751fb92SStefano Zampini     } else {
38*d751fb92SStefano Zampini       ierr = MatMult(U,Na->work1,y);CHKERRQ(ierr);
39*d751fb92SStefano Zampini     }
40*d751fb92SStefano Zampini   } else {
41*d751fb92SStefano Zampini     Mat               Uloc,Vloc;
42*d751fb92SStefano Zampini     Vec               yl,xl;
43*d751fb92SStefano Zampini     const PetscScalar *w1;
44*d751fb92SStefano Zampini     PetscScalar       *w2;
45*d751fb92SStefano Zampini     PetscInt          nwork;
46*d751fb92SStefano Zampini     PetscMPIInt       mpinwork;
47*d751fb92SStefano Zampini 
48*d751fb92SStefano Zampini     xl = transpose ? Na->yl : Na->xl;
49*d751fb92SStefano Zampini     yl = transpose ? Na->xl : Na->yl;
50*d751fb92SStefano Zampini     ierr = VecGetLocalVector(y,yl);CHKERRQ(ierr);
51*d751fb92SStefano Zampini     ierr = MatDenseGetLocalMatrix(U,&Uloc);CHKERRQ(ierr);
52*d751fb92SStefano Zampini     ierr = MatDenseGetLocalMatrix(V,&Vloc);CHKERRQ(ierr);
530575051bSJose E. Roman 
54ab92ecdeSBarry Smith     /* multiply the local part of V with the local part of x */
55*d751fb92SStefano Zampini     ierr = VecGetLocalVectorRead(x,xl);CHKERRQ(ierr);
56*d751fb92SStefano Zampini     ierr = MatMultHermitianTranspose(Vloc,xl,Na->work1);CHKERRQ(ierr);
57*d751fb92SStefano Zampini     ierr = VecRestoreLocalVectorRead(x,xl);CHKERRQ(ierr);
58ab92ecdeSBarry Smith 
5986eed97dSJose E. Roman     /* form the sum of all the local multiplies: this is work2 = V'*x =
60ab92ecdeSBarry Smith        sum_{all processors} work1 */
61*d751fb92SStefano Zampini     ierr = VecGetArrayRead(Na->work1,&w1);CHKERRQ(ierr);
62*d751fb92SStefano Zampini     ierr = VecGetArrayWrite(Na->work2,&w2);CHKERRQ(ierr);
63*d751fb92SStefano Zampini     ierr = VecGetLocalSize(Na->work1,&nwork);CHKERRQ(ierr);
64*d751fb92SStefano Zampini     ierr = PetscMPIIntCast(nwork,&mpinwork);CHKERRQ(ierr);
65*d751fb92SStefano Zampini     ierr = MPIU_Allreduce(w1,w2,mpinwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRMPI(ierr);
66*d751fb92SStefano Zampini     ierr = VecRestoreArrayRead(Na->work1,&w1);CHKERRQ(ierr);
67*d751fb92SStefano Zampini     ierr = VecRestoreArrayWrite(Na->work2,&w2);CHKERRQ(ierr);
68ab92ecdeSBarry Smith 
69986c4d72SJose E. Roman     if (Na->c) {  /* work2 = C*work2 */
70986c4d72SJose E. Roman       ierr = VecPointwiseMult(Na->work2,Na->c,Na->work2);CHKERRQ(ierr);
71986c4d72SJose E. Roman     }
72986c4d72SJose E. Roman 
7386eed97dSJose E. Roman     if (Na->A) {
74*d751fb92SStefano Zampini       /* form y = A*x or A^t*x */
75*d751fb92SStefano Zampini       if (transpose) {
76*d751fb92SStefano Zampini         ierr = MatMultTranspose(Na->A,x,y);CHKERRQ(ierr);
77*d751fb92SStefano Zampini       } else {
7886eed97dSJose E. Roman         ierr = MatMult(Na->A,x,y);CHKERRQ(ierr);
79*d751fb92SStefano Zampini       }
8086eed97dSJose E. Roman       /* multiply-add y = y + U*work2 */
81*d751fb92SStefano Zampini       ierr = MatMultAdd(Uloc,Na->work2,yl,yl);CHKERRQ(ierr);
8286eed97dSJose E. Roman     } else {
8386eed97dSJose E. Roman       /* multiply y = U*work2 */
84*d751fb92SStefano Zampini       ierr = MatMult(Uloc,Na->work2,yl);CHKERRQ(ierr);
8586eed97dSJose E. Roman     }
860575051bSJose E. Roman 
87*d751fb92SStefano Zampini     ierr = VecRestoreLocalVector(y,yl);CHKERRQ(ierr);
88*d751fb92SStefano Zampini   }
89ceb9aaf7SBarry Smith   PetscFunctionReturn(0);
90ceb9aaf7SBarry Smith }
91ceb9aaf7SBarry Smith 
92*d751fb92SStefano Zampini static PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y)
93*d751fb92SStefano Zampini {
94*d751fb92SStefano Zampini   PetscErrorCode ierr;
95*d751fb92SStefano Zampini 
96*d751fb92SStefano Zampini   PetscFunctionBegin;
97*d751fb92SStefano Zampini   ierr = MatMult_LRC_kernel(N,x,y,PETSC_FALSE);CHKERRQ(ierr);
98*d751fb92SStefano Zampini   PetscFunctionReturn(0);
99*d751fb92SStefano Zampini }
100*d751fb92SStefano Zampini 
101*d751fb92SStefano Zampini static PetscErrorCode MatMultTranspose_LRC(Mat N,Vec x,Vec y)
102*d751fb92SStefano Zampini {
103*d751fb92SStefano Zampini   PetscErrorCode ierr;
104*d751fb92SStefano Zampini 
105*d751fb92SStefano Zampini   PetscFunctionBegin;
106*d751fb92SStefano Zampini   ierr = MatMult_LRC_kernel(N,x,y,PETSC_TRUE);CHKERRQ(ierr);
107*d751fb92SStefano Zampini   PetscFunctionReturn(0);
108*d751fb92SStefano Zampini }
109*d751fb92SStefano Zampini 
110*d751fb92SStefano Zampini static PetscErrorCode MatDestroy_LRC(Mat N)
111ceb9aaf7SBarry Smith {
112ab92ecdeSBarry Smith   Mat_LRC        *Na = (Mat_LRC*)N->data;
113ceb9aaf7SBarry Smith   PetscErrorCode ierr;
114ceb9aaf7SBarry Smith 
115ceb9aaf7SBarry Smith   PetscFunctionBegin;
116bf0cc555SLisandro Dalcin   ierr = MatDestroy(&Na->A);CHKERRQ(ierr);
117bf0cc555SLisandro Dalcin   ierr = MatDestroy(&Na->U);CHKERRQ(ierr);
118bf0cc555SLisandro Dalcin   ierr = MatDestroy(&Na->V);CHKERRQ(ierr);
119986c4d72SJose E. Roman   ierr = VecDestroy(&Na->c);CHKERRQ(ierr);
1206bf464f9SBarry Smith   ierr = VecDestroy(&Na->work1);CHKERRQ(ierr);
1216bf464f9SBarry Smith   ierr = VecDestroy(&Na->work2);CHKERRQ(ierr);
1220575051bSJose E. Roman   ierr = VecDestroy(&Na->xl);CHKERRQ(ierr);
1230575051bSJose E. Roman   ierr = VecDestroy(&Na->yl);CHKERRQ(ierr);
124bf0cc555SLisandro Dalcin   ierr = PetscFree(N->data);CHKERRQ(ierr);
125*d751fb92SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)N,"MatLRCGetMats_C",NULL);CHKERRQ(ierr);
126267ca84cSJose E. Roman   PetscFunctionReturn(0);
127267ca84cSJose E. Roman }
128267ca84cSJose E. Roman 
129*d751fb92SStefano Zampini static PetscErrorCode MatLRCGetMats_LRC(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
130267ca84cSJose E. Roman {
131267ca84cSJose E. Roman   Mat_LRC *Na = (Mat_LRC*)N->data;
132267ca84cSJose E. Roman 
133267ca84cSJose E. Roman   PetscFunctionBegin;
134267ca84cSJose E. Roman   if (A) *A = Na->A;
135267ca84cSJose E. Roman   if (U) *U = Na->U;
136267ca84cSJose E. Roman   if (c) *c = Na->c;
137267ca84cSJose E. Roman   if (V) *V = Na->V;
138267ca84cSJose E. Roman   PetscFunctionReturn(0);
139267ca84cSJose E. Roman }
140267ca84cSJose E. Roman 
141267ca84cSJose E. Roman /*@
142267ca84cSJose E. Roman    MatLRCGetMats - Returns the constituents of an LRC matrix
143267ca84cSJose E. Roman 
144267ca84cSJose E. Roman    Collective on Mat
145267ca84cSJose E. Roman 
146267ca84cSJose E. Roman    Input Parameter:
147267ca84cSJose E. Roman .  N - matrix of type LRC
148267ca84cSJose E. Roman 
149267ca84cSJose E. Roman    Output Parameters:
150267ca84cSJose E. Roman +  A - the (sparse) matrix
1516b867d5aSJose E. Roman .  U - first dense rectangular (tall and skinny) matrix
1526b867d5aSJose E. Roman .  c - a sequential vector containing the diagonal of C
1536b867d5aSJose E. Roman -  V - second dense rectangular (tall and skinny) matrix
154267ca84cSJose E. Roman 
155267ca84cSJose E. Roman    Note:
156267ca84cSJose E. Roman    The returned matrices need not be destroyed by the caller.
157267ca84cSJose E. Roman 
158267ca84cSJose E. Roman    Level: intermediate
159267ca84cSJose E. Roman 
160267ca84cSJose E. Roman .seealso: MatCreateLRC()
161267ca84cSJose E. Roman @*/
162267ca84cSJose E. Roman PetscErrorCode MatLRCGetMats(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
163267ca84cSJose E. Roman {
164267ca84cSJose E. Roman   PetscErrorCode ierr;
165267ca84cSJose E. Roman 
166267ca84cSJose E. Roman   PetscFunctionBegin;
167267ca84cSJose E. Roman   ierr = PetscUseMethod(N,"MatLRCGetMats_C",(Mat,Mat*,Mat*,Vec*,Mat*),(N,A,U,c,V));CHKERRQ(ierr);
168ceb9aaf7SBarry Smith   PetscFunctionReturn(0);
169ceb9aaf7SBarry Smith }
170ceb9aaf7SBarry Smith 
171ceb9aaf7SBarry Smith /*@
172986c4d72SJose E. Roman    MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V'
173ceb9aaf7SBarry Smith 
174ceb9aaf7SBarry Smith    Collective on Mat
175ceb9aaf7SBarry Smith 
1769d9032efSJose E. Roman    Input Parameters:
17786eed97dSJose E. Roman +  A    - the (sparse) matrix (can be NULL)
178986c4d72SJose E. Roman .  U, V - two dense rectangular (tall and skinny) matrices
179*d751fb92SStefano Zampini -  c    - a vector containing the diagonal of C (can be NULL)
180ceb9aaf7SBarry Smith 
181ceb9aaf7SBarry Smith    Output Parameter:
182986c4d72SJose E. Roman .  N    - the matrix that represents A + U*C*V'
183ceb9aaf7SBarry Smith 
1849d9032efSJose E. Roman    Notes:
185986c4d72SJose E. Roman    The matrix A + U*C*V' is not formed! Rather the new matrix
186ceb9aaf7SBarry Smith    object performs the matrix-vector product by first multiplying by
1879d9032efSJose E. Roman    A and then adding the other term.
1889d9032efSJose E. Roman 
189986c4d72SJose E. Roman    C is a diagonal matrix (represented as a vector) of order k,
190986c4d72SJose E. Roman    where k is the number of columns of both U and V.
19186eed97dSJose E. Roman 
192986c4d72SJose E. Roman    If A is NULL then the new object behaves like a low-rank matrix U*C*V'.
193986c4d72SJose E. Roman 
194986c4d72SJose E. Roman    Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'.
195986c4d72SJose E. Roman 
196986c4d72SJose E. Roman    If c is NULL then the low-rank correction is just U*V'.
197*d751fb92SStefano Zampini    If a sequential c vector is used for a parallel matrix,
198*d751fb92SStefano Zampini    PETSc assumes that the values of the vector are consistently set across processors.
19986eed97dSJose E. Roman 
2009d9032efSJose E. Roman    Level: intermediate
201267ca84cSJose E. Roman 
202267ca84cSJose E. Roman .seealso: MatLRCGetMats()
203ceb9aaf7SBarry Smith @*/
204986c4d72SJose E. Roman PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N)
205ceb9aaf7SBarry Smith {
206ceb9aaf7SBarry Smith   PetscErrorCode ierr;
207986c4d72SJose E. Roman   PetscBool      match;
2089d9032efSJose E. Roman   PetscInt       m,n,k,m1,n1,k1;
209ab92ecdeSBarry Smith   Mat_LRC        *Na;
210*d751fb92SStefano Zampini   Mat            Uloc;
211*d751fb92SStefano Zampini   PetscMPIInt    size, csize = 0;
212ceb9aaf7SBarry Smith 
213ceb9aaf7SBarry Smith   PetscFunctionBegin;
21486eed97dSJose E. Roman   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1);
2159d9032efSJose E. Roman   PetscValidHeaderSpecific(U,MAT_CLASSID,2);
216986c4d72SJose E. Roman   if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3);
217986c4d72SJose E. Roman   if (V) {
218*d751fb92SStefano Zampini     PetscValidHeaderSpecific(V,MAT_CLASSID,4);
219*d751fb92SStefano Zampini     PetscCheckSameComm(U,2,V,4);
220*d751fb92SStefano Zampini   }
221*d751fb92SStefano Zampini   if (A) PetscCheckSameComm(A,1,U,2);
222*d751fb92SStefano Zampini 
223*d751fb92SStefano Zampini   if (!V) V = U;
224*d751fb92SStefano Zampini   ierr = PetscObjectBaseTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
225*d751fb92SStefano Zampini   PetscCheckFalse(!match,PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense, found %s",((PetscObject)U)->type_name);
226*d751fb92SStefano Zampini   ierr = PetscObjectBaseTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
227*d751fb92SStefano Zampini   PetscCheckFalse(!match,PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix V must be of type dense, found %s",((PetscObject)V)->type_name);
228*d751fb92SStefano Zampini   ierr = PetscStrcmp(U->defaultvectype,V->defaultvectype,&match);CHKERRQ(ierr);
229*d751fb92SStefano Zampini   PetscCheckFalse(!match,PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_WRONG,"Matrix U and V must have the same VecType %s != %s",U->defaultvectype,V->defaultvectype);
230*d751fb92SStefano Zampini   if (A) {
231*d751fb92SStefano Zampini     ierr = PetscStrcmp(A->defaultvectype,U->defaultvectype,&match);CHKERRQ(ierr);
232*d751fb92SStefano Zampini     PetscCheckFalse(!match,PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_WRONG,"Matrix A and U must have the same VecType %s != %s",A->defaultvectype,U->defaultvectype);
233986c4d72SJose E. Roman   }
2349d9032efSJose E. Roman 
235*d751fb92SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)U),&size);CHKERRMPI(ierr);
2369d9032efSJose E. Roman   ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr);
2379d9032efSJose E. Roman   ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr);
2382c71b3e2SJacob Faibussowitsch   PetscCheckFalse(k != k1,PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns (%" PetscInt_FMT " vs %" PetscInt_FMT ")",k,k1);
2399d9032efSJose E. Roman   ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr);
2409d9032efSJose E. Roman   ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr);
24186eed97dSJose E. Roman   if (A) {
2429d9032efSJose E. Roman     ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr);
243*d751fb92SStefano Zampini     PetscCheckFalse(m != m1,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local dimensions of U %" PetscInt_FMT " and A %" PetscInt_FMT " do not match",m,m1);
244*d751fb92SStefano Zampini     PetscCheckFalse(n != n1,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local dimensions of V %" PetscInt_FMT " and A %" PetscInt_FMT " do not match",n,n1);
24586eed97dSJose E. Roman   }
246986c4d72SJose E. Roman   if (c) {
247*d751fb92SStefano Zampini     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)c),&csize);CHKERRMPI(ierr);
248986c4d72SJose E. Roman     ierr = VecGetSize(c,&k1);CHKERRQ(ierr);
2492c71b3e2SJacob Faibussowitsch     PetscCheckFalse(k != k1,PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"The length of c %" PetscInt_FMT " does not match the number of columns of U and V (%" PetscInt_FMT ")",k1,k);
250*d751fb92SStefano Zampini     PetscCheckFalse(csize != 1 && csize != size, PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"U and c must have the same communicator size %d != %d",size,csize);
251986c4d72SJose E. Roman   }
2529d9032efSJose E. Roman 
25386eed97dSJose E. Roman   ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr);
2549d9032efSJose E. Roman   ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
255*d751fb92SStefano Zampini   ierr = MatSetVecType(*N,U->defaultvectype);CHKERRQ(ierr);
256ab92ecdeSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr);
257*d751fb92SStefano Zampini   /* Flag matrix as symmetric if A is symmetric and U == V */
258*d751fb92SStefano Zampini   ierr = MatSetOption(*N,MAT_SYMMETRIC,(PetscBool)((A ? A->symmetric : PETSC_TRUE) && U == V));CHKERRQ(ierr);
259ceb9aaf7SBarry Smith 
260b00a9115SJed Brown   ierr       = PetscNewLog(*N,&Na);CHKERRQ(ierr);
26138f2d2fdSLisandro Dalcin   (*N)->data = (void*)Na;
262ceb9aaf7SBarry Smith   Na->A      = A;
263a639b4dcSJose E. Roman   Na->U      = U;
264986c4d72SJose E. Roman   Na->c      = c;
265a639b4dcSJose E. Roman   Na->V      = V;
266ab92ecdeSBarry Smith 
267*d751fb92SStefano Zampini   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
268ab92ecdeSBarry Smith   ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr);
269ab92ecdeSBarry Smith   ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr);
270*d751fb92SStefano Zampini   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
271ceb9aaf7SBarry Smith 
272*d751fb92SStefano Zampini   ierr = MatDenseGetLocalMatrix(Na->U,&Uloc);CHKERRQ(ierr);
273*d751fb92SStefano Zampini   ierr = MatCreateVecs(Uloc,&Na->work1,NULL);CHKERRQ(ierr);
274*d751fb92SStefano Zampini   if (size != 1) {
275*d751fb92SStefano Zampini     Mat Vloc;
276*d751fb92SStefano Zampini 
277*d751fb92SStefano Zampini     if (Na->c && csize != 1) { /* scatter parallel vector to sequential */
278*d751fb92SStefano Zampini       VecScatter sct;
279*d751fb92SStefano Zampini 
280*d751fb92SStefano Zampini       ierr = VecScatterCreateToAll(Na->c,&sct,&c);CHKERRQ(ierr);
281*d751fb92SStefano Zampini       ierr = VecScatterBegin(sct,Na->c,c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
282*d751fb92SStefano Zampini       ierr = VecScatterEnd(sct,Na->c,c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
283*d751fb92SStefano Zampini       ierr = VecScatterDestroy(&sct);CHKERRQ(ierr);
284*d751fb92SStefano Zampini       ierr = VecDestroy(&Na->c);CHKERRQ(ierr);
285*d751fb92SStefano Zampini       ierr = PetscLogObjectParent((PetscObject)*N,(PetscObject)c);CHKERRQ(ierr);
286*d751fb92SStefano Zampini       Na->c = c;
287*d751fb92SStefano Zampini     }
288*d751fb92SStefano Zampini     ierr = MatDenseGetLocalMatrix(Na->V,&Vloc);CHKERRQ(ierr);
289ab92ecdeSBarry Smith     ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr);
290*d751fb92SStefano Zampini     ierr = MatCreateVecs(Vloc,NULL,&Na->xl);CHKERRQ(ierr);
291*d751fb92SStefano Zampini     ierr = MatCreateVecs(Uloc,NULL,&Na->yl);CHKERRQ(ierr);
292*d751fb92SStefano Zampini   }
293*d751fb92SStefano Zampini   ierr = PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->work1);CHKERRQ(ierr);
294*d751fb92SStefano Zampini   ierr = PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->work1);CHKERRQ(ierr);
295*d751fb92SStefano Zampini   ierr = PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->xl);CHKERRQ(ierr);
296*d751fb92SStefano Zampini   ierr = PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->yl);CHKERRQ(ierr);
297ab92ecdeSBarry Smith 
298*d751fb92SStefano Zampini   /* Internally create a scaling vector if roottypes do not match */
299*d751fb92SStefano Zampini   if (Na->c) {
300*d751fb92SStefano Zampini     VecType rt1,rt2;
301*d751fb92SStefano Zampini 
302*d751fb92SStefano Zampini     ierr = VecGetRootType_Private(Na->work1,&rt1);CHKERRQ(ierr);
303*d751fb92SStefano Zampini     ierr = VecGetRootType_Private(Na->c,&rt2);CHKERRQ(ierr);
304*d751fb92SStefano Zampini     ierr = PetscStrcmp(rt1,rt2,&match);CHKERRQ(ierr);
305*d751fb92SStefano Zampini     if (!match) {
306*d751fb92SStefano Zampini       ierr = VecDuplicate(Na->c,&c);CHKERRQ(ierr);
307*d751fb92SStefano Zampini       ierr = VecCopy(Na->c,c);CHKERRQ(ierr);
308*d751fb92SStefano Zampini       ierr = VecDestroy(&Na->c);CHKERRQ(ierr);
309*d751fb92SStefano Zampini       ierr = PetscLogObjectParent((PetscObject)*N,(PetscObject)c);CHKERRQ(ierr);
310*d751fb92SStefano Zampini       Na->c = c;
311*d751fb92SStefano Zampini     }
312*d751fb92SStefano Zampini   }
3130575051bSJose E. Roman 
314ab92ecdeSBarry Smith   (*N)->ops->destroy       = MatDestroy_LRC;
315ab92ecdeSBarry Smith   (*N)->ops->mult          = MatMult_LRC;
316*d751fb92SStefano Zampini   (*N)->ops->multtranspose = MatMultTranspose_LRC;
317*d751fb92SStefano Zampini 
318ceb9aaf7SBarry Smith   (*N)->assembled    = PETSC_TRUE;
31926c5da1bSJose E. Roman   (*N)->preallocated = PETSC_TRUE;
320267ca84cSJose E. Roman 
321267ca84cSJose E. Roman   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatLRCGetMats_C",MatLRCGetMats_LRC);CHKERRQ(ierr);
322*d751fb92SStefano Zampini   ierr = MatSetUp(*N);CHKERRQ(ierr);
323ceb9aaf7SBarry Smith   PetscFunctionReturn(0);
324ceb9aaf7SBarry Smith }
325