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 } 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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)); 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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)); 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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)); 118bcf2175cSMatthew Knepley PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCSetMats_C", NULL)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120267ca84cSJose E. Roman } 121267ca84cSJose E. Roman 122d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLRCGetMats_LRC(Mat N, Mat *A, Mat *U, Vec *c, Mat *V) 123d71ae5a4SJacob Faibussowitsch { 124267ca84cSJose E. Roman Mat_LRC *Na = (Mat_LRC *)N->data; 125267ca84cSJose E. Roman 126267ca84cSJose E. Roman PetscFunctionBegin; 127267ca84cSJose E. Roman if (A) *A = Na->A; 128267ca84cSJose E. Roman if (U) *U = Na->U; 129267ca84cSJose E. Roman if (c) *c = Na->c; 130267ca84cSJose E. Roman if (V) *V = Na->V; 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132267ca84cSJose E. Roman } 133267ca84cSJose E. Roman 134bcf2175cSMatthew Knepley static PetscErrorCode MatLRCSetMats_LRC(Mat N, Mat A, Mat U, Vec c, Mat V) 135bcf2175cSMatthew Knepley { 136bcf2175cSMatthew Knepley Mat_LRC *Na = (Mat_LRC *)N->data; 137bcf2175cSMatthew Knepley 138bcf2175cSMatthew Knepley PetscFunctionBegin; 139bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)A)); 140bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)U)); 141bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)V)); 142bcf2175cSMatthew Knepley PetscCall(PetscObjectReference((PetscObject)c)); 143bcf2175cSMatthew Knepley PetscCall(MatDestroy(&Na->A)); 144bcf2175cSMatthew Knepley PetscCall(MatDestroy(&Na->U)); 145bcf2175cSMatthew Knepley PetscCall(MatDestroy(&Na->V)); 146bcf2175cSMatthew Knepley PetscCall(VecDestroy(&Na->c)); 147bcf2175cSMatthew Knepley Na->A = A; 148bcf2175cSMatthew Knepley Na->U = U; 149bcf2175cSMatthew Knepley Na->c = c; 150bcf2175cSMatthew Knepley Na->V = V; 151bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 152bcf2175cSMatthew Knepley } 153bcf2175cSMatthew Knepley 154267ca84cSJose E. Roman /*@ 155267ca84cSJose E. Roman MatLRCGetMats - Returns the constituents of an LRC matrix 156267ca84cSJose E. Roman 157bcf2175cSMatthew Knepley Not collective 158267ca84cSJose E. Roman 159267ca84cSJose E. Roman Input Parameter: 16011a5261eSBarry Smith . N - matrix of type `MATLRC` 161267ca84cSJose E. Roman 162267ca84cSJose E. Roman Output Parameters: 163267ca84cSJose E. Roman + A - the (sparse) matrix 1646b867d5aSJose E. Roman . U - first dense rectangular (tall and skinny) matrix 1656b867d5aSJose E. Roman . c - a sequential vector containing the diagonal of C 1666b867d5aSJose E. Roman - V - second dense rectangular (tall and skinny) matrix 167267ca84cSJose E. Roman 168267ca84cSJose E. Roman Level: intermediate 169267ca84cSJose E. Roman 1702ef1f0ffSBarry Smith Notes: 171bcf2175cSMatthew Knepley The returned matrices should not be destroyed by the caller. 1722ef1f0ffSBarry Smith 1732ef1f0ffSBarry Smith `U`, `c`, `V` may be `NULL` if not needed 1742ef1f0ffSBarry Smith 175bcf2175cSMatthew Knepley .seealso: [](ch_matrices), `MatLRCSetMats()`, `Mat`, `MATLRC`, `MatCreateLRC()` 176267ca84cSJose E. Roman @*/ 177d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLRCGetMats(Mat N, Mat *A, Mat *U, Vec *c, Mat *V) 178d71ae5a4SJacob Faibussowitsch { 179267ca84cSJose E. Roman PetscFunctionBegin; 180cac4c232SBarry Smith PetscUseMethod(N, "MatLRCGetMats_C", (Mat, Mat *, Mat *, Vec *, Mat *), (N, A, U, c, V)); 1813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 182ceb9aaf7SBarry Smith } 183ceb9aaf7SBarry Smith 184bcf2175cSMatthew Knepley /*@ 185bcf2175cSMatthew Knepley MatLRCSetMats - Sets the constituents of an LRC matrix 186bcf2175cSMatthew Knepley 187bcf2175cSMatthew Knepley Logically collective 188bcf2175cSMatthew Knepley 189bcf2175cSMatthew Knepley Input Parameters: 190bcf2175cSMatthew Knepley + N - matrix of type `MATLRC` 191bcf2175cSMatthew Knepley . A - the (sparse) matrix 192bcf2175cSMatthew Knepley . U - first dense rectangular (tall and skinny) matrix 193bcf2175cSMatthew Knepley . c - a sequential vector containing the diagonal of C 194bcf2175cSMatthew Knepley - V - second dense rectangular (tall and skinny) matrix 195bcf2175cSMatthew Knepley 196bcf2175cSMatthew Knepley Level: intermediate 197bcf2175cSMatthew Knepley 198bcf2175cSMatthew Knepley Note: 199bcf2175cSMatthew Knepley If `V` is `NULL`, then it is assumed to be identical to `U`. 200bcf2175cSMatthew Knepley 201bcf2175cSMatthew Knepley .seealso: [](ch_matrices), `MatLRCGetMats()`, `Mat`, `MATLRC`, `MatCreateLRC()` 202bcf2175cSMatthew Knepley @*/ 203bcf2175cSMatthew Knepley PetscErrorCode MatLRCSetMats(Mat N, Mat A, Mat U, Vec c, Mat V) 204bcf2175cSMatthew Knepley { 205bcf2175cSMatthew Knepley PetscInt k, k1, m, n, m1, n1; 206bcf2175cSMatthew Knepley PetscBool match; 207bcf2175cSMatthew Knepley 208bcf2175cSMatthew Knepley PetscFunctionBegin; 209*0d6f747bSJacob Faibussowitsch if (A) PetscValidHeaderSpecific(A, MAT_CLASSID, 2); 210*0d6f747bSJacob Faibussowitsch PetscValidHeaderSpecific(U, MAT_CLASSID, 3); 211*0d6f747bSJacob Faibussowitsch if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 4); 212bcf2175cSMatthew Knepley if (V) { 213*0d6f747bSJacob Faibussowitsch PetscValidHeaderSpecific(V, MAT_CLASSID, 5); 214*0d6f747bSJacob Faibussowitsch PetscCheckSameComm(U, 3, V, 5); 215bcf2175cSMatthew Knepley } 216*0d6f747bSJacob Faibussowitsch if (A) PetscCheckSameComm(A, 2, U, 3); 217bcf2175cSMatthew Knepley if (!V) V = U; 218bcf2175cSMatthew Knepley PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)U, &match, MATSEQDENSE, MATMPIDENSE, "")); 219bcf2175cSMatthew Knepley PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix U must be of type dense, found %s", ((PetscObject)U)->type_name); 220bcf2175cSMatthew Knepley PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)V, &match, MATSEQDENSE, MATMPIDENSE, "")); 221bcf2175cSMatthew Knepley PetscCheck(match, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Matrix V must be of type dense, found %s", ((PetscObject)V)->type_name); 222bcf2175cSMatthew Knepley PetscCall(PetscStrcmp(U->defaultvectype, V->defaultvectype, &match)); 223bcf2175cSMatthew 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); 224bcf2175cSMatthew Knepley if (A) { 225bcf2175cSMatthew Knepley PetscCall(PetscStrcmp(A->defaultvectype, U->defaultvectype, &match)); 226bcf2175cSMatthew 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); 227bcf2175cSMatthew Knepley } 228bcf2175cSMatthew Knepley PetscCall(MatGetSize(U, NULL, &k)); 229bcf2175cSMatthew Knepley PetscCall(MatGetSize(V, NULL, &k1)); 230bcf2175cSMatthew 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); 231bcf2175cSMatthew Knepley PetscCall(MatGetLocalSize(U, &m, NULL)); 232bcf2175cSMatthew Knepley PetscCall(MatGetLocalSize(V, &n, NULL)); 233bcf2175cSMatthew Knepley if (A) { 234bcf2175cSMatthew Knepley PetscCall(MatGetLocalSize(A, &m1, &n1)); 235bcf2175cSMatthew 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); 236bcf2175cSMatthew 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); 237bcf2175cSMatthew Knepley } 238bcf2175cSMatthew Knepley if (c) { 239bcf2175cSMatthew Knepley PetscMPIInt size, csize; 240bcf2175cSMatthew Knepley 241bcf2175cSMatthew Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)U), &size)); 242bcf2175cSMatthew Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)c), &csize)); 243bcf2175cSMatthew Knepley PetscCall(VecGetSize(c, &k1)); 244bcf2175cSMatthew 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); 245bcf2175cSMatthew 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); 246bcf2175cSMatthew Knepley } 247bcf2175cSMatthew Knepley PetscCall(MatSetSizes(N, m, n, PETSC_DECIDE, PETSC_DECIDE)); 248bcf2175cSMatthew Knepley 249bcf2175cSMatthew Knepley PetscUseMethod(N, "MatLRCSetMats_C", (Mat, Mat, Mat, Vec, Mat), (N, A, U, c, V)); 250bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 251bcf2175cSMatthew Knepley } 252bcf2175cSMatthew Knepley 253bcf2175cSMatthew Knepley static PetscErrorCode MatSetUp_LRC(Mat N) 254bcf2175cSMatthew Knepley { 255bcf2175cSMatthew Knepley Mat_LRC *Na = (Mat_LRC *)N->data; 256bcf2175cSMatthew Knepley Mat A = Na->A; 257bcf2175cSMatthew Knepley Mat U = Na->U; 258bcf2175cSMatthew Knepley Mat V = Na->V; 259bcf2175cSMatthew Knepley Vec c = Na->c; 260bcf2175cSMatthew Knepley Mat Uloc; 261bcf2175cSMatthew Knepley PetscMPIInt size, csize = 0; 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 266bcf2175cSMatthew Knepley PetscCall(MatSetOption(N, MAT_SYMMETRIC, (PetscBool)((A ? A->symmetric == PETSC_BOOL3_TRUE : PETSC_TRUE) && U == V))); 267bcf2175cSMatthew Knepley PetscCall(MatDenseGetLocalMatrix(Na->U, &Uloc)); 268bcf2175cSMatthew Knepley PetscCall(MatCreateVecs(Uloc, &Na->work1, NULL)); 269bcf2175cSMatthew Knepley 270bcf2175cSMatthew Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)U), &size)); 271bcf2175cSMatthew Knepley if (c) PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)c), &csize)); 272bcf2175cSMatthew Knepley if (size != 1) { 273bcf2175cSMatthew Knepley Mat Vloc; 274bcf2175cSMatthew Knepley 275bcf2175cSMatthew Knepley if (Na->c && csize != 1) { /* scatter parallel vector to sequential */ 276bcf2175cSMatthew Knepley VecScatter sct; 277bcf2175cSMatthew Knepley 278bcf2175cSMatthew Knepley PetscCall(VecScatterCreateToAll(Na->c, &sct, &c)); 279bcf2175cSMatthew Knepley PetscCall(VecScatterBegin(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 280bcf2175cSMatthew Knepley PetscCall(VecScatterEnd(sct, Na->c, c, INSERT_VALUES, SCATTER_FORWARD)); 281bcf2175cSMatthew Knepley PetscCall(VecScatterDestroy(&sct)); 282bcf2175cSMatthew Knepley PetscCall(VecDestroy(&Na->c)); 283bcf2175cSMatthew Knepley Na->c = c; 284bcf2175cSMatthew Knepley } 285bcf2175cSMatthew Knepley PetscCall(MatDenseGetLocalMatrix(Na->V, &Vloc)); 286bcf2175cSMatthew Knepley PetscCall(VecDuplicate(Na->work1, &Na->work2)); 287bcf2175cSMatthew Knepley PetscCall(MatCreateVecs(Vloc, NULL, &Na->xl)); 288bcf2175cSMatthew Knepley PetscCall(MatCreateVecs(Uloc, NULL, &Na->yl)); 289bcf2175cSMatthew Knepley } 290bcf2175cSMatthew Knepley // Internally create a scaling vector if roottypes do not match 291bcf2175cSMatthew Knepley if (Na->c) { 292bcf2175cSMatthew Knepley VecType rt1, rt2; 293bcf2175cSMatthew Knepley PetscBool match; 294bcf2175cSMatthew Knepley 295bcf2175cSMatthew Knepley PetscCall(VecGetRootType_Private(Na->work1, &rt1)); 296bcf2175cSMatthew Knepley PetscCall(VecGetRootType_Private(Na->c, &rt2)); 297bcf2175cSMatthew Knepley PetscCall(PetscStrcmp(rt1, rt2, &match)); 298bcf2175cSMatthew Knepley if (!match) { 299bcf2175cSMatthew Knepley PetscCall(VecDuplicate(Na->c, &c)); 300bcf2175cSMatthew Knepley PetscCall(VecCopy(Na->c, c)); 301bcf2175cSMatthew Knepley PetscCall(VecDestroy(&Na->c)); 302bcf2175cSMatthew Knepley Na->c = c; 303bcf2175cSMatthew Knepley } 304bcf2175cSMatthew Knepley } 305bcf2175cSMatthew Knepley N->assembled = PETSC_TRUE; 306bcf2175cSMatthew Knepley N->preallocated = PETSC_TRUE; 307bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 308bcf2175cSMatthew Knepley } 309bcf2175cSMatthew Knepley 310bcf2175cSMatthew Knepley PETSC_EXTERN PetscErrorCode MatCreate_LRC(Mat N) 311bcf2175cSMatthew Knepley { 312bcf2175cSMatthew Knepley Mat_LRC *Na; 313bcf2175cSMatthew Knepley 314bcf2175cSMatthew Knepley PetscFunctionBegin; 315bcf2175cSMatthew Knepley PetscCall(PetscObjectChangeTypeName((PetscObject)N, MATLRC)); 316bcf2175cSMatthew Knepley PetscCall(PetscNew(&Na)); 317bcf2175cSMatthew Knepley N->data = (void *)Na; 318bcf2175cSMatthew Knepley N->ops->destroy = MatDestroy_LRC; 319bcf2175cSMatthew Knepley N->ops->setup = MatSetUp_LRC; 320bcf2175cSMatthew Knepley N->ops->mult = MatMult_LRC; 321bcf2175cSMatthew Knepley N->ops->multtranspose = MatMultTranspose_LRC; 322bcf2175cSMatthew Knepley 323bcf2175cSMatthew Knepley PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCGetMats_C", MatLRCGetMats_LRC)); 324bcf2175cSMatthew Knepley PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatLRCSetMats_C", MatLRCSetMats_LRC)); 325bcf2175cSMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 326bcf2175cSMatthew Knepley } 327bcf2175cSMatthew Knepley 32811a5261eSBarry Smith /*MC 32911a5261eSBarry Smith MATLRC - "lrc" - a matrix object that behaves like A + U*C*V' 330ceb9aaf7SBarry Smith 33111a5261eSBarry Smith Note: 33211a5261eSBarry Smith The matrix A + U*C*V' is not formed! Rather the matrix object performs the matrix-vector product `MatMult()`, by first multiplying by 33311a5261eSBarry Smith A and then adding the other term. 33411a5261eSBarry Smith 33511a5261eSBarry Smith Level: advanced 33611a5261eSBarry Smith 337bcf2175cSMatthew Knepley .seealso: [](ch_matrices), `Mat`, `MatCreateLRC()`, `MatMult()`, `MatLRCGetMats()`, `MatLRCSetMats()` 33811a5261eSBarry Smith M*/ 33911a5261eSBarry Smith 34011a5261eSBarry Smith /*@ 34111a5261eSBarry Smith MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V' of type `MATLRC` 34211a5261eSBarry Smith 343c3339decSBarry Smith Collective 344ceb9aaf7SBarry Smith 3459d9032efSJose E. Roman Input Parameters: 3462ef1f0ffSBarry Smith + A - the (sparse) matrix (can be `NULL`) 3472ef1f0ffSBarry Smith . U - dense rectangular (tall and skinny) matrix 3482ef1f0ffSBarry Smith . V - dense rectangular (tall and skinny) matrix 3492ef1f0ffSBarry Smith - c - a vector containing the diagonal of C (can be `NULL`) 350ceb9aaf7SBarry Smith 351ceb9aaf7SBarry Smith Output Parameter: 352986c4d72SJose E. Roman . N - the matrix that represents A + U*C*V' 353ceb9aaf7SBarry Smith 3542ef1f0ffSBarry Smith Level: intermediate 3552ef1f0ffSBarry Smith 3569d9032efSJose E. Roman Notes: 357986c4d72SJose E. Roman The matrix A + U*C*V' is not formed! Rather the new matrix 35811a5261eSBarry Smith object performs the matrix-vector product `MatMult()`, by first multiplying by 3599d9032efSJose E. Roman A and then adding the other term. 3609d9032efSJose E. Roman 3612ef1f0ffSBarry Smith `C` is a diagonal matrix (represented as a vector) of order k, 3622ef1f0ffSBarry Smith where k is the number of columns of both `U` and `V`. 36386eed97dSJose E. Roman 3642ef1f0ffSBarry Smith If `A` is `NULL` then the new object behaves like a low-rank matrix U*C*V'. 365986c4d72SJose E. Roman 3662ef1f0ffSBarry Smith Use `V`=`U` (or `V`=`NULL`) for a symmetric low-rank correction, A + U*C*U'. 367986c4d72SJose E. Roman 3682ef1f0ffSBarry Smith If `c` is `NULL` then the low-rank correction is just U*V'. 3692ef1f0ffSBarry Smith If a sequential `c` vector is used for a parallel matrix, 370d751fb92SStefano Zampini PETSc assumes that the values of the vector are consistently set across processors. 37186eed97dSJose E. Roman 3721cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATLRC`, `MatLRCGetMats()` 373ceb9aaf7SBarry Smith @*/ 374d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateLRC(Mat A, Mat U, Vec c, Mat V, Mat *N) 375d71ae5a4SJacob Faibussowitsch { 376ceb9aaf7SBarry Smith PetscFunctionBegin; 3779566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)U), N)); 378bcf2175cSMatthew Knepley PetscCall(MatSetType(*N, MATLRC)); 379bcf2175cSMatthew Knepley PetscCall(MatLRCSetMats(*N, A, U, c, V)); 3809566063dSJacob Faibussowitsch PetscCall(MatSetUp(*N)); 3813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 382ceb9aaf7SBarry Smith } 383