xref: /petsc/src/ksp/pc/impls/lmvm/lmvmpc.c (revision 0e02354e4efb6ae7920bb0e7d191735ccbbd1e1e)
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 /*@
17b9ac7092SAlp Dener    PCLMVMSetMatLMVM - Replaces the LMVM matrix inside the preconditioner with
18b9ac7092SAlp Dener    the one provided by the user.
19b9ac7092SAlp Dener 
20b9ac7092SAlp Dener    Input Parameters:
21b9ac7092SAlp Dener +  pc - An LMVM preconditioner
22b9ac7092SAlp Dener -  B  - An LMVM-type matrix (MATLDFP, MATLBFGS, MATLSR1, MATLBRDN, MATLMBRDN, MATLSBRDN)
23b9ac7092SAlp Dener 
24b9ac7092SAlp Dener    Level: intermediate
25b9ac7092SAlp Dener @*/
26b9ac7092SAlp Dener PetscErrorCode PCLMVMSetMatLMVM(PC pc, Mat B)
27b9ac7092SAlp Dener {
28b9ac7092SAlp Dener   PC_LMVM          *ctx = (PC_LMVM*)pc->data;
29b9ac7092SAlp Dener   PetscBool        same;
30b9ac7092SAlp Dener 
31b9ac7092SAlp Dener   PetscFunctionBegin;
32b9ac7092SAlp Dener   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
33b9ac7092SAlp Dener   PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
3528b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type.");
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &same));
3728b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix must be an LMVM-type.");
389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->B));
399566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)B));
40b9ac7092SAlp Dener   ctx->B = B;
41b9ac7092SAlp Dener   PetscFunctionReturn(0);
42b9ac7092SAlp Dener }
43b9ac7092SAlp Dener 
44b9ac7092SAlp Dener /*@
45b9ac7092SAlp Dener    PCLMVMGetMatLMVM - Returns a pointer to the underlying LMVM matrix.
46b9ac7092SAlp Dener 
47b9ac7092SAlp Dener    Input Parameters:
48b9ac7092SAlp Dener .  pc - An LMVM preconditioner
49b9ac7092SAlp Dener 
50b9ac7092SAlp Dener    Output Parameters:
51b9ac7092SAlp Dener .  B - LMVM matrix inside the preconditioner
52b9ac7092SAlp Dener 
53b9ac7092SAlp Dener    Level: intermediate
54b9ac7092SAlp Dener @*/
55b9ac7092SAlp Dener PetscErrorCode PCLMVMGetMatLMVM(PC pc, Mat *B)
56b9ac7092SAlp Dener {
57b9ac7092SAlp Dener   PC_LMVM          *ctx = (PC_LMVM*)pc->data;
58b9ac7092SAlp Dener   PetscBool        same;
59b9ac7092SAlp Dener 
60b9ac7092SAlp Dener   PetscFunctionBegin;
61b9ac7092SAlp Dener   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
6328b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type.");
64b9ac7092SAlp Dener   *B = ctx->B;
65b9ac7092SAlp Dener   PetscFunctionReturn(0);
66b9ac7092SAlp Dener }
67b9ac7092SAlp Dener 
68b9ac7092SAlp Dener /*@
69b9ac7092SAlp Dener    PCLMVMSetIS - Sets the index sets that reduce the PC application.
70b9ac7092SAlp Dener 
71b9ac7092SAlp Dener    Input Parameters:
72b9ac7092SAlp Dener +  pc - An LMVM preconditioner
73b9ac7092SAlp Dener -  inactive - Index set defining the variables removed from the problem
74b9ac7092SAlp Dener 
75b9ac7092SAlp Dener    Level: intermediate
76b9ac7092SAlp Dener 
77db781477SPatrick Sanan .seealso: `MatLMVMUpdate()`
78b9ac7092SAlp Dener @*/
79b9ac7092SAlp Dener PetscErrorCode PCLMVMSetIS(PC pc, IS inactive)
80b9ac7092SAlp Dener {
81b9ac7092SAlp Dener   PC_LMVM          *ctx = (PC_LMVM*)pc->data;
82b9ac7092SAlp Dener   PetscBool        same;
83b9ac7092SAlp Dener 
84b9ac7092SAlp Dener   PetscFunctionBegin;
85b9ac7092SAlp Dener   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
86b2d8c577SAlp Dener   PetscValidHeaderSpecific(inactive, IS_CLASSID, 2);
879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
8828b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type.");
899566063dSJacob Faibussowitsch   PetscCall(PCLMVMClearIS(pc));
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)inactive));
91b9ac7092SAlp Dener   ctx->inactive = inactive;
92b9ac7092SAlp Dener   PetscFunctionReturn(0);
93b9ac7092SAlp Dener }
94b9ac7092SAlp Dener 
95b9ac7092SAlp Dener /*@
96b9ac7092SAlp Dener    PCLMVMClearIS - Removes the inactive variable index set.
97b9ac7092SAlp Dener 
98b9ac7092SAlp Dener    Input Parameters:
99b9ac7092SAlp Dener .  pc - An LMVM preconditioner
100b9ac7092SAlp Dener 
101b9ac7092SAlp Dener    Level: intermediate
102b9ac7092SAlp Dener 
103db781477SPatrick Sanan .seealso: `MatLMVMUpdate()`
104b9ac7092SAlp Dener @*/
105b9ac7092SAlp Dener PetscErrorCode PCLMVMClearIS(PC pc)
106b9ac7092SAlp Dener {
107b9ac7092SAlp Dener   PC_LMVM          *ctx = (PC_LMVM*)pc->data;
108b9ac7092SAlp Dener   PetscBool        same;
109b9ac7092SAlp Dener 
110b9ac7092SAlp Dener   PetscFunctionBegin;
111b9ac7092SAlp Dener   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
11328b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type.");
114b9ac7092SAlp Dener   if (ctx->inactive) {
1159566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ctx->inactive));
116b9ac7092SAlp Dener   }
117b9ac7092SAlp Dener   PetscFunctionReturn(0);
118b9ac7092SAlp Dener }
119b9ac7092SAlp Dener 
120b9ac7092SAlp Dener static PetscErrorCode PCApply_LMVM(PC pc,Vec x,Vec y)
121b9ac7092SAlp Dener {
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 
145b9ac7092SAlp Dener static PetscErrorCode PCReset_LMVM(PC pc)
146b9ac7092SAlp Dener {
147b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
148b9ac7092SAlp Dener 
149b9ac7092SAlp Dener   PetscFunctionBegin;
150b9ac7092SAlp Dener   if (ctx->xwork) {
1519566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ctx->xwork));
152b9ac7092SAlp Dener   }
153b9ac7092SAlp Dener   if (ctx->ywork) {
1549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ctx->ywork));
155b9ac7092SAlp Dener   }
156b9ac7092SAlp Dener   PetscFunctionReturn(0);
157b9ac7092SAlp Dener }
158b9ac7092SAlp Dener 
159b9ac7092SAlp Dener static PetscErrorCode PCSetUp_LMVM(PC pc)
160b9ac7092SAlp Dener {
161b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
162b9ac7092SAlp Dener   PetscInt       n, N;
163b2d8c577SAlp Dener   PetscBool      allocated;
164*0e02354eSStefano Zampini   const char     *prefix;
165b9ac7092SAlp Dener 
166b9ac7092SAlp Dener   PetscFunctionBegin;
1679566063dSJacob Faibussowitsch   PetscCall(MatLMVMIsAllocated(ctx->B, &allocated));
168b2d8c577SAlp Dener   if (!allocated) {
1699566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->mat, &ctx->xwork, &ctx->ywork));
1709566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ctx->xwork, &n));
1719566063dSJacob Faibussowitsch     PetscCall(VecGetSize(ctx->xwork, &N));
1729566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(ctx->B, n, n, N, N));
1739566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(ctx->B, ctx->xwork, ctx->ywork));
174b2d8c577SAlp Dener   } else {
1759566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(ctx->B, &ctx->xwork, &ctx->ywork));
176b2d8c577SAlp Dener   }
177*0e02354eSStefano Zampini   PetscCall(PCGetOptionsPrefix(pc, &prefix));
178*0e02354eSStefano Zampini   PetscCall(MatSetOptionsPrefix(ctx->B, prefix));
179*0e02354eSStefano Zampini   PetscCall(MatAppendOptionsPrefix(ctx->B, "pc_lmvm_"));
180b9ac7092SAlp Dener   PetscFunctionReturn(0);
181b9ac7092SAlp Dener }
182b9ac7092SAlp Dener 
183b9ac7092SAlp Dener static PetscErrorCode PCView_LMVM(PC pc,PetscViewer viewer)
184b9ac7092SAlp Dener {
185b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
186b9ac7092SAlp Dener   PetscBool      iascii;
187b9ac7092SAlp Dener 
188b9ac7092SAlp Dener   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
190dbf6bb8dSprj-   if (iascii && ctx->B->assembled) {
1919566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1929566063dSJacob Faibussowitsch     PetscCall(MatView(ctx->B, viewer));
1939566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
194b9ac7092SAlp Dener   }
195b9ac7092SAlp Dener   PetscFunctionReturn(0);
196b9ac7092SAlp Dener }
197b9ac7092SAlp Dener 
198b9ac7092SAlp Dener static PetscErrorCode PCSetFromOptions_LMVM(PetscOptionItems* PetscOptionsObject, PC pc)
199b9ac7092SAlp Dener {
200b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
201*0e02354eSStefano Zampini   const char     *prefix;
202b9ac7092SAlp Dener 
203b9ac7092SAlp Dener   PetscFunctionBegin;
204*0e02354eSStefano Zampini   PetscCall(PCGetOptionsPrefix(pc, &prefix));
205*0e02354eSStefano Zampini   PetscCall(MatSetOptionsPrefix(ctx->B, prefix));
206*0e02354eSStefano Zampini   PetscCall(MatAppendOptionsPrefix(ctx->B, "pc_lmvm_"));
2079566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(ctx->B));
208b9ac7092SAlp Dener   PetscFunctionReturn(0);
209b9ac7092SAlp Dener }
210b9ac7092SAlp Dener 
211b9ac7092SAlp Dener static PetscErrorCode PCDestroy_LMVM(PC pc)
212b9ac7092SAlp Dener {
213b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
214b9ac7092SAlp Dener 
215b9ac7092SAlp Dener   PetscFunctionBegin;
216b9ac7092SAlp Dener   if (ctx->inactive) {
2179566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ctx->inactive));
218b9ac7092SAlp Dener   }
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));
225b9ac7092SAlp Dener   PetscFunctionReturn(0);
226b9ac7092SAlp Dener }
227b9ac7092SAlp Dener 
228b9ac7092SAlp Dener /*MC
229b9ac7092SAlp Dener    PCLMVM - Creates a preconditioner around an LMVM matrix. Options for the
230b9ac7092SAlp Dener             underlying LMVM matrix can be access with the "-pc_lmvm_" prefix.
231b9ac7092SAlp Dener 
23206dd6b0eSSatish Balay    Level: intermediate
23306dd6b0eSSatish Balay 
234db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`,
235db781477SPatrick Sanan           `PC`, `MATLMVM`, `PCLMVMUpdate()`, `PCLMVMSetMatLMVM()`, `PCLMVMGetMatLMVM()`
236b9ac7092SAlp Dener M*/
237b9ac7092SAlp Dener PETSC_EXTERN PetscErrorCode PCCreate_LMVM(PC pc)
238b9ac7092SAlp Dener {
239b9ac7092SAlp Dener   PC_LMVM        *ctx;
240b9ac7092SAlp Dener 
241b9ac7092SAlp Dener   PetscFunctionBegin;
2429566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&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(PCSetReusePreconditioner(pc, PETSC_TRUE));
259b9ac7092SAlp Dener 
2609566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &ctx->B));
2619566063dSJacob Faibussowitsch   PetscCall(MatSetType(ctx->B, MATLMVMBFGS));
2629566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->B, (PetscObject)pc, 1));
263b9ac7092SAlp Dener   PetscFunctionReturn(0);
264b9ac7092SAlp Dener }
265