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 149371c9d4SSatish Balay static PetscErrorCode MatMult_LRC_kernel(Mat N, Vec x, Vec y, PetscBool transpose) { 15ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC *)N->data; 16d751fb92SStefano Zampini PetscMPIInt size; 17d751fb92SStefano Zampini Mat U, V; 18ceb9aaf7SBarry Smith 19ceb9aaf7SBarry Smith PetscFunctionBegin; 20d751fb92SStefano Zampini U = transpose ? Na->V : Na->U; 21d751fb92SStefano Zampini V = transpose ? Na->U : Na->V; 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)N), &size)); 23d751fb92SStefano Zampini if (size == 1) { 249566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(V, x, Na->work1)); 251baa6e33SBarry Smith if (Na->c) PetscCall(VecPointwiseMult(Na->work1, Na->c, Na->work1)); 26d751fb92SStefano Zampini if (Na->A) { 27d751fb92SStefano Zampini if (transpose) { 289566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A, x, y)); 29d751fb92SStefano Zampini } else { 309566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, x, y)); 31d751fb92SStefano Zampini } 329566063dSJacob Faibussowitsch PetscCall(MatMultAdd(U, Na->work1, y, y)); 33d751fb92SStefano Zampini } else { 349566063dSJacob Faibussowitsch PetscCall(MatMult(U, Na->work1, y)); 35d751fb92SStefano Zampini } 36d751fb92SStefano Zampini } else { 37d751fb92SStefano Zampini Mat Uloc, Vloc; 38d751fb92SStefano Zampini Vec yl, xl; 39d751fb92SStefano Zampini const PetscScalar *w1; 40d751fb92SStefano Zampini PetscScalar *w2; 41d751fb92SStefano Zampini PetscInt nwork; 42d751fb92SStefano Zampini PetscMPIInt mpinwork; 43d751fb92SStefano Zampini 44d751fb92SStefano Zampini xl = transpose ? Na->yl : Na->xl; 45d751fb92SStefano Zampini yl = transpose ? Na->xl : Na->yl; 469566063dSJacob Faibussowitsch PetscCall(VecGetLocalVector(y, yl)); 479566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(U, &Uloc)); 489566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(V, &Vloc)); 490575051bSJose E. Roman 50ab92ecdeSBarry Smith /* multiply the local part of V with the local part of x */ 519566063dSJacob Faibussowitsch PetscCall(VecGetLocalVectorRead(x, xl)); 529566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Vloc, xl, Na->work1)); 539566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVectorRead(x, xl)); 54ab92ecdeSBarry Smith 5586eed97dSJose E. Roman /* form the sum of all the local multiplies: this is work2 = V'*x = 56ab92ecdeSBarry Smith sum_{all processors} work1 */ 579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Na->work1, &w1)); 589566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(Na->work2, &w2)); 599566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Na->work1, &nwork)); 609566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nwork, &mpinwork)); 611c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(w1, w2, mpinwork, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N))); 629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Na->work1, &w1)); 639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(Na->work2, &w2)); 64ab92ecdeSBarry Smith 65986c4d72SJose E. Roman if (Na->c) { /* work2 = C*work2 */ 669566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->work2, Na->c, Na->work2)); 67986c4d72SJose E. Roman } 68986c4d72SJose E. Roman 6986eed97dSJose E. Roman if (Na->A) { 70d751fb92SStefano Zampini /* form y = A*x or A^t*x */ 71d751fb92SStefano Zampini if (transpose) { 729566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A, x, y)); 73d751fb92SStefano Zampini } else { 749566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, x, y)); 75d751fb92SStefano Zampini } 7686eed97dSJose E. Roman /* multiply-add y = y + U*work2 */ 779566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Uloc, Na->work2, yl, yl)); 7886eed97dSJose E. Roman } else { 7986eed97dSJose E. Roman /* multiply y = U*work2 */ 809566063dSJacob Faibussowitsch PetscCall(MatMult(Uloc, Na->work2, yl)); 8186eed97dSJose E. Roman } 820575051bSJose E. Roman 839566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVector(y, yl)); 84d751fb92SStefano Zampini } 85ceb9aaf7SBarry Smith PetscFunctionReturn(0); 86ceb9aaf7SBarry Smith } 87ceb9aaf7SBarry Smith 889371c9d4SSatish Balay static PetscErrorCode MatMult_LRC(Mat N, Vec x, Vec y) { 89d751fb92SStefano Zampini PetscFunctionBegin; 909566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N, x, y, PETSC_FALSE)); 91d751fb92SStefano Zampini PetscFunctionReturn(0); 92d751fb92SStefano Zampini } 93d751fb92SStefano Zampini 949371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_LRC(Mat N, Vec x, Vec y) { 95d751fb92SStefano Zampini PetscFunctionBegin; 969566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N, x, y, PETSC_TRUE)); 97d751fb92SStefano Zampini PetscFunctionReturn(0); 98d751fb92SStefano Zampini } 99d751fb92SStefano Zampini 1009371c9d4SSatish Balay static PetscErrorCode MatDestroy_LRC(Mat N) { 101ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC *)N->data; 102ceb9aaf7SBarry Smith 103ceb9aaf7SBarry Smith PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->U)); 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->V)); 1079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work1)); 1099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work2)); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->xl)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->yl)); 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCGetMats_C", NULL)); 114267ca84cSJose E. Roman PetscFunctionReturn(0); 115267ca84cSJose E. Roman } 116267ca84cSJose E. Roman 1179371c9d4SSatish Balay static PetscErrorCode MatLRCGetMats_LRC(Mat N, Mat *A, Mat *U, Vec *c, Mat *V) { 118267ca84cSJose E. Roman Mat_LRC *Na = (Mat_LRC *)N->data; 119267ca84cSJose E. Roman 120267ca84cSJose E. Roman PetscFunctionBegin; 121267ca84cSJose E. Roman if (A) *A = Na->A; 122267ca84cSJose E. Roman if (U) *U = Na->U; 123267ca84cSJose E. Roman if (c) *c = Na->c; 124267ca84cSJose E. Roman if (V) *V = Na->V; 125267ca84cSJose E. Roman PetscFunctionReturn(0); 126267ca84cSJose E. Roman } 127267ca84cSJose E. Roman 128267ca84cSJose E. Roman /*@ 129267ca84cSJose E. Roman MatLRCGetMats - Returns the constituents of an LRC matrix 130267ca84cSJose E. Roman 13111a5261eSBarry Smith Collective on N 132267ca84cSJose E. Roman 133267ca84cSJose E. Roman Input Parameter: 13411a5261eSBarry Smith . N - matrix of type `MATLRC` 135267ca84cSJose E. Roman 136267ca84cSJose E. Roman Output Parameters: 137267ca84cSJose E. Roman + A - the (sparse) matrix 1386b867d5aSJose E. Roman . U - first dense rectangular (tall and skinny) matrix 1396b867d5aSJose E. Roman . c - a sequential vector containing the diagonal of C 1406b867d5aSJose E. Roman - V - second dense rectangular (tall and skinny) matrix 141267ca84cSJose E. Roman 142267ca84cSJose E. Roman Note: 143267ca84cSJose E. Roman The returned matrices need not be destroyed by the caller. 144267ca84cSJose E. Roman 145267ca84cSJose E. Roman Level: intermediate 146267ca84cSJose E. Roman 14711a5261eSBarry Smith .seealso: `MATLRC`, `MatCreateLRC()` 148267ca84cSJose E. Roman @*/ 1499371c9d4SSatish Balay PetscErrorCode MatLRCGetMats(Mat N, Mat *A, Mat *U, Vec *c, Mat *V) { 150267ca84cSJose E. Roman PetscFunctionBegin; 151cac4c232SBarry Smith PetscUseMethod(N, "MatLRCGetMats_C", (Mat, Mat *, Mat *, Vec *, Mat *), (N, A, U, c, V)); 152ceb9aaf7SBarry Smith PetscFunctionReturn(0); 153ceb9aaf7SBarry Smith } 154ceb9aaf7SBarry Smith 15511a5261eSBarry Smith /*MC 15611a5261eSBarry Smith MATLRC - "lrc" - a matrix object that behaves like A + U*C*V' 157ceb9aaf7SBarry Smith 15811a5261eSBarry Smith Note: 15911a5261eSBarry Smith The matrix A + U*C*V' is not formed! Rather the matrix object performs the matrix-vector product `MatMult()`, by first multiplying by 16011a5261eSBarry Smith A and then adding the other term. 16111a5261eSBarry Smith 16211a5261eSBarry Smith Level: advanced 16311a5261eSBarry Smith 16411a5261eSBarry Smith .seealso: `MatCreateLRC` 16511a5261eSBarry Smith M*/ 16611a5261eSBarry Smith 16711a5261eSBarry Smith /*@ 16811a5261eSBarry Smith MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V' of type `MATLRC` 16911a5261eSBarry Smith 17011a5261eSBarry Smith Collective on A 171ceb9aaf7SBarry Smith 1729d9032efSJose E. Roman Input Parameters: 17386eed97dSJose E. Roman + A - the (sparse) matrix (can be NULL) 174986c4d72SJose E. Roman . U, V - two dense rectangular (tall and skinny) matrices 175d751fb92SStefano Zampini - c - a vector containing the diagonal of C (can be NULL) 176ceb9aaf7SBarry Smith 177ceb9aaf7SBarry Smith Output Parameter: 178986c4d72SJose E. Roman . N - the matrix that represents A + U*C*V' 179ceb9aaf7SBarry Smith 1809d9032efSJose E. Roman Notes: 181986c4d72SJose E. Roman The matrix A + U*C*V' is not formed! Rather the new matrix 18211a5261eSBarry Smith object performs the matrix-vector product `MatMult()`, by first multiplying by 1839d9032efSJose E. Roman A and then adding the other term. 1849d9032efSJose E. Roman 185986c4d72SJose E. Roman C is a diagonal matrix (represented as a vector) of order k, 186986c4d72SJose E. Roman where k is the number of columns of both U and V. 18786eed97dSJose E. Roman 188986c4d72SJose E. Roman If A is NULL then the new object behaves like a low-rank matrix U*C*V'. 189986c4d72SJose E. Roman 190986c4d72SJose E. Roman Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'. 191986c4d72SJose E. Roman 192986c4d72SJose E. Roman If c is NULL then the low-rank correction is just U*V'. 193d751fb92SStefano Zampini If a sequential c vector is used for a parallel matrix, 194d751fb92SStefano Zampini PETSc assumes that the values of the vector are consistently set across processors. 19586eed97dSJose E. Roman 1969d9032efSJose E. Roman Level: intermediate 197267ca84cSJose E. Roman 19811a5261eSBarry Smith .seealso: `MATLRC`, `MatLRCGetMats()` 199ceb9aaf7SBarry Smith @*/ 2009371c9d4SSatish Balay PetscErrorCode MatCreateLRC(Mat A, Mat U, Vec c, Mat V, Mat *N) { 201986c4d72SJose E. Roman PetscBool match; 2029d9032efSJose E. Roman PetscInt m, n, k, m1, n1, k1; 203ab92ecdeSBarry Smith Mat_LRC *Na; 204d751fb92SStefano Zampini Mat Uloc; 205d751fb92SStefano Zampini PetscMPIInt size, csize = 0; 206ceb9aaf7SBarry Smith 207ceb9aaf7SBarry Smith PetscFunctionBegin; 20886eed97dSJose E. Roman if (A) PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2099d9032efSJose E. Roman PetscValidHeaderSpecific(U, MAT_CLASSID, 2); 210986c4d72SJose E. Roman if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 3); 211986c4d72SJose E. Roman if (V) { 212d751fb92SStefano Zampini PetscValidHeaderSpecific(V, MAT_CLASSID, 4); 213d751fb92SStefano Zampini PetscCheckSameComm(U, 2, V, 4); 214d751fb92SStefano Zampini } 215d751fb92SStefano Zampini if (A) PetscCheckSameComm(A, 1, U, 2); 216d751fb92SStefano Zampini 217d751fb92SStefano Zampini if (!V) V = U; 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)U, &match, MATSEQDENSE, MATMPIDENSE, "")); 21928b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix U must be of type dense, found %s", ((PetscObject)U)->type_name); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)V, &match, MATSEQDENSE, MATMPIDENSE, "")); 22128b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix V must be of type dense, found %s", ((PetscObject)V)->type_name); 2229566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(U->defaultvectype, V->defaultvectype, &match)); 22328b400f6SJacob 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); 224d751fb92SStefano Zampini if (A) { 2259566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(A->defaultvectype, U->defaultvectype, &match)); 22628b400f6SJacob 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); 227986c4d72SJose E. Roman } 2289d9032efSJose E. Roman 2299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)U), &size)); 2309566063dSJacob Faibussowitsch PetscCall(MatGetSize(U, NULL, &k)); 2319566063dSJacob Faibussowitsch PetscCall(MatGetSize(V, NULL, &k1)); 23208401ef6SPierre 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); 2339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(U, &m, NULL)); 2349566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(V, &n, NULL)); 23586eed97dSJose E. Roman if (A) { 2369566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m1, &n1)); 23708401ef6SPierre 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); 23808401ef6SPierre 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); 23986eed97dSJose E. Roman } 240986c4d72SJose E. Roman if (c) { 2419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)c), &csize)); 2429566063dSJacob Faibussowitsch PetscCall(VecGetSize(c, &k1)); 24308401ef6SPierre 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); 244aed4548fSBarry 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); 245986c4d72SJose E. Roman } 2469d9032efSJose E. Roman 2479566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)U), N)); 2489566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N, m, n, PETSC_DECIDE, PETSC_DECIDE)); 2499566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N, U->defaultvectype)); 2509566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATLRC)); 251d751fb92SStefano Zampini /* Flag matrix as symmetric if A is symmetric and U == V */ 252b94d7dedSBarry Smith PetscCall(MatSetOption(*N, MAT_SYMMETRIC, (PetscBool)((A ? A->symmetric == PETSC_BOOL3_TRUE : PETSC_TRUE) && U == V))); 253ceb9aaf7SBarry Smith 254*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na)); 25538f2d2fdSLisandro Dalcin (*N)->data = (void *)Na; 256ceb9aaf7SBarry Smith Na->A = A; 257a639b4dcSJose E. Roman Na->U = U; 258986c4d72SJose E. Roman Na->c = c; 259a639b4dcSJose E. Roman Na->V = V; 260ab92ecdeSBarry Smith 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Na->U)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Na->V)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 265ceb9aaf7SBarry Smith 2669566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Na->U, &Uloc)); 2679566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Uloc, &Na->work1, NULL)); 268d751fb92SStefano Zampini if (size != 1) { 269d751fb92SStefano Zampini Mat Vloc; 270d751fb92SStefano Zampini 271d751fb92SStefano Zampini if (Na->c && csize != 1) { /* scatter parallel vector to sequential */ 272d751fb92SStefano Zampini VecScatter sct; 273d751fb92SStefano Zampini 2749566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToAll(Na->c, &sct, &c)); 2759566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 2769566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 2779566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sct)); 2789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 279d751fb92SStefano Zampini Na->c = c; 280d751fb92SStefano Zampini } 2819566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Na->V, &Vloc)); 2829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->work1, &Na->work2)); 2839566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Vloc, NULL, &Na->xl)); 2849566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Uloc, NULL, &Na->yl)); 285d751fb92SStefano Zampini } 286ab92ecdeSBarry Smith 287d751fb92SStefano Zampini /* Internally create a scaling vector if roottypes do not match */ 288d751fb92SStefano Zampini if (Na->c) { 289d751fb92SStefano Zampini VecType rt1, rt2; 290d751fb92SStefano Zampini 2919566063dSJacob Faibussowitsch PetscCall(VecGetRootType_Private(Na->work1, &rt1)); 2929566063dSJacob Faibussowitsch PetscCall(VecGetRootType_Private(Na->c, &rt2)); 2939566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(rt1, rt2, &match)); 294d751fb92SStefano Zampini if (!match) { 2959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->c, &c)); 2969566063dSJacob Faibussowitsch PetscCall(VecCopy(Na->c, c)); 2979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 298d751fb92SStefano Zampini Na->c = c; 299d751fb92SStefano Zampini } 300d751fb92SStefano Zampini } 3010575051bSJose E. Roman 302ab92ecdeSBarry Smith (*N)->ops->destroy = MatDestroy_LRC; 303ab92ecdeSBarry Smith (*N)->ops->mult = MatMult_LRC; 304d751fb92SStefano Zampini (*N)->ops->multtranspose = MatMultTranspose_LRC; 305d751fb92SStefano Zampini 306ceb9aaf7SBarry Smith (*N)->assembled = PETSC_TRUE; 30726c5da1bSJose E. Roman (*N)->preallocated = PETSC_TRUE; 308267ca84cSJose E. Roman 3099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatLRCGetMats_C", MatLRCGetMats_LRC)); 3109566063dSJacob Faibussowitsch PetscCall(MatSetUp(*N)); 311ceb9aaf7SBarry Smith PetscFunctionReturn(0); 312ceb9aaf7SBarry Smith } 313