xref: /petsc/src/ksp/pc/impls/lmvm/lmvmpc.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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);
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
35*28b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type.");
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &same));
37*28b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix must be an LMVM-type.");
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&ctx->B));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(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);
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
63*28b400f6SJacob 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 
77b9ac7092SAlp Dener .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);
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
88*28b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type.");
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PCLMVMClearIS(pc));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
103b9ac7092SAlp Dener .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);
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same));
113*28b400f6SJacob Faibussowitsch   PetscCheck(same,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type.");
114b9ac7092SAlp Dener   if (ctx->inactive) {
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(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) {
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(ctx->xwork));
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetSubVector(ctx->xwork, ctx->inactive, &xsub));
1295f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(x, xsub));
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreSubVector(ctx->xwork, ctx->inactive, &xsub));
131b9ac7092SAlp Dener   } else {
1325f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(x, ctx->xwork));
133b9ac7092SAlp Dener   }
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolve(ctx->B, ctx->xwork, ctx->ywork));
135b9ac7092SAlp Dener   if (ctx->inactive) {
1365f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetSubVector(ctx->ywork, ctx->inactive, &ysub));
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(ysub, y));
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreSubVector(ctx->ywork, ctx->inactive, &ysub));
139b9ac7092SAlp Dener   } else {
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(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) {
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&ctx->xwork));
152b9ac7092SAlp Dener   }
153b9ac7092SAlp Dener   if (ctx->ywork) {
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(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;
164b9ac7092SAlp Dener 
165b9ac7092SAlp Dener   PetscFunctionBegin;
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLMVMIsAllocated(ctx->B, &allocated));
167b2d8c577SAlp Dener   if (!allocated) {
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(pc->mat, &ctx->xwork, &ctx->ywork));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(ctx->xwork, &n));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetSize(ctx->xwork, &N));
1715f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(ctx->B, n, n, N, N));
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLMVMAllocate(ctx->B, ctx->xwork, ctx->ywork));
173b2d8c577SAlp Dener   } else {
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(ctx->B, &ctx->xwork, &ctx->ywork));
175b2d8c577SAlp Dener   }
176b9ac7092SAlp Dener   PetscFunctionReturn(0);
177b9ac7092SAlp Dener }
178b9ac7092SAlp Dener 
179b9ac7092SAlp Dener static PetscErrorCode PCView_LMVM(PC pc,PetscViewer viewer)
180b9ac7092SAlp Dener {
181b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
182b9ac7092SAlp Dener   PetscBool      iascii;
183b9ac7092SAlp Dener 
184b9ac7092SAlp Dener   PetscFunctionBegin;
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
186dbf6bb8dSprj-   if (iascii && ctx->B->assembled) {
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(ctx->B, viewer));
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(viewer));
190b9ac7092SAlp Dener   }
191b9ac7092SAlp Dener   PetscFunctionReturn(0);
192b9ac7092SAlp Dener }
193b9ac7092SAlp Dener 
194b9ac7092SAlp Dener static PetscErrorCode PCSetFromOptions_LMVM(PetscOptionItems* PetscOptionsObject, PC pc)
195b9ac7092SAlp Dener {
196b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
197b9ac7092SAlp Dener 
198b9ac7092SAlp Dener   PetscFunctionBegin;
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(ctx->B));
200b9ac7092SAlp Dener   PetscFunctionReturn(0);
201b9ac7092SAlp Dener }
202b9ac7092SAlp Dener 
203b9ac7092SAlp Dener static PetscErrorCode PCDestroy_LMVM(PC pc)
204b9ac7092SAlp Dener {
205b9ac7092SAlp Dener   PC_LMVM        *ctx = (PC_LMVM*)pc->data;
206b9ac7092SAlp Dener 
207b9ac7092SAlp Dener   PetscFunctionBegin;
208b9ac7092SAlp Dener   if (ctx->inactive) {
2095f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&ctx->inactive));
210b9ac7092SAlp Dener   }
211b9ac7092SAlp Dener   if (pc->setupcalled) {
2125f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&ctx->xwork));
2135f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&ctx->ywork));
214b9ac7092SAlp Dener   }
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&ctx->B));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pc->data));
217b9ac7092SAlp Dener   PetscFunctionReturn(0);
218b9ac7092SAlp Dener }
219b9ac7092SAlp Dener 
220b9ac7092SAlp Dener /*MC
221b9ac7092SAlp Dener    PCLMVM - Creates a preconditioner around an LMVM matrix. Options for the
222b9ac7092SAlp Dener             underlying LMVM matrix can be access with the "-pc_lmvm_" prefix.
223b9ac7092SAlp Dener 
22406dd6b0eSSatish Balay    Level: intermediate
22506dd6b0eSSatish Balay 
226b9ac7092SAlp Dener .seealso:  PCCreate(), PCSetType(), PCType (for list of available types),
227b9ac7092SAlp Dener            PC, MATLMVM, PCLMVMUpdate(), PCLMVMSetMatLMVM(), PCLMVMGetMatLMVM()
228b9ac7092SAlp Dener M*/
229b9ac7092SAlp Dener PETSC_EXTERN PetscErrorCode PCCreate_LMVM(PC pc)
230b9ac7092SAlp Dener {
231b9ac7092SAlp Dener   PC_LMVM        *ctx;
232b9ac7092SAlp Dener 
233b9ac7092SAlp Dener   PetscFunctionBegin;
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(pc,&ctx));
235b9ac7092SAlp Dener   pc->data = (void*)ctx;
236b9ac7092SAlp Dener 
237b9ac7092SAlp Dener   pc->ops->reset           = PCReset_LMVM;
238b9ac7092SAlp Dener   pc->ops->setup           = PCSetUp_LMVM;
239b9ac7092SAlp Dener   pc->ops->destroy         = PCDestroy_LMVM;
240b9ac7092SAlp Dener   pc->ops->view            = PCView_LMVM;
241b9ac7092SAlp Dener   pc->ops->apply           = PCApply_LMVM;
242b9ac7092SAlp Dener   pc->ops->setfromoptions  = PCSetFromOptions_LMVM;
2430a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
2440a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
2450a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
2460a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
2470a545947SLisandro Dalcin   pc->ops->presolve        = NULL;
2480a545947SLisandro Dalcin   pc->ops->postsolve       = NULL;
249b9ac7092SAlp Dener 
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetReusePreconditioner(pc, PETSC_TRUE));
251b9ac7092SAlp Dener 
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject)pc), &ctx->B));
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(ctx->B, MATLMVMBFGS));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectIncrementTabLevel((PetscObject)ctx->B, (PetscObject)pc, 1));
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(ctx->B, "pc_lmvm_"));
256b9ac7092SAlp Dener   PetscFunctionReturn(0);
257b9ac7092SAlp Dener }
258