1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 2ceb9aaf7SBarry Smith 3d751fb92SStefano Zampini PETSC_EXTERN PetscErrorCode VecGetRootType_Private(Vec, VecType *); 4d751fb92SStefano Zampini 5ceb9aaf7SBarry Smith typedef struct { 6986c4d72SJose E. Roman Mat A; /* sparse matrix */ 7986c4d72SJose E. Roman Mat U, V; /* dense tall-skinny matrices */ 8986c4d72SJose E. Roman Vec c; /* sequential vector containing the diagonal of C */ 9986c4d72SJose E. Roman Vec work1, work2; /* sequential vectors that hold partial products */ 100575051bSJose E. Roman Vec xl, yl; /* auxiliary sequential vectors for matmult operation */ 11ab92ecdeSBarry Smith } Mat_LRC; 12ab92ecdeSBarry Smith 13d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_LRC_kernel(Mat N, Vec x, Vec y, PetscBool transpose) 14d71ae5a4SJacob Faibussowitsch { 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)); 61462c564dSBarry Smith PetscCallMPI(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 } 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86ceb9aaf7SBarry Smith } 87ceb9aaf7SBarry Smith 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_LRC(Mat N, Vec x, Vec y) 89d71ae5a4SJacob Faibussowitsch { 90d751fb92SStefano Zampini PetscFunctionBegin; 919566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N, x, y, PETSC_FALSE)); 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 93d751fb92SStefano Zampini } 94d751fb92SStefano Zampini 95d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_LRC(Mat N, Vec x, Vec y) 96d71ae5a4SJacob Faibussowitsch { 97d751fb92SStefano Zampini PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(MatMult_LRC_kernel(N, x, y, PETSC_TRUE)); 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 100d751fb92SStefano Zampini } 101d751fb92SStefano Zampini 102d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_LRC(Mat N) 103d71ae5a4SJacob Faibussowitsch { 104ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC *)N->data; 105ceb9aaf7SBarry Smith 106ceb9aaf7SBarry Smith PetscFunctionBegin; 1079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->U)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->V)); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->c)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work1)); 1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->work2)); 1139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->xl)); 1149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->yl)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCGetMats_C", NULL)); 117bcf2175cSMatthew Knepley PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCSetMats_C", NULL)); 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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; 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 131267ca84cSJose E. Roman } 132267ca84cSJose E. Roman 133bcf2175cSMatthew Knepley static PetscErrorCode MatLRCSetMats_LRC(Mat N, Mat A, Mat U, Vec c, Mat V) 134bcf2175cSMatthew Knepley { 135bcf2175cSMatthew Knepley Mat_LRC *Na = (Mat_LRC *)N->data; 136bcf2175cSMatthew Knepley 137bcf2175cSMatthew Knepley PetscFunctionBegin; 138bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)A)); 139bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)U)); 140bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)V)); 141bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)c)); 142bcf2175cSMatthew Knepley PetscCall(MatDestroy(&Na->A)); 143bcf2175cSMatthew Knepley PetscCall(MatDestroy(&Na->U)); 144bcf2175cSMatthew Knepley PetscCall(MatDestroy(&Na->V)); 145bcf2175cSMatthew Knepley PetscCall(VecDestroy(&Na->c)); 146bcf2175cSMatthew Knepley Na->A = A; 147bcf2175cSMatthew Knepley Na->U = U; 148bcf2175cSMatthew Knepley Na->c = c; 149bcf2175cSMatthew Knepley Na->V = V; 150bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 151bcf2175cSMatthew Knepley } 152bcf2175cSMatthew Knepley 153267ca84cSJose E. Roman /*@ 154267ca84cSJose E. Roman MatLRCGetMats - Returns the constituents of an LRC matrix 155267ca84cSJose E. Roman 156bcf2175cSMatthew Knepley Not collective 157267ca84cSJose E. Roman 158267ca84cSJose E. Roman Input Parameter: 15911a5261eSBarry Smith . N - matrix of type `MATLRC` 160267ca84cSJose E. Roman 161267ca84cSJose E. Roman Output Parameters: 162267ca84cSJose E. Roman + A - the (sparse) matrix 1636b867d5aSJose E. Roman . U - first dense rectangular (tall and skinny) matrix 1646b867d5aSJose E. Roman . c - a sequential vector containing the diagonal of C 1656b867d5aSJose E. Roman - V - second dense rectangular (tall and skinny) matrix 166267ca84cSJose E. Roman 167267ca84cSJose E. Roman Level: intermediate 168267ca84cSJose E. Roman 1692ef1f0ffSBarry Smith Notes: 170bcf2175cSMatthew Knepley The returned matrices should not be destroyed by the caller. 1712ef1f0ffSBarry Smith 1722ef1f0ffSBarry Smith `U`, `c`, `V` may be `NULL` if not needed 1732ef1f0ffSBarry Smith 174bcf2175cSMatthew Knepley .seealso: [](ch_matrices), `MatLRCSetMats()`, `Mat`, `MATLRC`, `MatCreateLRC()` 175267ca84cSJose E. Roman @*/ 176d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLRCGetMats(Mat N, Mat *A, Mat *U, Vec *c, Mat *V) 177d71ae5a4SJacob Faibussowitsch { 178267ca84cSJose E. Roman PetscFunctionBegin; 179cac4c232SBarry Smith PetscUseMethod(N, "MatLRCGetMats_C", (Mat, Mat *, Mat *, Vec *, Mat *), (N, A, U, c, V)); 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 181ceb9aaf7SBarry Smith } 182ceb9aaf7SBarry Smith 183bcf2175cSMatthew Knepley /*@ 184bcf2175cSMatthew Knepley MatLRCSetMats - Sets the constituents of an LRC matrix 185bcf2175cSMatthew Knepley 186bcf2175cSMatthew Knepley Logically collective 187bcf2175cSMatthew Knepley 188bcf2175cSMatthew Knepley Input Parameters: 189bcf2175cSMatthew Knepley + N - matrix of type `MATLRC` 190bcf2175cSMatthew Knepley . A - the (sparse) matrix 191bcf2175cSMatthew Knepley . U - first dense rectangular (tall and skinny) matrix 192bcf2175cSMatthew Knepley . c - a sequential vector containing the diagonal of C 193bcf2175cSMatthew Knepley - V - second dense rectangular (tall and skinny) matrix 194bcf2175cSMatthew Knepley 195bcf2175cSMatthew Knepley Level: intermediate 196bcf2175cSMatthew Knepley 197bcf2175cSMatthew Knepley Note: 198bcf2175cSMatthew Knepley If `V` is `NULL`, then it is assumed to be identical to `U`. 199bcf2175cSMatthew Knepley 200bcf2175cSMatthew Knepley .seealso: [](ch_matrices), `MatLRCGetMats()`, `Mat`, `MATLRC`, `MatCreateLRC()` 201bcf2175cSMatthew Knepley @*/ 202bcf2175cSMatthew Knepley PetscErrorCode MatLRCSetMats(Mat N, Mat A, Mat U, Vec c, Mat V) 203bcf2175cSMatthew Knepley { 204bcf2175cSMatthew Knepley PetscInt k, k1, m, n, m1, n1; 205bcf2175cSMatthew Knepley PetscBool match; 206bcf2175cSMatthew Knepley 207bcf2175cSMatthew Knepley PetscFunctionBegin; 2080d6f747bSJacob Faibussowitsch if (A) PetscValidHeaderSpecific(A, MAT_CLASSID, 2); 2090d6f747bSJacob Faibussowitsch PetscValidHeaderSpecific(U, MAT_CLASSID, 3); 2100d6f747bSJacob Faibussowitsch if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 4); 211bcf2175cSMatthew Knepley if (V) { 2120d6f747bSJacob Faibussowitsch PetscValidHeaderSpecific(V, MAT_CLASSID, 5); 2130d6f747bSJacob Faibussowitsch PetscCheckSameComm(U, 3, V, 5); 214bcf2175cSMatthew Knepley } 2150d6f747bSJacob Faibussowitsch if (A) PetscCheckSameComm(A, 2, U, 3); 216bcf2175cSMatthew Knepley if (!V) V = U; 217bcf2175cSMatthew Knepley PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)U, &match, MATSEQDENSE, MATMPIDENSE, "")); 218bcf2175cSMatthew Knepley PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix U must be of type dense, found %s", ((PetscObject)U)->type_name); 219bcf2175cSMatthew Knepley PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)V, &match, MATSEQDENSE, MATMPIDENSE, "")); 220bcf2175cSMatthew Knepley PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix V must be of type dense, found %s", ((PetscObject)V)->type_name); 221bcf2175cSMatthew Knepley PetscCall(PetscStrcmp(U->defaultvectype, V->defaultvectype, &match)); 222bcf2175cSMatthew Knepley PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_ARG_WRONG, "Matrix U and V must have the same VecType %s != %s", U->defaultvectype, V->defaultvectype); 223bcf2175cSMatthew Knepley if (A) { 224bcf2175cSMatthew Knepley PetscCall(PetscStrcmp(A->defaultvectype, U->defaultvectype, &match)); 225bcf2175cSMatthew Knepley PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_ARG_WRONG, "Matrix A and U must have the same VecType %s != %s", A->defaultvectype, U->defaultvectype); 226bcf2175cSMatthew Knepley } 227bcf2175cSMatthew Knepley PetscCall(MatGetSize(U, NULL, &k)); 228bcf2175cSMatthew Knepley PetscCall(MatGetSize(V, NULL, &k1)); 229bcf2175cSMatthew Knepley 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); 230bcf2175cSMatthew Knepley PetscCall(MatGetLocalSize(U, &m, NULL)); 231bcf2175cSMatthew Knepley PetscCall(MatGetLocalSize(V, &n, NULL)); 232bcf2175cSMatthew Knepley if (A) { 233bcf2175cSMatthew Knepley PetscCall(MatGetLocalSize(A, &m1, &n1)); 234bcf2175cSMatthew Knepley 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); 235bcf2175cSMatthew Knepley 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); 236bcf2175cSMatthew Knepley } 237bcf2175cSMatthew Knepley if (c) { 238bcf2175cSMatthew Knepley PetscMPIInt size, csize; 239bcf2175cSMatthew Knepley 240bcf2175cSMatthew Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)U), &size)); 241bcf2175cSMatthew Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)c), &csize)); 242bcf2175cSMatthew Knepley PetscCall(VecGetSize(c, &k1)); 243bcf2175cSMatthew Knepley 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); 244bcf2175cSMatthew Knepley 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); 245bcf2175cSMatthew Knepley } 246bcf2175cSMatthew Knepley PetscCall(MatSetSizes(N, m, n, PETSC_DECIDE, PETSC_DECIDE)); 247bcf2175cSMatthew Knepley 248bcf2175cSMatthew Knepley PetscUseMethod(N, "MatLRCSetMats_C", (Mat, Mat, Mat, Vec, Mat), (N, A, U, c, V)); 249bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 250bcf2175cSMatthew Knepley } 251bcf2175cSMatthew Knepley 252bcf2175cSMatthew Knepley static PetscErrorCode MatSetUp_LRC(Mat N) 253bcf2175cSMatthew Knepley { 254bcf2175cSMatthew Knepley Mat_LRC *Na = (Mat_LRC *)N->data; 255bcf2175cSMatthew Knepley Mat A = Na->A; 256bcf2175cSMatthew Knepley Mat U = Na->U; 257bcf2175cSMatthew Knepley Mat V = Na->V; 258bcf2175cSMatthew Knepley Vec c = Na->c; 259bcf2175cSMatthew Knepley Mat Uloc; 260bcf2175cSMatthew Knepley PetscMPIInt size, csize = 0; 261*d070fe2eSStefano Zampini PetscBool sym = (PetscBool)(U == V), dummy; 262bcf2175cSMatthew Knepley 263bcf2175cSMatthew Knepley PetscFunctionBegin; 264bcf2175cSMatthew Knepley PetscCall(MatSetVecType(N, U->defaultvectype)); 265bcf2175cSMatthew Knepley // Flag matrix as symmetric if A is symmetric and U == V 266*d070fe2eSStefano Zampini if (A && sym) PetscCall(MatIsSymmetricKnown(A, &dummy, &sym)); 267*d070fe2eSStefano Zampini PetscCall(MatSetOption(N, MAT_SYMMETRIC, sym)); 268bcf2175cSMatthew Knepley PetscCall(MatDenseGetLocalMatrix(Na->U, &Uloc)); 269bcf2175cSMatthew Knepley PetscCall(MatCreateVecs(Uloc, &Na->work1, NULL)); 270bcf2175cSMatthew Knepley 271bcf2175cSMatthew Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)U), &size)); 272bcf2175cSMatthew Knepley if (c) PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)c), &csize)); 273bcf2175cSMatthew Knepley if (size != 1) { 274bcf2175cSMatthew Knepley Mat Vloc; 275bcf2175cSMatthew Knepley 276bcf2175cSMatthew Knepley if (Na->c && csize != 1) { /* scatter parallel vector to sequential */ 277bcf2175cSMatthew Knepley VecScatter sct; 278bcf2175cSMatthew Knepley 279bcf2175cSMatthew Knepley PetscCall(VecScatterCreateToAll(Na->c, &sct, &c)); 280bcf2175cSMatthew Knepley PetscCall(VecScatterBegin(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 281bcf2175cSMatthew Knepley PetscCall(VecScatterEnd(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 282bcf2175cSMatthew Knepley PetscCall(VecScatterDestroy(&sct)); 283bcf2175cSMatthew Knepley PetscCall(VecDestroy(&Na->c)); 284bcf2175cSMatthew Knepley Na->c = c; 285bcf2175cSMatthew Knepley } 286bcf2175cSMatthew Knepley PetscCall(MatDenseGetLocalMatrix(Na->V, &Vloc)); 287bcf2175cSMatthew Knepley PetscCall(VecDuplicate(Na->work1, &Na->work2)); 288bcf2175cSMatthew Knepley PetscCall(MatCreateVecs(Vloc, NULL, &Na->xl)); 289bcf2175cSMatthew Knepley PetscCall(MatCreateVecs(Uloc, NULL, &Na->yl)); 290bcf2175cSMatthew Knepley } 291bcf2175cSMatthew Knepley // Internally create a scaling vector if roottypes do not match 292bcf2175cSMatthew Knepley if (Na->c) { 293bcf2175cSMatthew Knepley VecType rt1, rt2; 294bcf2175cSMatthew Knepley PetscBool match; 295bcf2175cSMatthew Knepley 296bcf2175cSMatthew Knepley PetscCall(VecGetRootType_Private(Na->work1, &rt1)); 297bcf2175cSMatthew Knepley PetscCall(VecGetRootType_Private(Na->c, &rt2)); 298bcf2175cSMatthew Knepley PetscCall(PetscStrcmp(rt1, rt2, &match)); 299bcf2175cSMatthew Knepley if (!match) { 300bcf2175cSMatthew Knepley PetscCall(VecDuplicate(Na->c, &c)); 301bcf2175cSMatthew Knepley PetscCall(VecCopy(Na->c, c)); 302bcf2175cSMatthew Knepley PetscCall(VecDestroy(&Na->c)); 303bcf2175cSMatthew Knepley Na->c = c; 304bcf2175cSMatthew Knepley } 305bcf2175cSMatthew Knepley } 306bcf2175cSMatthew Knepley N->assembled = PETSC_TRUE; 307bcf2175cSMatthew Knepley N->preallocated = PETSC_TRUE; 308bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 309bcf2175cSMatthew Knepley } 310bcf2175cSMatthew Knepley 311bcf2175cSMatthew Knepley PETSC_EXTERN PetscErrorCode MatCreate_LRC(Mat N) 312bcf2175cSMatthew Knepley { 313bcf2175cSMatthew Knepley Mat_LRC *Na; 314bcf2175cSMatthew Knepley 315bcf2175cSMatthew Knepley PetscFunctionBegin; 316bcf2175cSMatthew Knepley PetscCall(PetscObjectChangeTypeName((PetscObject)N, MATLRC)); 317bcf2175cSMatthew Knepley PetscCall(PetscNew(&Na)); 318bcf2175cSMatthew Knepley N->data = (void *)Na; 319bcf2175cSMatthew Knepley N->ops->destroy = MatDestroy_LRC; 320bcf2175cSMatthew Knepley N->ops->setup = MatSetUp_LRC; 321bcf2175cSMatthew Knepley N->ops->mult = MatMult_LRC; 322bcf2175cSMatthew Knepley N->ops->multtranspose = MatMultTranspose_LRC; 323bcf2175cSMatthew Knepley 324bcf2175cSMatthew Knepley PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCGetMats_C", MatLRCGetMats_LRC)); 325bcf2175cSMatthew Knepley PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCSetMats_C", MatLRCSetMats_LRC)); 326bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 327bcf2175cSMatthew Knepley } 328bcf2175cSMatthew Knepley 32911a5261eSBarry Smith /*MC 33011a5261eSBarry Smith MATLRC - "lrc" - a matrix object that behaves like A + U*C*V' 331ceb9aaf7SBarry Smith 33211a5261eSBarry Smith Note: 33311a5261eSBarry Smith The matrix A + U*C*V' is not formed! Rather the matrix object performs the matrix-vector product `MatMult()`, by first multiplying by 33411a5261eSBarry Smith A and then adding the other term. 33511a5261eSBarry Smith 33611a5261eSBarry Smith Level: advanced 33711a5261eSBarry Smith 338bcf2175cSMatthew Knepley .seealso: [](ch_matrices), `Mat`, `MatCreateLRC()`, `MatMult()`, `MatLRCGetMats()`, `MatLRCSetMats()` 33911a5261eSBarry Smith M*/ 34011a5261eSBarry Smith 34111a5261eSBarry Smith /*@ 34211a5261eSBarry Smith MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V' of type `MATLRC` 34311a5261eSBarry Smith 344c3339decSBarry Smith Collective 345ceb9aaf7SBarry Smith 3469d9032efSJose E. Roman Input Parameters: 3472ef1f0ffSBarry Smith + A - the (sparse) matrix (can be `NULL`) 3482ef1f0ffSBarry Smith . U - dense rectangular (tall and skinny) matrix 3492ef1f0ffSBarry Smith . V - dense rectangular (tall and skinny) matrix 3502ef1f0ffSBarry Smith - c - a vector containing the diagonal of C (can be `NULL`) 351ceb9aaf7SBarry Smith 352ceb9aaf7SBarry Smith Output Parameter: 353986c4d72SJose E. Roman . N - the matrix that represents A + U*C*V' 354ceb9aaf7SBarry Smith 3552ef1f0ffSBarry Smith Level: intermediate 3562ef1f0ffSBarry Smith 3579d9032efSJose E. Roman Notes: 358986c4d72SJose E. Roman The matrix A + U*C*V' is not formed! Rather the new matrix 35911a5261eSBarry Smith object performs the matrix-vector product `MatMult()`, by first multiplying by 3609d9032efSJose E. Roman A and then adding the other term. 3619d9032efSJose E. Roman 3622ef1f0ffSBarry Smith `C` is a diagonal matrix (represented as a vector) of order k, 3632ef1f0ffSBarry Smith where k is the number of columns of both `U` and `V`. 36486eed97dSJose E. Roman 3652ef1f0ffSBarry Smith If `A` is `NULL` then the new object behaves like a low-rank matrix U*C*V'. 366986c4d72SJose E. Roman 3672ef1f0ffSBarry Smith Use `V`=`U` (or `V`=`NULL`) for a symmetric low-rank correction, A + U*C*U'. 368986c4d72SJose E. Roman 3692ef1f0ffSBarry Smith If `c` is `NULL` then the low-rank correction is just U*V'. 3702ef1f0ffSBarry Smith If a sequential `c` vector is used for a parallel matrix, 371d751fb92SStefano Zampini PETSc assumes that the values of the vector are consistently set across processors. 37286eed97dSJose E. Roman 3731cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATLRC`, `MatLRCGetMats()` 374ceb9aaf7SBarry Smith @*/ 375d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateLRC(Mat A, Mat U, Vec c, Mat V, Mat *N) 376d71ae5a4SJacob Faibussowitsch { 377ceb9aaf7SBarry Smith PetscFunctionBegin; 3789566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)U), N)); 379bcf2175cSMatthew Knepley PetscCall(MatSetType(*N, MATLRC)); 380bcf2175cSMatthew Knepley PetscCall(MatLRCSetMats(*N, A, U, c, V)); 3819566063dSJacob Faibussowitsch PetscCall(MatSetUp(*N)); 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 383ceb9aaf7SBarry Smith } 384