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 14d751fb92SStefano 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; 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)); 26d751fb92SStefano Zampini if (Na->c) { 279566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->work1,Na->c,Na->work1)); 28d751fb92SStefano Zampini } 29d751fb92SStefano Zampini if (Na->A) { 30d751fb92SStefano Zampini if (transpose) { 319566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A,x,y)); 32d751fb92SStefano Zampini } else { 339566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,x,y)); 34d751fb92SStefano Zampini } 359566063dSJacob Faibussowitsch PetscCall(MatMultAdd(U,Na->work1,y,y)); 36d751fb92SStefano Zampini } else { 379566063dSJacob Faibussowitsch PetscCall(MatMult(U,Na->work1,y)); 38d751fb92SStefano Zampini } 39d751fb92SStefano Zampini } else { 40d751fb92SStefano Zampini Mat Uloc,Vloc; 41d751fb92SStefano Zampini Vec yl,xl; 42d751fb92SStefano Zampini const PetscScalar *w1; 43d751fb92SStefano Zampini PetscScalar *w2; 44d751fb92SStefano Zampini PetscInt nwork; 45d751fb92SStefano Zampini PetscMPIInt mpinwork; 46d751fb92SStefano Zampini 47d751fb92SStefano Zampini xl = transpose ? Na->yl : Na->xl; 48d751fb92SStefano Zampini yl = transpose ? Na->xl : Na->yl; 499566063dSJacob Faibussowitsch PetscCall(VecGetLocalVector(y,yl)); 509566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(U,&Uloc)); 519566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(V,&Vloc)); 520575051bSJose E. Roman 53ab92ecdeSBarry Smith /* multiply the local part of V with the local part of x */ 549566063dSJacob Faibussowitsch PetscCall(VecGetLocalVectorRead(x,xl)); 559566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(Vloc,xl,Na->work1)); 569566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVectorRead(x,xl)); 57ab92ecdeSBarry Smith 5886eed97dSJose E. Roman /* form the sum of all the local multiplies: this is work2 = V'*x = 59ab92ecdeSBarry Smith sum_{all processors} work1 */ 609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Na->work1,&w1)); 619566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(Na->work2,&w2)); 629566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Na->work1,&nwork)); 639566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nwork,&mpinwork)); 649566063dSJacob Faibussowitsch PetscCallMPI(MPIU_Allreduce(w1,w2,mpinwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N))); 659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Na->work1,&w1)); 669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(Na->work2,&w2)); 67ab92ecdeSBarry Smith 68986c4d72SJose E. Roman if (Na->c) { /* work2 = C*work2 */ 699566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->work2,Na->c,Na->work2)); 70986c4d72SJose E. Roman } 71986c4d72SJose E. Roman 7286eed97dSJose E. Roman if (Na->A) { 73d751fb92SStefano Zampini /* form y = A*x or A^t*x */ 74d751fb92SStefano Zampini if (transpose) { 759566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A,x,y)); 76d751fb92SStefano Zampini } else { 779566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,x,y)); 78d751fb92SStefano Zampini } 7986eed97dSJose E. Roman /* multiply-add y = y + U*work2 */ 809566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Uloc,Na->work2,yl,yl)); 8186eed97dSJose E. Roman } else { 8286eed97dSJose E. Roman /* multiply y = U*work2 */ 839566063dSJacob Faibussowitsch PetscCall(MatMult(Uloc,Na->work2,yl)); 8486eed97dSJose E. Roman } 850575051bSJose E. Roman 869566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVector(y,yl)); 87d751fb92SStefano Zampini } 88ceb9aaf7SBarry Smith PetscFunctionReturn(0); 89ceb9aaf7SBarry Smith } 90ceb9aaf7SBarry Smith 91d751fb92SStefano Zampini static PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y) 92d751fb92SStefano Zampini { 93d751fb92SStefano Zampini PetscFunctionBegin; 949566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N,x,y,PETSC_FALSE)); 95d751fb92SStefano Zampini PetscFunctionReturn(0); 96d751fb92SStefano Zampini } 97d751fb92SStefano Zampini 98d751fb92SStefano Zampini static PetscErrorCode MatMultTranspose_LRC(Mat N,Vec x,Vec y) 99d751fb92SStefano Zampini { 100d751fb92SStefano Zampini PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N,x,y,PETSC_TRUE)); 102d751fb92SStefano Zampini PetscFunctionReturn(0); 103d751fb92SStefano Zampini } 104d751fb92SStefano Zampini 105d751fb92SStefano Zampini static PetscErrorCode MatDestroy_LRC(Mat N) 106ceb9aaf7SBarry Smith { 107ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 108ceb9aaf7SBarry Smith 109ceb9aaf7SBarry Smith PetscFunctionBegin; 1109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->U)); 1129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->V)); 1139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 1149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work1)); 1159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work2)); 1169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->xl)); 1179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->yl)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 1199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatLRCGetMats_C",NULL)); 120267ca84cSJose E. Roman PetscFunctionReturn(0); 121267ca84cSJose E. Roman } 122267ca84cSJose E. Roman 123d751fb92SStefano Zampini static PetscErrorCode MatLRCGetMats_LRC(Mat N,Mat *A,Mat *U,Vec *c,Mat *V) 124267ca84cSJose E. Roman { 125267ca84cSJose E. Roman Mat_LRC *Na = (Mat_LRC*)N->data; 126267ca84cSJose E. Roman 127267ca84cSJose E. Roman PetscFunctionBegin; 128267ca84cSJose E. Roman if (A) *A = Na->A; 129267ca84cSJose E. Roman if (U) *U = Na->U; 130267ca84cSJose E. Roman if (c) *c = Na->c; 131267ca84cSJose E. Roman if (V) *V = Na->V; 132267ca84cSJose E. Roman PetscFunctionReturn(0); 133267ca84cSJose E. Roman } 134267ca84cSJose E. Roman 135267ca84cSJose E. Roman /*@ 136267ca84cSJose E. Roman MatLRCGetMats - Returns the constituents of an LRC matrix 137267ca84cSJose E. Roman 138267ca84cSJose E. Roman Collective on Mat 139267ca84cSJose E. Roman 140267ca84cSJose E. Roman Input Parameter: 141267ca84cSJose E. Roman . N - matrix of type LRC 142267ca84cSJose E. Roman 143267ca84cSJose E. Roman Output Parameters: 144267ca84cSJose E. Roman + A - the (sparse) matrix 1456b867d5aSJose E. Roman . U - first dense rectangular (tall and skinny) matrix 1466b867d5aSJose E. Roman . c - a sequential vector containing the diagonal of C 1476b867d5aSJose E. Roman - V - second dense rectangular (tall and skinny) matrix 148267ca84cSJose E. Roman 149267ca84cSJose E. Roman Note: 150267ca84cSJose E. Roman The returned matrices need not be destroyed by the caller. 151267ca84cSJose E. Roman 152267ca84cSJose E. Roman Level: intermediate 153267ca84cSJose E. Roman 154267ca84cSJose E. Roman .seealso: MatCreateLRC() 155267ca84cSJose E. Roman @*/ 156267ca84cSJose E. Roman PetscErrorCode MatLRCGetMats(Mat N,Mat *A,Mat *U,Vec *c,Mat *V) 157267ca84cSJose E. Roman { 158267ca84cSJose E. Roman PetscFunctionBegin; 159*cac4c232SBarry Smith PetscUseMethod(N,"MatLRCGetMats_C",(Mat,Mat*,Mat*,Vec*,Mat*),(N,A,U,c,V)); 160ceb9aaf7SBarry Smith PetscFunctionReturn(0); 161ceb9aaf7SBarry Smith } 162ceb9aaf7SBarry Smith 163ceb9aaf7SBarry Smith /*@ 164986c4d72SJose E. Roman MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V' 165ceb9aaf7SBarry Smith 166ceb9aaf7SBarry Smith Collective on Mat 167ceb9aaf7SBarry Smith 1689d9032efSJose E. Roman Input Parameters: 16986eed97dSJose E. Roman + A - the (sparse) matrix (can be NULL) 170986c4d72SJose E. Roman . U, V - two dense rectangular (tall and skinny) matrices 171d751fb92SStefano Zampini - c - a vector containing the diagonal of C (can be NULL) 172ceb9aaf7SBarry Smith 173ceb9aaf7SBarry Smith Output Parameter: 174986c4d72SJose E. Roman . N - the matrix that represents A + U*C*V' 175ceb9aaf7SBarry Smith 1769d9032efSJose E. Roman Notes: 177986c4d72SJose E. Roman The matrix A + U*C*V' is not formed! Rather the new matrix 178ceb9aaf7SBarry Smith object performs the matrix-vector product by first multiplying by 1799d9032efSJose E. Roman A and then adding the other term. 1809d9032efSJose E. Roman 181986c4d72SJose E. Roman C is a diagonal matrix (represented as a vector) of order k, 182986c4d72SJose E. Roman where k is the number of columns of both U and V. 18386eed97dSJose E. Roman 184986c4d72SJose E. Roman If A is NULL then the new object behaves like a low-rank matrix U*C*V'. 185986c4d72SJose E. Roman 186986c4d72SJose E. Roman Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'. 187986c4d72SJose E. Roman 188986c4d72SJose E. Roman If c is NULL then the low-rank correction is just U*V'. 189d751fb92SStefano Zampini If a sequential c vector is used for a parallel matrix, 190d751fb92SStefano Zampini PETSc assumes that the values of the vector are consistently set across processors. 19186eed97dSJose E. Roman 1929d9032efSJose E. Roman Level: intermediate 193267ca84cSJose E. Roman 194267ca84cSJose E. Roman .seealso: MatLRCGetMats() 195ceb9aaf7SBarry Smith @*/ 196986c4d72SJose E. Roman PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N) 197ceb9aaf7SBarry Smith { 198986c4d72SJose E. Roman PetscBool match; 1999d9032efSJose E. Roman PetscInt m,n,k,m1,n1,k1; 200ab92ecdeSBarry Smith Mat_LRC *Na; 201d751fb92SStefano Zampini Mat Uloc; 202d751fb92SStefano Zampini PetscMPIInt size, csize = 0; 203ceb9aaf7SBarry Smith 204ceb9aaf7SBarry Smith PetscFunctionBegin; 20586eed97dSJose E. Roman if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2069d9032efSJose E. Roman PetscValidHeaderSpecific(U,MAT_CLASSID,2); 207986c4d72SJose E. Roman if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3); 208986c4d72SJose E. Roman if (V) { 209d751fb92SStefano Zampini PetscValidHeaderSpecific(V,MAT_CLASSID,4); 210d751fb92SStefano Zampini PetscCheckSameComm(U,2,V,4); 211d751fb92SStefano Zampini } 212d751fb92SStefano Zampini if (A) PetscCheckSameComm(A,1,U,2); 213d751fb92SStefano Zampini 214d751fb92SStefano Zampini if (!V) V = U; 2159566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"")); 21628b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense, found %s",((PetscObject)U)->type_name); 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"")); 21828b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix V must be of type dense, found %s",((PetscObject)V)->type_name); 2199566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(U->defaultvectype,V->defaultvectype,&match)); 22028b400f6SJacob 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); 221d751fb92SStefano Zampini if (A) { 2229566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(A->defaultvectype,U->defaultvectype,&match)); 22328b400f6SJacob 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); 224986c4d72SJose E. Roman } 2259d9032efSJose E. Roman 2269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)U),&size)); 2279566063dSJacob Faibussowitsch PetscCall(MatGetSize(U,NULL,&k)); 2289566063dSJacob Faibussowitsch PetscCall(MatGetSize(V,NULL,&k1)); 2292c71b3e2SJacob 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); 2309566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(U,&m,NULL)); 2319566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(V,&n,NULL)); 23286eed97dSJose E. Roman if (A) { 2339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m1,&n1)); 234d751fb92SStefano 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); 235d751fb92SStefano 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); 23686eed97dSJose E. Roman } 237986c4d72SJose E. Roman if (c) { 2389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)c),&csize)); 2399566063dSJacob Faibussowitsch PetscCall(VecGetSize(c,&k1)); 2402c71b3e2SJacob 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); 241d751fb92SStefano 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); 242986c4d72SJose E. Roman } 2439d9032efSJose E. Roman 2449566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)U),N)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N,U->defaultvectype)); 2479566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N,MATLRC)); 248d751fb92SStefano Zampini /* Flag matrix as symmetric if A is symmetric and U == V */ 2499566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N,MAT_SYMMETRIC,(PetscBool)((A ? A->symmetric : PETSC_TRUE) && U == V))); 250ceb9aaf7SBarry Smith 2519566063dSJacob Faibussowitsch PetscCall(PetscNewLog(*N,&Na)); 25238f2d2fdSLisandro Dalcin (*N)->data = (void*)Na; 253ceb9aaf7SBarry Smith Na->A = A; 254a639b4dcSJose E. Roman Na->U = U; 255986c4d72SJose E. Roman Na->c = c; 256a639b4dcSJose E. Roman Na->V = V; 257ab92ecdeSBarry Smith 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Na->U)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Na->V)); 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 262ceb9aaf7SBarry Smith 2639566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Na->U,&Uloc)); 2649566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Uloc,&Na->work1,NULL)); 265d751fb92SStefano Zampini if (size != 1) { 266d751fb92SStefano Zampini Mat Vloc; 267d751fb92SStefano Zampini 268d751fb92SStefano Zampini if (Na->c && csize != 1) { /* scatter parallel vector to sequential */ 269d751fb92SStefano Zampini VecScatter sct; 270d751fb92SStefano Zampini 2719566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToAll(Na->c,&sct,&c)); 2729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sct,Na->c,c,INSERT_VALUES,SCATTER_FORWARD)); 2739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sct,Na->c,c,INSERT_VALUES,SCATTER_FORWARD)); 2749566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sct)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 2769566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)*N,(PetscObject)c)); 277d751fb92SStefano Zampini Na->c = c; 278d751fb92SStefano Zampini } 2799566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Na->V,&Vloc)); 2809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->work1,&Na->work2)); 2819566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Vloc,NULL,&Na->xl)); 2829566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Uloc,NULL,&Na->yl)); 283d751fb92SStefano Zampini } 2849566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->work1)); 2859566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->work1)); 2869566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->xl)); 2879566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)*N,(PetscObject)Na->yl)); 288ab92ecdeSBarry Smith 289d751fb92SStefano Zampini /* Internally create a scaling vector if roottypes do not match */ 290d751fb92SStefano Zampini if (Na->c) { 291d751fb92SStefano Zampini VecType rt1,rt2; 292d751fb92SStefano Zampini 2939566063dSJacob Faibussowitsch PetscCall(VecGetRootType_Private(Na->work1,&rt1)); 2949566063dSJacob Faibussowitsch PetscCall(VecGetRootType_Private(Na->c,&rt2)); 2959566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(rt1,rt2,&match)); 296d751fb92SStefano Zampini if (!match) { 2979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->c,&c)); 2989566063dSJacob Faibussowitsch PetscCall(VecCopy(Na->c,c)); 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 3009566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)*N,(PetscObject)c)); 301d751fb92SStefano Zampini Na->c = c; 302d751fb92SStefano Zampini } 303d751fb92SStefano Zampini } 3040575051bSJose E. Roman 305ab92ecdeSBarry Smith (*N)->ops->destroy = MatDestroy_LRC; 306ab92ecdeSBarry Smith (*N)->ops->mult = MatMult_LRC; 307d751fb92SStefano Zampini (*N)->ops->multtranspose = MatMultTranspose_LRC; 308d751fb92SStefano Zampini 309ceb9aaf7SBarry Smith (*N)->assembled = PETSC_TRUE; 31026c5da1bSJose E. Roman (*N)->preallocated = PETSC_TRUE; 311267ca84cSJose E. Roman 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatLRCGetMats_C",MatLRCGetMats_LRC)); 3139566063dSJacob Faibussowitsch PetscCall(MatSetUp(*N)); 314ceb9aaf7SBarry Smith PetscFunctionReturn(0); 315ceb9aaf7SBarry Smith } 316