xref: /petsc/src/snes/utils/dmsnes.c (revision 670f3ff9a12f3667ff7d4cebfc50ebc8579c9e33)
1 #include <private/snesimpl.h>   /*I "petscsnes.h" I*/
2 #include <petscdm.h>            /*I "petscdm.h" I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "PetscContainerDestroy_SNESDM"
6 static PetscErrorCode PetscContainerDestroy_SNESDM(void *ctx)
7 {
8   PetscErrorCode ierr;
9   SNESDM sdm = (SNESDM)ctx;
10 
11   PetscFunctionBegin;
12   if (sdm->dmloopref) {         /* Reconnect old loop before we drop our reference */
13     DM dminner;
14     ierr = PetscObjectQuery((PetscObject)sdm->vec_sol,"DM",(PetscObject*)&dminner);CHKERRQ(ierr);
15     if (dminner) {ierr = PetscObjectReference((PetscObject)dminner);CHKERRQ(ierr);}
16   }
17   ierr = VecDestroy(&sdm->vec_sol);CHKERRQ(ierr);
18   if (sdm->destroy) {ierr = (*sdm->destroy)(sdm);CHKERRQ(ierr);}
19   ierr = PetscFree(sdm);CHKERRQ(ierr);
20   PetscFunctionReturn(0);
21 }
22 
23 #undef __FUNCT__
24 #define __FUNCT__ "DMCoarsenHook_SNESDM"
25 /* Attaches the SNESDM to the coarse level.
26  * Under what conditions should we copy versus duplicate?
27  */
28 static PetscErrorCode DMCoarsenHook_SNESDM(DM dm,DM dmc,void *ctx)
29 {
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   ierr = DMSNESCopyContext(dm,dmc);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 #undef __FUNCT__
38 #define __FUNCT__ "DMRestrictHook_SNESDM"
39 /* Restrict state vector from finest level
40  */
41 static PetscErrorCode DMRestrictHook_SNESDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
42 {
43   PetscErrorCode ierr;
44   SNESDM sdm,sdmc;
45   Vec Xcoarse;
46 
47   PetscFunctionBegin;
48   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
49   ierr = DMSNESGetContextWrite(dmc,&sdmc);CHKERRQ(ierr);
50   if (sdm->vec_sol) {
51     ierr = DMSNESGetSolution(dmc,&Xcoarse);CHKERRQ(ierr);
52     ierr = MatRestrict(Restrict,sdm->vec_sol,Xcoarse);CHKERRQ(ierr);
53     ierr = VecPointwiseMult(Xcoarse,Xcoarse,rscale);CHKERRQ(ierr);
54   }
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "DMSNESGetContext"
60 /*@C
61    DMSNESGetContext - get read-only private SNESDM context from a DM
62 
63    Not Collective
64 
65    Input Argument:
66 .  dm - DM to be used with SNES
67 
68    Output Argument:
69 .  snesdm - private SNESDM context
70 
71    Level: developer
72 
73    Notes:
74    Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
75 
76 .seealso: DMSNESGetContextWrite()
77 @*/
78 PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm)
79 {
80   PetscErrorCode ierr;
81   PetscContainer container;
82   SNESDM         sdm;
83 
84   PetscFunctionBegin;
85   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
86   ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
87   if (container) {
88     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
89   } else {
90     ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr);
91     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
92     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
93     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
94     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
95     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
96     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
97     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
98     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
99   }
100   PetscFunctionReturn(0);
101 }
102 
103 #undef __FUNCT__
104 #define __FUNCT__ "DMSNESGetContextWrite"
105 /*@C
106    DMSNESGetContextWrite - get write access to private SNESDM context from a DM
107 
108    Not Collective
109 
110    Input Argument:
111 .  dm - DM to be used with SNES
112 
113    Output Argument:
114 .  snesdm - private SNESDM context
115 
116    Level: developer
117 
118 .seealso: DMSNESGetContext()
119 @*/
120 PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm)
121 {
122   PetscErrorCode ierr;
123   SNESDM         sdm;
124 
125   PetscFunctionBegin;
126   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
127   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
128   if (!sdm->originaldm) sdm->originaldm = dm;
129   if (sdm->originaldm != dm) {  /* Copy on write */
130     PetscContainer container;
131     SNESDM         oldsdm = sdm;
132     ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr);
133     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
134     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
135     ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr);
136     sdm->vec_sol = PETSC_NULL;
137     sdm->dmloopref = PETSC_FALSE;
138     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
139     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
140     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
141     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
142   }
143   *snesdm = sdm;
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNCT__
148 #define __FUNCT__ "DMSNESCopyContext"
149 /*@C
150    DMSNESCopyContext - copies a DM context to a new DM
151 
152    Logically Collective
153 
154    Input Arguments:
155 +  dmsrc - DM to obtain context from
156 -  dmdest - DM to add context to
157 
158    Level: developer
159 
160    Note:
161    The context is copied by reference. This function does not ensure that a context exists.
162 
163 .seealso: DMSNESGetContext(), SNESSetDM()
164 @*/
165 PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest)
166 {
167   PetscErrorCode ierr;
168   PetscContainer container;
169 
170   PetscFunctionBegin;
171   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
172   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
173   ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
174   if (container) {
175     ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
176     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
177   }
178   PetscFunctionReturn(0);
179 }
180 
181 #undef __FUNCT__
182 #define __FUNCT__ "DMSNESSetFunction"
183 /*@C
184    DMSNESSetFunction - set SNES residual evaluation function
185 
186    Not Collective
187 
188    Input Arguments:
189 +  dm - DM to be used with SNES
190 .  func - residual evaluation function, see SNESSetFunction() for calling sequence
191 -  ctx - context for residual evaluation
192 
193    Level: advanced
194 
195    Note:
196    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
197    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
198    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
199 
200 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
201 @*/
202 PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
203 {
204   PetscErrorCode ierr;
205   SNESDM sdm;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
209   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
210   if (func) sdm->computefunction = func;
211   if (ctx)  sdm->functionctx = ctx;
212   PetscFunctionReturn(0);
213 }
214 
215 #undef __FUNCT__
216 #define __FUNCT__ "DMSNESGetFunction"
217 /*@C
218    DMSNESGetFunction - get SNES residual evaluation function
219 
220    Not Collective
221 
222    Input Argument:
223 .  dm - DM to be used with SNES
224 
225    Output Arguments:
226 +  func - residual evaluation function, see SNESSetFunction() for calling sequence
227 -  ctx - context for residual evaluation
228 
229    Level: advanced
230 
231    Note:
232    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
233    associated with the DM.
234 
235 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction()
236 @*/
237 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
238 {
239   PetscErrorCode ierr;
240   SNESDM sdm;
241 
242   PetscFunctionBegin;
243   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
244   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
245   if (func) *func = sdm->computefunction;
246   if (ctx)  *ctx = sdm->functionctx;
247   PetscFunctionReturn(0);
248 }
249 
250 #undef __FUNCT__
251 #define __FUNCT__ "DMSNESSetGS"
252 /*@C
253    DMSNESSetGS - set SNES Gauss-Seidel relaxation function
254 
255    Not Collective
256 
257    Input Argument:
258 +  dm - DM to be used with SNES
259 .  func - relaxation function, see SNESSetGS() for calling sequence
260 -  ctx - context for residual evaluation
261 
262    Level: advanced
263 
264    Note:
265    SNESSetGS() is normally used, but it calls this function internally because the user context is actually
266    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
267    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
268 
269 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction()
270 @*/
271 PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
272 {
273   PetscErrorCode ierr;
274   SNESDM sdm;
275 
276   PetscFunctionBegin;
277   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
278   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
279   if (func) sdm->computegs = func;
280   if (ctx)  sdm->gsctx = ctx;
281   PetscFunctionReturn(0);
282 }
283 
284 #undef __FUNCT__
285 #define __FUNCT__ "DMSNESGetGS"
286 /*@C
287    DMSNESGetGS - get SNES Gauss-Seidel relaxation function
288 
289    Not Collective
290 
291    Input Argument:
292 .  dm - DM to be used with SNES
293 
294    Output Arguments:
295 +  func - relaxation function, see SNESSetGS() for calling sequence
296 -  ctx - context for residual evaluation
297 
298    Level: advanced
299 
300    Note:
301    SNESGetGS() is normally used, but it calls this function internally because the user context is actually
302    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
303    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
304 
305 .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction()
306 @*/
307 PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
308 {
309   PetscErrorCode ierr;
310   SNESDM sdm;
311 
312   PetscFunctionBegin;
313   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
314   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
315   if (func) *func = sdm->computegs;
316   if (ctx)  *ctx = sdm->gsctx;
317   PetscFunctionReturn(0);
318 }
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "DMSNESSetJacobian"
322 /*@C
323    DMSNESSetFunction - set SNES Jacobian evaluation function
324 
325    Not Collective
326 
327    Input Argument:
328 +  dm - DM to be used with SNES
329 .  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
330 -  ctx - context for residual evaluation
331 
332    Level: advanced
333 
334    Note:
335    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
336    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
337    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
338 
339 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian()
340 @*/
341 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
342 {
343   PetscErrorCode ierr;
344   SNESDM sdm;
345 
346   PetscFunctionBegin;
347   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
348   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
349   if (func) sdm->computejacobian = func;
350   if (ctx)  sdm->jacobianctx = ctx;
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "DMSNESGetJacobian"
356 /*@C
357    DMSNESGetFunction - get SNES Jacobian evaluation function
358 
359    Not Collective
360 
361    Input Argument:
362 .  dm - DM to be used with SNES
363 
364    Output Arguments:
365 +  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
366 -  ctx - context for residual evaluation
367 
368    Level: advanced
369 
370    Note:
371    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
372    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
373    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
374 
375 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
376 @*/
377 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
378 {
379   PetscErrorCode ierr;
380   SNESDM sdm;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
384   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
385   if (func) *func = sdm->computejacobian;
386   if (ctx)  *ctx = sdm->jacobianctx;
387   PetscFunctionReturn(0);
388 }
389 
390 #undef __FUNCT__
391 #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy"
392 static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx)
393 {
394   PetscErrorCode ierr;
395   DM             dm;
396 
397   PetscFunctionBegin;
398   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
399   ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr);
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy"
405 static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
406 {
407   PetscErrorCode ierr;
408   DM             dm;
409 
410   PetscFunctionBegin;
411   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
412   ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr);
413   PetscFunctionReturn(0);
414 }
415 
416 #undef __FUNCT__
417 #define __FUNCT__ "DMSNESSetUpLegacy"
418 /* Sets up calling of legacy DM routines */
419 PetscErrorCode DMSNESSetUpLegacy(DM dm)
420 {
421   PetscErrorCode ierr;
422   SNESDM         sdm;
423 
424   PetscFunctionBegin;
425   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
426   if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
427   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
428   if (!sdm->computejacobian) {ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
429   PetscFunctionReturn(0);
430 }
431 
432 #undef __FUNCT__
433 #define __FUNCT__ "DMSNESGetSolution"
434 PetscErrorCode DMSNESGetSolution(DM dm,Vec *vec_sol)
435 {
436   PetscErrorCode ierr;
437   SNESDM sdm;
438 
439   PetscFunctionBegin;
440   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
441   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
442   if (!sdm->vec_sol) {
443     ierr = DMCreateGlobalVector(dm,&sdm->vec_sol);CHKERRQ(ierr);
444     /* This vector holds an extra reference to the DM, decrement the DM reference count so that it remains constant */
445     sdm->dmloopref = PETSC_TRUE;
446     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
447   }
448   *vec_sol = sdm->vec_sol;
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "DMSNESSetSolution"
454 PetscErrorCode DMSNESSetSolution(DM dm,Vec vec_sol)
455 {
456   PetscErrorCode ierr;
457   DM dminner;
458   SNESDM sdm;
459 
460   PetscFunctionBegin;
461   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
462   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
463   if (sdm->vec_sol == vec_sol) PetscFunctionReturn(0);
464   if (sdm->dmloopref) {         /* Reconnect old loop before we drop our reference */
465     ierr = PetscObjectQuery((PetscObject)sdm->vec_sol,"DM",(PetscObject*)&dminner);CHKERRQ(ierr);
466     if (dminner) {ierr = PetscObjectReference((PetscObject)dminner);CHKERRQ(ierr);}
467   }
468   sdm->dmloopref = PETSC_FALSE;
469 
470   if (vec_sol) {ierr = PetscObjectReference((PetscObject)vec_sol);CHKERRQ(ierr);}
471   ierr = VecDestroy(&sdm->vec_sol);CHKERRQ(ierr);
472   sdm->vec_sol = vec_sol;
473 
474   /* Break new loop if we just created one */
475   if (sdm->vec_sol) {
476     ierr = PetscObjectQuery((PetscObject)sdm->vec_sol,"DM",(PetscObject*)&dminner);CHKERRQ(ierr);
477     if (dminner == dm) {
478       sdm->dmloopref = PETSC_TRUE;
479       if (((PetscObject)dm)->refct <= 1) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONG,"The incoming Vec holds the last reference to the DM");
480       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
481     }
482   }
483   PetscFunctionReturn(0);
484 }
485