1ceb9aaf7SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3ceb9aaf7SBarry Smith 4d751fb92SStefano Zampini PETSC_EXTERN PetscErrorCode VecGetRootType_Private(Vec, VecType *); 5d751fb92SStefano 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 14d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_LRC_kernel(Mat N, Vec x, Vec y, PetscBool transpose) 15d71ae5a4SJacob Faibussowitsch { 16ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC *)N->data; 17d751fb92SStefano Zampini PetscMPIInt size; 18d751fb92SStefano Zampini Mat U, V; 19ceb9aaf7SBarry Smith 20ceb9aaf7SBarry Smith PetscFunctionBegin; 21d751fb92SStefano Zampini U = transpose ? Na->V : Na->U; 22d751fb92SStefano Zampini V = transpose ? Na->U : Na->V; 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)N), &size)); 24d751fb92SStefano Zampini if (size == 1) { 259566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(V, x, Na->work1)); 261baa6e33SBarry Smith if (Na->c) PetscCall(VecPointwiseMult(Na->work1, Na->c, Na->work1)); 27d751fb92SStefano Zampini if (Na->A) { 28d751fb92SStefano Zampini if (transpose) { 299566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A, x, y)); 30d751fb92SStefano Zampini } else { 319566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, x, y)); 32d751fb92SStefano Zampini } 339566063dSJacob Faibussowitsch PetscCall(MatMultAdd(U, Na->work1, y, y)); 34d751fb92SStefano Zampini } else { 359566063dSJacob Faibussowitsch PetscCall(MatMult(U, Na->work1, y)); 36d751fb92SStefano Zampini } 37d751fb92SStefano Zampini } else { 38d751fb92SStefano Zampini Mat Uloc, Vloc; 39d751fb92SStefano Zampini Vec yl, xl; 40d751fb92SStefano Zampini const PetscScalar *w1; 41d751fb92SStefano Zampini PetscScalar *w2; 42d751fb92SStefano Zampini PetscInt nwork; 43d751fb92SStefano Zampini PetscMPIInt mpinwork; 44d751fb92SStefano Zampini 45d751fb92SStefano Zampini xl = transpose ? Na->yl : Na->xl; 46d751fb92SStefano Zampini yl = transpose ? Na->xl : Na->yl; 479566063dSJacob Faibussowitsch PetscCall(VecGetLocalVector(y, yl)); 489566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(U, &Uloc)); 499566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(V, &Vloc)); 500575051bSJose E. Roman 51ab92ecdeSBarry Smith /* multiply the local part of V with the local part of x */ 529566063dSJacob Faibussowitsch PetscCall(VecGetLocalVectorRead(x, xl)); 539566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Vloc, xl, Na->work1)); 549566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVectorRead(x, xl)); 55ab92ecdeSBarry Smith 5686eed97dSJose E. Roman /* form the sum of all the local multiplies: this is work2 = V'*x = 57ab92ecdeSBarry Smith sum_{all processors} work1 */ 589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Na->work1, &w1)); 599566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(Na->work2, &w2)); 609566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Na->work1, &nwork)); 619566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nwork, &mpinwork)); 621c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(w1, w2, mpinwork, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N))); 639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Na->work1, &w1)); 649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(Na->work2, &w2)); 65ab92ecdeSBarry Smith 66986c4d72SJose E. Roman if (Na->c) { /* work2 = C*work2 */ 679566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->work2, Na->c, Na->work2)); 68986c4d72SJose E. Roman } 69986c4d72SJose E. Roman 7086eed97dSJose E. Roman if (Na->A) { 71d751fb92SStefano Zampini /* form y = A*x or A^t*x */ 72d751fb92SStefano Zampini if (transpose) { 739566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A, x, y)); 74d751fb92SStefano Zampini } else { 759566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, x, y)); 76d751fb92SStefano Zampini } 7786eed97dSJose E. Roman /* multiply-add y = y + U*work2 */ 789566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Uloc, Na->work2, yl, yl)); 7986eed97dSJose E. Roman } else { 8086eed97dSJose E. Roman /* multiply y = U*work2 */ 819566063dSJacob Faibussowitsch PetscCall(MatMult(Uloc, Na->work2, yl)); 8286eed97dSJose E. Roman } 830575051bSJose E. Roman 849566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVector(y, yl)); 85d751fb92SStefano Zampini } 86ceb9aaf7SBarry Smith PetscFunctionReturn(0); 87ceb9aaf7SBarry Smith } 88ceb9aaf7SBarry Smith 89d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_LRC(Mat N, Vec x, Vec y) 90d71ae5a4SJacob Faibussowitsch { 91d751fb92SStefano Zampini PetscFunctionBegin; 929566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N, x, y, PETSC_FALSE)); 93d751fb92SStefano Zampini PetscFunctionReturn(0); 94d751fb92SStefano Zampini } 95d751fb92SStefano Zampini 96d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_LRC(Mat N, Vec x, Vec y) 97d71ae5a4SJacob Faibussowitsch { 98d751fb92SStefano Zampini PetscFunctionBegin; 999566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N, x, y, PETSC_TRUE)); 100d751fb92SStefano Zampini PetscFunctionReturn(0); 101d751fb92SStefano Zampini } 102d751fb92SStefano Zampini 103d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_LRC(Mat N) 104d71ae5a4SJacob Faibussowitsch { 105ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC *)N->data; 106ceb9aaf7SBarry Smith 107ceb9aaf7SBarry Smith PetscFunctionBegin; 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->U)); 1109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->V)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work1)); 1139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work2)); 1149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->xl)); 1159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->yl)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 1179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCGetMats_C", NULL)); 118267ca84cSJose E. Roman PetscFunctionReturn(0); 119267ca84cSJose E. Roman } 120267ca84cSJose E. Roman 121d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLRCGetMats_LRC(Mat N, Mat *A, Mat *U, Vec *c, Mat *V) 122d71ae5a4SJacob Faibussowitsch { 123267ca84cSJose E. Roman Mat_LRC *Na = (Mat_LRC *)N->data; 124267ca84cSJose E. Roman 125267ca84cSJose E. Roman PetscFunctionBegin; 126267ca84cSJose E. Roman if (A) *A = Na->A; 127267ca84cSJose E. Roman if (U) *U = Na->U; 128267ca84cSJose E. Roman if (c) *c = Na->c; 129267ca84cSJose E. Roman if (V) *V = Na->V; 130267ca84cSJose E. Roman PetscFunctionReturn(0); 131267ca84cSJose E. Roman } 132267ca84cSJose E. Roman 133267ca84cSJose E. Roman /*@ 134267ca84cSJose E. Roman MatLRCGetMats - Returns the constituents of an LRC matrix 135267ca84cSJose E. Roman 136*c3339decSBarry Smith Collective 137267ca84cSJose E. Roman 138267ca84cSJose E. Roman Input Parameter: 13911a5261eSBarry Smith . N - matrix of type `MATLRC` 140267ca84cSJose E. Roman 141267ca84cSJose E. Roman Output Parameters: 142267ca84cSJose E. Roman + A - the (sparse) matrix 1436b867d5aSJose E. Roman . U - first dense rectangular (tall and skinny) matrix 1446b867d5aSJose E. Roman . c - a sequential vector containing the diagonal of C 1456b867d5aSJose E. Roman - V - second dense rectangular (tall and skinny) matrix 146267ca84cSJose E. Roman 147267ca84cSJose E. Roman Note: 148267ca84cSJose E. Roman The returned matrices need not be destroyed by the caller. 149267ca84cSJose E. Roman 150267ca84cSJose E. Roman Level: intermediate 151267ca84cSJose E. Roman 15211a5261eSBarry Smith .seealso: `MATLRC`, `MatCreateLRC()` 153267ca84cSJose E. Roman @*/ 154d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLRCGetMats(Mat N, Mat *A, Mat *U, Vec *c, Mat *V) 155d71ae5a4SJacob Faibussowitsch { 156267ca84cSJose E. Roman PetscFunctionBegin; 157cac4c232SBarry Smith PetscUseMethod(N, "MatLRCGetMats_C", (Mat, Mat *, Mat *, Vec *, Mat *), (N, A, U, c, V)); 158ceb9aaf7SBarry Smith PetscFunctionReturn(0); 159ceb9aaf7SBarry Smith } 160ceb9aaf7SBarry Smith 16111a5261eSBarry Smith /*MC 16211a5261eSBarry Smith MATLRC - "lrc" - a matrix object that behaves like A + U*C*V' 163ceb9aaf7SBarry Smith 16411a5261eSBarry Smith Note: 16511a5261eSBarry Smith The matrix A + U*C*V' is not formed! Rather the matrix object performs the matrix-vector product `MatMult()`, by first multiplying by 16611a5261eSBarry Smith A and then adding the other term. 16711a5261eSBarry Smith 16811a5261eSBarry Smith Level: advanced 16911a5261eSBarry Smith 17011a5261eSBarry Smith .seealso: `MatCreateLRC` 17111a5261eSBarry Smith M*/ 17211a5261eSBarry Smith 17311a5261eSBarry Smith /*@ 17411a5261eSBarry Smith MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V' of type `MATLRC` 17511a5261eSBarry Smith 176*c3339decSBarry Smith Collective 177ceb9aaf7SBarry Smith 1789d9032efSJose E. Roman Input Parameters: 17986eed97dSJose E. Roman + A - the (sparse) matrix (can be NULL) 180986c4d72SJose E. Roman . U, V - two dense rectangular (tall and skinny) matrices 181d751fb92SStefano Zampini - c - a vector containing the diagonal of C (can be NULL) 182ceb9aaf7SBarry Smith 183ceb9aaf7SBarry Smith Output Parameter: 184986c4d72SJose E. Roman . N - the matrix that represents A + U*C*V' 185ceb9aaf7SBarry Smith 1869d9032efSJose E. Roman Notes: 187986c4d72SJose E. Roman The matrix A + U*C*V' is not formed! Rather the new matrix 18811a5261eSBarry Smith object performs the matrix-vector product `MatMult()`, by first multiplying by 1899d9032efSJose E. Roman A and then adding the other term. 1909d9032efSJose E. Roman 191986c4d72SJose E. Roman C is a diagonal matrix (represented as a vector) of order k, 192986c4d72SJose E. Roman where k is the number of columns of both U and V. 19386eed97dSJose E. Roman 194986c4d72SJose E. Roman If A is NULL then the new object behaves like a low-rank matrix U*C*V'. 195986c4d72SJose E. Roman 196986c4d72SJose E. Roman Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'. 197986c4d72SJose E. Roman 198986c4d72SJose E. Roman If c is NULL then the low-rank correction is just U*V'. 199d751fb92SStefano Zampini If a sequential c vector is used for a parallel matrix, 200d751fb92SStefano Zampini PETSc assumes that the values of the vector are consistently set across processors. 20186eed97dSJose E. Roman 2029d9032efSJose E. Roman Level: intermediate 203267ca84cSJose E. Roman 20411a5261eSBarry Smith .seealso: `MATLRC`, `MatLRCGetMats()` 205ceb9aaf7SBarry Smith @*/ 206d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateLRC(Mat A, Mat U, Vec c, Mat V, Mat *N) 207d71ae5a4SJacob Faibussowitsch { 208986c4d72SJose E. Roman PetscBool match; 2099d9032efSJose E. Roman PetscInt m, n, k, m1, n1, k1; 210ab92ecdeSBarry Smith Mat_LRC *Na; 211d751fb92SStefano Zampini Mat Uloc; 212d751fb92SStefano Zampini PetscMPIInt size, csize = 0; 213ceb9aaf7SBarry Smith 214ceb9aaf7SBarry Smith PetscFunctionBegin; 21586eed97dSJose E. Roman if (A) PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2169d9032efSJose E. Roman PetscValidHeaderSpecific(U, MAT_CLASSID, 2); 217986c4d72SJose E. Roman if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 3); 218986c4d72SJose E. Roman if (V) { 219d751fb92SStefano Zampini PetscValidHeaderSpecific(V, MAT_CLASSID, 4); 220d751fb92SStefano Zampini PetscCheckSameComm(U, 2, V, 4); 221d751fb92SStefano Zampini } 222d751fb92SStefano Zampini if (A) PetscCheckSameComm(A, 1, U, 2); 223d751fb92SStefano Zampini 224d751fb92SStefano Zampini if (!V) V = U; 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)U, &match, MATSEQDENSE, MATMPIDENSE, "")); 22628b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix U must be of type dense, found %s", ((PetscObject)U)->type_name); 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)V, &match, MATSEQDENSE, MATMPIDENSE, "")); 22828b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix V must be of type dense, found %s", ((PetscObject)V)->type_name); 2299566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(U->defaultvectype, V->defaultvectype, &match)); 23028b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_ARG_WRONG, "Matrix U and V must have the same VecType %s != %s", U->defaultvectype, V->defaultvectype); 231d751fb92SStefano Zampini if (A) { 2329566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(A->defaultvectype, U->defaultvectype, &match)); 23328b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_ARG_WRONG, "Matrix A and U must have the same VecType %s != %s", A->defaultvectype, U->defaultvectype); 234986c4d72SJose E. Roman } 2359d9032efSJose E. Roman 2369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)U), &size)); 2379566063dSJacob Faibussowitsch PetscCall(MatGetSize(U, NULL, &k)); 2389566063dSJacob Faibussowitsch PetscCall(MatGetSize(V, NULL, &k1)); 23908401ef6SPierre Jolivet PetscCheck(k == k1, PetscObjectComm((PetscObject)U), PETSC_ERR_ARG_INCOMP, "U and V have different number of columns (%" PetscInt_FMT " vs %" PetscInt_FMT ")", k, k1); 2409566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(U, &m, NULL)); 2419566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(V, &n, NULL)); 24286eed97dSJose E. Roman if (A) { 2439566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m1, &n1)); 24408401ef6SPierre Jolivet PetscCheck(m == m1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local dimensions of U %" PetscInt_FMT " and A %" PetscInt_FMT " do not match", m, m1); 24508401ef6SPierre Jolivet PetscCheck(n == n1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local dimensions of V %" PetscInt_FMT " and A %" PetscInt_FMT " do not match", n, n1); 24686eed97dSJose E. Roman } 247986c4d72SJose E. Roman if (c) { 2489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)c), &csize)); 2499566063dSJacob Faibussowitsch PetscCall(VecGetSize(c, &k1)); 25008401ef6SPierre Jolivet PetscCheck(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); 251aed4548fSBarry Smith PetscCheck(csize == 1 || csize == size, PetscObjectComm((PetscObject)c), PETSC_ERR_ARG_INCOMP, "U and c must have the same communicator size %d != %d", size, csize); 252986c4d72SJose E. Roman } 2539d9032efSJose E. Roman 2549566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)U), N)); 2559566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N, m, n, PETSC_DECIDE, PETSC_DECIDE)); 2569566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N, U->defaultvectype)); 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATLRC)); 258d751fb92SStefano Zampini /* Flag matrix as symmetric if A is symmetric and U == V */ 259b94d7dedSBarry Smith PetscCall(MatSetOption(*N, MAT_SYMMETRIC, (PetscBool)((A ? A->symmetric == PETSC_BOOL3_TRUE : PETSC_TRUE) && U == V))); 260ceb9aaf7SBarry Smith 2614dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na)); 26238f2d2fdSLisandro Dalcin (*N)->data = (void *)Na; 263ceb9aaf7SBarry Smith Na->A = A; 264a639b4dcSJose E. Roman Na->U = U; 265986c4d72SJose E. Roman Na->c = c; 266a639b4dcSJose E. Roman Na->V = V; 267ab92ecdeSBarry Smith 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Na->U)); 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Na->V)); 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 272ceb9aaf7SBarry Smith 2739566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Na->U, &Uloc)); 2749566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Uloc, &Na->work1, NULL)); 275d751fb92SStefano Zampini if (size != 1) { 276d751fb92SStefano Zampini Mat Vloc; 277d751fb92SStefano Zampini 278d751fb92SStefano Zampini if (Na->c && csize != 1) { /* scatter parallel vector to sequential */ 279d751fb92SStefano Zampini VecScatter sct; 280d751fb92SStefano Zampini 2819566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToAll(Na->c, &sct, &c)); 2829566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 2839566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 2849566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sct)); 2859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 286d751fb92SStefano Zampini Na->c = c; 287d751fb92SStefano Zampini } 2889566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Na->V, &Vloc)); 2899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->work1, &Na->work2)); 2909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Vloc, NULL, &Na->xl)); 2919566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Uloc, NULL, &Na->yl)); 292d751fb92SStefano Zampini } 293ab92ecdeSBarry Smith 294d751fb92SStefano Zampini /* Internally create a scaling vector if roottypes do not match */ 295d751fb92SStefano Zampini if (Na->c) { 296d751fb92SStefano Zampini VecType rt1, rt2; 297d751fb92SStefano Zampini 2989566063dSJacob Faibussowitsch PetscCall(VecGetRootType_Private(Na->work1, &rt1)); 2999566063dSJacob Faibussowitsch PetscCall(VecGetRootType_Private(Na->c, &rt2)); 3009566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(rt1, rt2, &match)); 301d751fb92SStefano Zampini if (!match) { 3029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->c, &c)); 3039566063dSJacob Faibussowitsch PetscCall(VecCopy(Na->c, c)); 3049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 305d751fb92SStefano Zampini Na->c = c; 306d751fb92SStefano Zampini } 307d751fb92SStefano Zampini } 3080575051bSJose E. Roman 309ab92ecdeSBarry Smith (*N)->ops->destroy = MatDestroy_LRC; 310ab92ecdeSBarry Smith (*N)->ops->mult = MatMult_LRC; 311d751fb92SStefano Zampini (*N)->ops->multtranspose = MatMultTranspose_LRC; 312d751fb92SStefano Zampini 313ceb9aaf7SBarry Smith (*N)->assembled = PETSC_TRUE; 31426c5da1bSJose E. Roman (*N)->preallocated = PETSC_TRUE; 315267ca84cSJose E. Roman 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatLRCGetMats_C", MatLRCGetMats_LRC)); 3179566063dSJacob Faibussowitsch PetscCall(MatSetUp(*N)); 318ceb9aaf7SBarry Smith PetscFunctionReturn(0); 319ceb9aaf7SBarry Smith } 320