1b9ac7092SAlp Dener /* 2f7bf01afSAlp Dener This provides a thin wrapper around LMVM matrices in order to use their MatLMVMSolve 3b9ac7092SAlp Dener methods as preconditioner applications in KSP solves. 4b9ac7092SAlp Dener */ 5b9ac7092SAlp Dener 6b9ac7092SAlp Dener #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 7dbf6bb8dSprj- #include <petsc/private/matimpl.h> 8b9ac7092SAlp Dener 9b9ac7092SAlp Dener typedef struct { 10b9ac7092SAlp Dener Vec xwork, ywork; 11b9ac7092SAlp Dener IS inactive; 12b9ac7092SAlp Dener Mat B; 13b2d8c577SAlp Dener PetscBool allocated; 14b9ac7092SAlp Dener } PC_LMVM; 15b9ac7092SAlp Dener 16b9ac7092SAlp Dener /*@ 17f1580f4eSBarry Smith PCLMVMSetMatLMVM - Replaces the `MATLMVM` matrix inside the preconditioner with 18b9ac7092SAlp Dener the one provided by the user. 19b9ac7092SAlp Dener 20b9ac7092SAlp Dener Input Parameters: 21f1580f4eSBarry Smith + pc - An `PCLMVM` preconditioner 22f1580f4eSBarry Smith - B - An LMVM-type matrix (`MATLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`) 23b9ac7092SAlp Dener 24b9ac7092SAlp Dener Level: intermediate 25f1580f4eSBarry Smith 26f1580f4eSBarry Smith .seealso: `PCLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`, `PCLMVMGetMatLMVM()` 27b9ac7092SAlp Dener @*/ 28d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMSetMatLMVM(PC pc, Mat B) 29d71ae5a4SJacob Faibussowitsch { 30b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 31b9ac7092SAlp Dener PetscBool same; 32b9ac7092SAlp Dener 33b9ac7092SAlp Dener PetscFunctionBegin; 34b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 35b9ac7092SAlp Dener PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 3728b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 389566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &same)); 3928b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix must be an LMVM-type."); 409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->B)); 419566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 42b9ac7092SAlp Dener ctx->B = B; 433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44b9ac7092SAlp Dener } 45b9ac7092SAlp Dener 46b9ac7092SAlp Dener /*@ 47f1580f4eSBarry Smith PCLMVMGetMatLMVM - Returns a pointer to the underlying `MATLMVM` matrix. 48b9ac7092SAlp Dener 49f1580f4eSBarry Smith Input Parameter: 50f1580f4eSBarry Smith . pc - An `PCLMVM` preconditioner 51b9ac7092SAlp Dener 52f1580f4eSBarry Smith Output Parameter: 53f1580f4eSBarry Smith . B - matrix inside the preconditioner, one of type `MATLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN` 54b9ac7092SAlp Dener 55b9ac7092SAlp Dener Level: intermediate 56f1580f4eSBarry Smith 57f1580f4eSBarry Smith .seealso: `PCLMVM`, `MATLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`, `PCLMVMSetMatLMVM()` 58b9ac7092SAlp Dener @*/ 59d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMGetMatLMVM(PC pc, Mat *B) 60d71ae5a4SJacob Faibussowitsch { 61b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 62b9ac7092SAlp Dener PetscBool same; 63b9ac7092SAlp Dener 64b9ac7092SAlp Dener PetscFunctionBegin; 65b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 669566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 6728b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 68b9ac7092SAlp Dener *B = ctx->B; 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70b9ac7092SAlp Dener } 71b9ac7092SAlp Dener 72b9ac7092SAlp Dener /*@ 73f1580f4eSBarry Smith PCLMVMSetIS - Sets the index sets that reduce the `PC` application. 74b9ac7092SAlp Dener 75b9ac7092SAlp Dener Input Parameters: 76f1580f4eSBarry Smith + pc - An `PCLMVM` preconditioner 77b9ac7092SAlp Dener - inactive - Index set defining the variables removed from the problem 78b9ac7092SAlp Dener 79b9ac7092SAlp Dener Level: intermediate 80b9ac7092SAlp Dener 81*feefa0e1SJacob Faibussowitsch Developer Notes: 82f1580f4eSBarry Smith Need to explain the purpose of this `IS` 83f1580f4eSBarry Smith 84f1580f4eSBarry Smith .seealso: `PCLMVM`, `MatLMVMUpdate()` 85b9ac7092SAlp Dener @*/ 86d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMSetIS(PC pc, IS inactive) 87d71ae5a4SJacob Faibussowitsch { 88b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 89b9ac7092SAlp Dener PetscBool same; 90b9ac7092SAlp Dener 91b9ac7092SAlp Dener PetscFunctionBegin; 92b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 93b2d8c577SAlp Dener PetscValidHeaderSpecific(inactive, IS_CLASSID, 2); 949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 9528b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 969566063dSJacob Faibussowitsch PetscCall(PCLMVMClearIS(pc)); 979566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)inactive)); 98b9ac7092SAlp Dener ctx->inactive = inactive; 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 100b9ac7092SAlp Dener } 101b9ac7092SAlp Dener 102b9ac7092SAlp Dener /*@ 103f1580f4eSBarry Smith PCLMVMClearIS - Removes the inactive variable index set from a `PCLMVM` 104b9ac7092SAlp Dener 105f1580f4eSBarry Smith Input Parameter: 106f1580f4eSBarry Smith . pc - An `PCLMVM` preconditioner 107b9ac7092SAlp Dener 108b9ac7092SAlp Dener Level: intermediate 109b9ac7092SAlp Dener 110f1580f4eSBarry Smith .seealso: `PCLMVM`, `MatLMVMUpdate()` 111b9ac7092SAlp Dener @*/ 112d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMClearIS(PC pc) 113d71ae5a4SJacob Faibussowitsch { 114b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 115b9ac7092SAlp Dener PetscBool same; 116b9ac7092SAlp Dener 117b9ac7092SAlp Dener PetscFunctionBegin; 118b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1199566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 12028b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 12148a46eb9SPierre Jolivet if (ctx->inactive) PetscCall(ISDestroy(&ctx->inactive)); 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123b9ac7092SAlp Dener } 124b9ac7092SAlp Dener 125d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LMVM(PC pc, Vec x, Vec y) 126d71ae5a4SJacob Faibussowitsch { 127b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 128b9ac7092SAlp Dener Vec xsub, ysub; 129b9ac7092SAlp Dener 130b9ac7092SAlp Dener PetscFunctionBegin; 131b9ac7092SAlp Dener if (ctx->inactive) { 1329566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(ctx->xwork)); 1339566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ctx->xwork, ctx->inactive, &xsub)); 1349566063dSJacob Faibussowitsch PetscCall(VecCopy(x, xsub)); 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ctx->xwork, ctx->inactive, &xsub)); 136b9ac7092SAlp Dener } else { 1379566063dSJacob Faibussowitsch PetscCall(VecCopy(x, ctx->xwork)); 138b9ac7092SAlp Dener } 1399566063dSJacob Faibussowitsch PetscCall(MatSolve(ctx->B, ctx->xwork, ctx->ywork)); 140b9ac7092SAlp Dener if (ctx->inactive) { 1419566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ctx->ywork, ctx->inactive, &ysub)); 1429566063dSJacob Faibussowitsch PetscCall(VecCopy(ysub, y)); 1439566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ctx->ywork, ctx->inactive, &ysub)); 144b9ac7092SAlp Dener } else { 1459566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->ywork, y)); 146b9ac7092SAlp Dener } 1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148b9ac7092SAlp Dener } 149b9ac7092SAlp Dener 150d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LMVM(PC pc) 151d71ae5a4SJacob Faibussowitsch { 152b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 153b9ac7092SAlp Dener 154b9ac7092SAlp Dener PetscFunctionBegin; 15548a46eb9SPierre Jolivet if (ctx->xwork) PetscCall(VecDestroy(&ctx->xwork)); 15648a46eb9SPierre Jolivet if (ctx->ywork) PetscCall(VecDestroy(&ctx->ywork)); 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 158b9ac7092SAlp Dener } 159b9ac7092SAlp Dener 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LMVM(PC pc) 161d71ae5a4SJacob Faibussowitsch { 162b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 163b9ac7092SAlp Dener PetscInt n, N; 164b2d8c577SAlp Dener PetscBool allocated; 1650e02354eSStefano Zampini const char *prefix; 166b9ac7092SAlp Dener 167b9ac7092SAlp Dener PetscFunctionBegin; 1683ba16761SJacob Faibussowitsch if (pc->setupcalled) PetscFunctionReturn(PETSC_SUCCESS); 1699566063dSJacob Faibussowitsch PetscCall(MatLMVMIsAllocated(ctx->B, &allocated)); 170b2d8c577SAlp Dener if (!allocated) { 1719566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->mat, &ctx->xwork, &ctx->ywork)); 1729566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->xwork, &n)); 1739566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->xwork, &N)); 1749566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx->B, n, n, N, N)); 1759566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(ctx->B, ctx->xwork, ctx->ywork)); 176b2d8c577SAlp Dener } else { 1779566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx->B, &ctx->xwork, &ctx->ywork)); 178b2d8c577SAlp Dener } 1790e02354eSStefano Zampini PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1800e02354eSStefano Zampini PetscCall(MatSetOptionsPrefix(ctx->B, prefix)); 1810e02354eSStefano Zampini PetscCall(MatAppendOptionsPrefix(ctx->B, "pc_lmvm_")); 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183b9ac7092SAlp Dener } 184b9ac7092SAlp Dener 185d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_LMVM(PC pc, PetscViewer viewer) 186d71ae5a4SJacob Faibussowitsch { 187b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 188b9ac7092SAlp Dener PetscBool iascii; 189b9ac7092SAlp Dener 190b9ac7092SAlp Dener PetscFunctionBegin; 1919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 192dbf6bb8dSprj- if (iascii && ctx->B->assembled) { 1939566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 1949566063dSJacob Faibussowitsch PetscCall(MatView(ctx->B, viewer)); 1959566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 196b9ac7092SAlp Dener } 1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 198b9ac7092SAlp Dener } 199b9ac7092SAlp Dener 200d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_LMVM(PC pc, PetscOptionItems *PetscOptionsObject) 201d71ae5a4SJacob Faibussowitsch { 202b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 2030e02354eSStefano Zampini const char *prefix; 204b9ac7092SAlp Dener 205b9ac7092SAlp Dener PetscFunctionBegin; 2060e02354eSStefano Zampini PetscCall(PCGetOptionsPrefix(pc, &prefix)); 2070e02354eSStefano Zampini PetscCall(MatSetOptionsPrefix(ctx->B, prefix)); 2080e02354eSStefano Zampini PetscCall(MatAppendOptionsPrefix(ctx->B, "pc_lmvm_")); 2099566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ctx->B)); 2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 211b9ac7092SAlp Dener } 212b9ac7092SAlp Dener 213d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LMVM(PC pc) 214d71ae5a4SJacob Faibussowitsch { 215b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 216b9ac7092SAlp Dener 217b9ac7092SAlp Dener PetscFunctionBegin; 21848a46eb9SPierre Jolivet if (ctx->inactive) PetscCall(ISDestroy(&ctx->inactive)); 219b9ac7092SAlp Dener if (pc->setupcalled) { 2209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->xwork)); 2219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->ywork)); 222b9ac7092SAlp Dener } 2239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->B)); 2249566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 226b9ac7092SAlp Dener } 227b9ac7092SAlp Dener 228b9ac7092SAlp Dener /*MC 229f1580f4eSBarry Smith PCLMVM - A a preconditioner constructed from a `MATLMVM` matrix. Options for the 230f1580f4eSBarry Smith underlying `MATLMVM` matrix can be access with the -pc_lmvm_ prefix. 231b9ac7092SAlp Dener 23206dd6b0eSSatish Balay Level: intermediate 23306dd6b0eSSatish Balay 234f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PCLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`, 235db781477SPatrick Sanan `PC`, `MATLMVM`, `PCLMVMUpdate()`, `PCLMVMSetMatLMVM()`, `PCLMVMGetMatLMVM()` 236b9ac7092SAlp Dener M*/ 237d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LMVM(PC pc) 238d71ae5a4SJacob Faibussowitsch { 239b9ac7092SAlp Dener PC_LMVM *ctx; 240b9ac7092SAlp Dener 241b9ac7092SAlp Dener PetscFunctionBegin; 2424dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 243b9ac7092SAlp Dener pc->data = (void *)ctx; 244b9ac7092SAlp Dener 245b9ac7092SAlp Dener pc->ops->reset = PCReset_LMVM; 246b9ac7092SAlp Dener pc->ops->setup = PCSetUp_LMVM; 247b9ac7092SAlp Dener pc->ops->destroy = PCDestroy_LMVM; 248b9ac7092SAlp Dener pc->ops->view = PCView_LMVM; 249b9ac7092SAlp Dener pc->ops->apply = PCApply_LMVM; 250b9ac7092SAlp Dener pc->ops->setfromoptions = PCSetFromOptions_LMVM; 2510a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 2520a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 2530a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 2540a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2550a545947SLisandro Dalcin pc->ops->presolve = NULL; 2560a545947SLisandro Dalcin pc->ops->postsolve = NULL; 257b9ac7092SAlp Dener 2589566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &ctx->B)); 2599566063dSJacob Faibussowitsch PetscCall(MatSetType(ctx->B, MATLMVMBFGS)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->B, (PetscObject)pc, 1)); 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 262b9ac7092SAlp Dener } 263