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 /*@ 17*f1580f4eSBarry Smith PCLMVMSetMatLMVM - Replaces the `MATLMVM` matrix inside the preconditioner with 18b9ac7092SAlp Dener the one provided by the user. 19b9ac7092SAlp Dener 20b9ac7092SAlp Dener Input Parameters: 21*f1580f4eSBarry Smith + pc - An `PCLMVM` preconditioner 22*f1580f4eSBarry Smith - B - An LMVM-type matrix (`MATLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`) 23b9ac7092SAlp Dener 24b9ac7092SAlp Dener Level: intermediate 25*f1580f4eSBarry Smith 26*f1580f4eSBarry Smith .seealso: `PCLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`, `PCLMVMGetMatLMVM()` 27b9ac7092SAlp Dener @*/ 289371c9d4SSatish Balay PetscErrorCode PCLMVMSetMatLMVM(PC pc, Mat B) { 29b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 30b9ac7092SAlp Dener PetscBool same; 31b9ac7092SAlp Dener 32b9ac7092SAlp Dener PetscFunctionBegin; 33b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 34b9ac7092SAlp Dener PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 3628b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 379566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &same)); 3828b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix must be an LMVM-type."); 399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->B)); 409566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 41b9ac7092SAlp Dener ctx->B = B; 42b9ac7092SAlp Dener PetscFunctionReturn(0); 43b9ac7092SAlp Dener } 44b9ac7092SAlp Dener 45b9ac7092SAlp Dener /*@ 46*f1580f4eSBarry Smith PCLMVMGetMatLMVM - Returns a pointer to the underlying `MATLMVM` matrix. 47b9ac7092SAlp Dener 48*f1580f4eSBarry Smith Input Parameter: 49*f1580f4eSBarry Smith . pc - An `PCLMVM` preconditioner 50b9ac7092SAlp Dener 51*f1580f4eSBarry Smith Output Parameter: 52*f1580f4eSBarry Smith . B - matrix inside the preconditioner, one of type `MATLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN` 53b9ac7092SAlp Dener 54b9ac7092SAlp Dener Level: intermediate 55*f1580f4eSBarry Smith 56*f1580f4eSBarry Smith .seealso: `PCLMVM`, `MATLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`, `PCLMVMSetMatLMVM()` 57b9ac7092SAlp Dener @*/ 589371c9d4SSatish Balay PetscErrorCode PCLMVMGetMatLMVM(PC pc, Mat *B) { 59b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 60b9ac7092SAlp Dener PetscBool same; 61b9ac7092SAlp Dener 62b9ac7092SAlp Dener PetscFunctionBegin; 63b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 6528b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 66b9ac7092SAlp Dener *B = ctx->B; 67b9ac7092SAlp Dener PetscFunctionReturn(0); 68b9ac7092SAlp Dener } 69b9ac7092SAlp Dener 70b9ac7092SAlp Dener /*@ 71*f1580f4eSBarry Smith PCLMVMSetIS - Sets the index sets that reduce the `PC` application. 72b9ac7092SAlp Dener 73b9ac7092SAlp Dener Input Parameters: 74*f1580f4eSBarry Smith + pc - An `PCLMVM` preconditioner 75b9ac7092SAlp Dener - inactive - Index set defining the variables removed from the problem 76b9ac7092SAlp Dener 77b9ac7092SAlp Dener Level: intermediate 78b9ac7092SAlp Dener 79*f1580f4eSBarry Smith Developer Note: 80*f1580f4eSBarry Smith Need to explain the purpose of this `IS` 81*f1580f4eSBarry Smith 82*f1580f4eSBarry Smith .seealso: `PCLMVM`, `MatLMVMUpdate()` 83b9ac7092SAlp Dener @*/ 849371c9d4SSatish Balay PetscErrorCode PCLMVMSetIS(PC pc, IS inactive) { 85b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 86b9ac7092SAlp Dener PetscBool same; 87b9ac7092SAlp Dener 88b9ac7092SAlp Dener PetscFunctionBegin; 89b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 90b2d8c577SAlp Dener PetscValidHeaderSpecific(inactive, IS_CLASSID, 2); 919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 9228b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 939566063dSJacob Faibussowitsch PetscCall(PCLMVMClearIS(pc)); 949566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)inactive)); 95b9ac7092SAlp Dener ctx->inactive = inactive; 96b9ac7092SAlp Dener PetscFunctionReturn(0); 97b9ac7092SAlp Dener } 98b9ac7092SAlp Dener 99b9ac7092SAlp Dener /*@ 100*f1580f4eSBarry Smith PCLMVMClearIS - Removes the inactive variable index set from a `PCLMVM` 101b9ac7092SAlp Dener 102*f1580f4eSBarry Smith Input Parameter: 103*f1580f4eSBarry Smith . pc - An `PCLMVM` preconditioner 104b9ac7092SAlp Dener 105b9ac7092SAlp Dener Level: intermediate 106b9ac7092SAlp Dener 107*f1580f4eSBarry Smith .seealso: `PCLMVM`, `MatLMVMUpdate()` 108b9ac7092SAlp Dener @*/ 1099371c9d4SSatish Balay PetscErrorCode PCLMVMClearIS(PC pc) { 110b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 111b9ac7092SAlp Dener PetscBool same; 112b9ac7092SAlp Dener 113b9ac7092SAlp Dener PetscFunctionBegin; 114b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 11628b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 11748a46eb9SPierre Jolivet if (ctx->inactive) PetscCall(ISDestroy(&ctx->inactive)); 118b9ac7092SAlp Dener PetscFunctionReturn(0); 119b9ac7092SAlp Dener } 120b9ac7092SAlp Dener 1219371c9d4SSatish Balay static PetscErrorCode PCApply_LMVM(PC pc, Vec x, Vec y) { 122b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 123b9ac7092SAlp Dener Vec xsub, ysub; 124b9ac7092SAlp Dener 125b9ac7092SAlp Dener PetscFunctionBegin; 126b9ac7092SAlp Dener if (ctx->inactive) { 1279566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(ctx->xwork)); 1289566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ctx->xwork, ctx->inactive, &xsub)); 1299566063dSJacob Faibussowitsch PetscCall(VecCopy(x, xsub)); 1309566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ctx->xwork, ctx->inactive, &xsub)); 131b9ac7092SAlp Dener } else { 1329566063dSJacob Faibussowitsch PetscCall(VecCopy(x, ctx->xwork)); 133b9ac7092SAlp Dener } 1349566063dSJacob Faibussowitsch PetscCall(MatSolve(ctx->B, ctx->xwork, ctx->ywork)); 135b9ac7092SAlp Dener if (ctx->inactive) { 1369566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ctx->ywork, ctx->inactive, &ysub)); 1379566063dSJacob Faibussowitsch PetscCall(VecCopy(ysub, y)); 1389566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ctx->ywork, ctx->inactive, &ysub)); 139b9ac7092SAlp Dener } else { 1409566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->ywork, y)); 141b9ac7092SAlp Dener } 142b9ac7092SAlp Dener PetscFunctionReturn(0); 143b9ac7092SAlp Dener } 144b9ac7092SAlp Dener 1459371c9d4SSatish Balay static PetscErrorCode PCReset_LMVM(PC pc) { 146b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 147b9ac7092SAlp Dener 148b9ac7092SAlp Dener PetscFunctionBegin; 14948a46eb9SPierre Jolivet if (ctx->xwork) PetscCall(VecDestroy(&ctx->xwork)); 15048a46eb9SPierre Jolivet if (ctx->ywork) PetscCall(VecDestroy(&ctx->ywork)); 151b9ac7092SAlp Dener PetscFunctionReturn(0); 152b9ac7092SAlp Dener } 153b9ac7092SAlp Dener 1549371c9d4SSatish Balay static PetscErrorCode PCSetUp_LMVM(PC pc) { 155b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 156b9ac7092SAlp Dener PetscInt n, N; 157b2d8c577SAlp Dener PetscBool allocated; 1580e02354eSStefano Zampini const char *prefix; 159b9ac7092SAlp Dener 160b9ac7092SAlp Dener PetscFunctionBegin; 1617ce701aeSStefano Zampini if (pc->setupcalled) PetscFunctionReturn(0); 1629566063dSJacob Faibussowitsch PetscCall(MatLMVMIsAllocated(ctx->B, &allocated)); 163b2d8c577SAlp Dener if (!allocated) { 1649566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->mat, &ctx->xwork, &ctx->ywork)); 1659566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->xwork, &n)); 1669566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->xwork, &N)); 1679566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx->B, n, n, N, N)); 1689566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(ctx->B, ctx->xwork, ctx->ywork)); 169b2d8c577SAlp Dener } else { 1709566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx->B, &ctx->xwork, &ctx->ywork)); 171b2d8c577SAlp Dener } 1720e02354eSStefano Zampini PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1730e02354eSStefano Zampini PetscCall(MatSetOptionsPrefix(ctx->B, prefix)); 1740e02354eSStefano Zampini PetscCall(MatAppendOptionsPrefix(ctx->B, "pc_lmvm_")); 175b9ac7092SAlp Dener PetscFunctionReturn(0); 176b9ac7092SAlp Dener } 177b9ac7092SAlp Dener 1789371c9d4SSatish Balay static PetscErrorCode PCView_LMVM(PC pc, PetscViewer viewer) { 179b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 180b9ac7092SAlp Dener PetscBool iascii; 181b9ac7092SAlp Dener 182b9ac7092SAlp Dener PetscFunctionBegin; 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 184dbf6bb8dSprj- if (iascii && ctx->B->assembled) { 1859566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 1869566063dSJacob Faibussowitsch PetscCall(MatView(ctx->B, viewer)); 1879566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 188b9ac7092SAlp Dener } 189b9ac7092SAlp Dener PetscFunctionReturn(0); 190b9ac7092SAlp Dener } 191b9ac7092SAlp Dener 1929371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_LMVM(PC pc, PetscOptionItems *PetscOptionsObject) { 193b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 1940e02354eSStefano Zampini const char *prefix; 195b9ac7092SAlp Dener 196b9ac7092SAlp Dener PetscFunctionBegin; 1970e02354eSStefano Zampini PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1980e02354eSStefano Zampini PetscCall(MatSetOptionsPrefix(ctx->B, prefix)); 1990e02354eSStefano Zampini PetscCall(MatAppendOptionsPrefix(ctx->B, "pc_lmvm_")); 2009566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ctx->B)); 201b9ac7092SAlp Dener PetscFunctionReturn(0); 202b9ac7092SAlp Dener } 203b9ac7092SAlp Dener 2049371c9d4SSatish Balay static PetscErrorCode PCDestroy_LMVM(PC pc) { 205b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 206b9ac7092SAlp Dener 207b9ac7092SAlp Dener PetscFunctionBegin; 20848a46eb9SPierre Jolivet if (ctx->inactive) PetscCall(ISDestroy(&ctx->inactive)); 209b9ac7092SAlp Dener if (pc->setupcalled) { 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->xwork)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->ywork)); 212b9ac7092SAlp Dener } 2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->B)); 2149566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 215b9ac7092SAlp Dener PetscFunctionReturn(0); 216b9ac7092SAlp Dener } 217b9ac7092SAlp Dener 218b9ac7092SAlp Dener /*MC 219*f1580f4eSBarry Smith PCLMVM - A a preconditioner constructed from a `MATLMVM` matrix. Options for the 220*f1580f4eSBarry Smith underlying `MATLMVM` matrix can be access with the -pc_lmvm_ prefix. 221b9ac7092SAlp Dener 22206dd6b0eSSatish Balay Level: intermediate 22306dd6b0eSSatish Balay 224*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PCLMVM`, `MATLDFP`, `MATLBFGS`, `MATLSR1`, `MATLBRDN`, `MATLMBRDN`, `MATLSBRDN`, 225db781477SPatrick Sanan `PC`, `MATLMVM`, `PCLMVMUpdate()`, `PCLMVMSetMatLMVM()`, `PCLMVMGetMatLMVM()` 226b9ac7092SAlp Dener M*/ 2279371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_LMVM(PC pc) { 228b9ac7092SAlp Dener PC_LMVM *ctx; 229b9ac7092SAlp Dener 230b9ac7092SAlp Dener PetscFunctionBegin; 2319566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &ctx)); 232b9ac7092SAlp Dener pc->data = (void *)ctx; 233b9ac7092SAlp Dener 234b9ac7092SAlp Dener pc->ops->reset = PCReset_LMVM; 235b9ac7092SAlp Dener pc->ops->setup = PCSetUp_LMVM; 236b9ac7092SAlp Dener pc->ops->destroy = PCDestroy_LMVM; 237b9ac7092SAlp Dener pc->ops->view = PCView_LMVM; 238b9ac7092SAlp Dener pc->ops->apply = PCApply_LMVM; 239b9ac7092SAlp Dener pc->ops->setfromoptions = PCSetFromOptions_LMVM; 2400a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 2410a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 2420a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 2430a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2440a545947SLisandro Dalcin pc->ops->presolve = NULL; 2450a545947SLisandro Dalcin pc->ops->postsolve = NULL; 246b9ac7092SAlp Dener 2479566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &ctx->B)); 2489566063dSJacob Faibussowitsch PetscCall(MatSetType(ctx->B, MATLMVMBFGS)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->B, (PetscObject)pc, 1)); 250b9ac7092SAlp Dener PetscFunctionReturn(0); 251b9ac7092SAlp Dener } 252