xref: /petsc/src/ts/utils/dmts.c (revision 5edf65166e22df6feda2ed6bd8a6e348419dee17)
1 #include <petsc-private/tsimpl.h>     /*I "petscts.h" I*/
2 #include <petsc-private/dmimpl.h>     /*I "petscdm.h" I*/
3 
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMCoarsenHook_TSDM"
7 /* Attaches the SNESDM to the coarse level.
8  * Under what conditions should we copy versus duplicate?
9  */
10 static PetscErrorCode DMCoarsenHook_TSDM(DM dm,DM dmc,void *ctx)
11 {
12   PetscErrorCode ierr;
13 
14   PetscFunctionBegin;
15   ierr = DMTSCopyContext(dm,dmc);CHKERRQ(ierr);
16   PetscFunctionReturn(0);
17 }
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "DMRestrictHook_TSDM"
21 /* This could restrict auxiliary information to the coarse level.
22  */
23 static PetscErrorCode DMRestrictHook_TSDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
24 {
25 
26   PetscFunctionBegin;
27   PetscFunctionReturn(0);
28 }
29 
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "PetscContainerDestroy_TSDM"
33 static PetscErrorCode PetscContainerDestroy_TSDM(void *ctx)
34 {
35   PetscErrorCode ierr;
36   TSDM tsdm = (TSDM)ctx;
37 
38   PetscFunctionBegin;
39   if (tsdm->destroy) {ierr = (*tsdm->destroy)(tsdm);CHKERRQ(ierr);}
40   ierr = PetscFree(tsdm);CHKERRQ(ierr);
41   PetscFunctionReturn(0);
42 }
43 
44 #undef __FUNCT__
45 #define __FUNCT__ "DMTSGetContext"
46 /*@C
47    DMTSGetContext - get read-only private TSDM context from a DM
48 
49    Not Collective
50 
51    Input Argument:
52 .  dm - DM to be used with TS
53 
54    Output Argument:
55 .  tsdm - private TSDM context
56 
57    Level: developer
58 
59    Notes:
60    Use DMTSGetContextWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
61 
62 .seealso: DMTSGetContextWrite()
63 @*/
64 PetscErrorCode DMTSGetContext(DM dm,TSDM *tsdm)
65 {
66   PetscErrorCode ierr;
67   PetscContainer container;
68   TSDM           tsdmnew;
69 
70 
71   PetscFunctionBegin;
72   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
73   ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
74   if (container) {
75     ierr = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr);
76   } else {
77     ierr = PetscInfo(dm,"Creating new TSDM\n");CHKERRQ(ierr);
78     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
79     ierr = PetscNewLog(dm,struct _n_TSDM,&tsdmnew);CHKERRQ(ierr);
80     ierr = PetscContainerSetPointer(container,tsdmnew);CHKERRQ(ierr);
81     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_TSDM);CHKERRQ(ierr);
82     ierr = PetscObjectCompose((PetscObject)dm,"TSDM",(PetscObject)container);CHKERRQ(ierr);
83     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSDM,DMRestrictHook_TSDM,PETSC_NULL);CHKERRQ(ierr);
84     ierr = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr);
85     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
86   }
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "DMTSGetContextWrite"
92 /*@C
93    DMTSGetContextWrite - get write access to private TSDM context from a DM
94 
95    Not Collective
96 
97    Input Argument:
98 .  dm - DM to be used with TS
99 
100    Output Argument:
101 .  tsdm - private TSDM context
102 
103    Level: developer
104 
105 .seealso: DMTSGetContext()
106 @*/
107 PetscErrorCode DMTSGetContextWrite(DM dm,TSDM *tsdm)
108 {
109   PetscErrorCode ierr;
110   TSDM           sdm;
111 
112   PetscFunctionBegin;
113   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
114   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
115   if (!sdm->originaldm) sdm->originaldm = dm;
116   if (sdm->originaldm != dm) {  /* Copy on write */
117     PetscContainer container;
118     TSDM         oldsdm = sdm;
119     ierr = PetscInfo(dm,"Copying TSDM due to write\n");CHKERRQ(ierr);
120     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
121     ierr = PetscNewLog(dm,struct _n_TSDM,&sdm);CHKERRQ(ierr);
122     ierr = PetscMemcpy(sdm,oldsdm,sizeof(*sdm));CHKERRQ(ierr);
123     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
124     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_TSDM);CHKERRQ(ierr);
125     ierr = PetscObjectCompose((PetscObject)dm,"TSDM",(PetscObject)container);CHKERRQ(ierr);
126     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
127   }
128   *tsdm = sdm;
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "DMTSCopyContext"
134 /*@C
135    DMTSCopyContext - copies a DM context to a new DM
136 
137    Logically Collective
138 
139    Input Arguments:
140 +  dmsrc - DM to obtain context from
141 -  dmdest - DM to add context to
142 
143    Level: developer
144 
145    Note:
146    The context is copied by reference. This function does not ensure that a context exists.
147 
148 .seealso: DMTSGetContext(), TSSetDM()
149 @*/
150 PetscErrorCode DMTSCopyContext(DM dmsrc,DM dmdest)
151 {
152   PetscErrorCode ierr;
153   PetscContainer container;
154 
155   PetscFunctionBegin;
156   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
157   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
158   ierr = PetscObjectQuery((PetscObject)dmsrc,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
159   if (container) {
160     ierr = PetscObjectCompose((PetscObject)dmdest,"TSDM",(PetscObject)container);CHKERRQ(ierr);
161     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_TSDM,DMRestrictHook_TSDM,PETSC_NULL);CHKERRQ(ierr);
162   }
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "DMTSSetIFunction"
168 /*@C
169    DMTSSetIFunction - set TS implicit function evaluation function
170 
171    Not Collective
172 
173    Input Arguments:
174 +  dm - DM to be used with TS
175 .  func - function evaluation function, see TSSetIFunction() for calling sequence
176 -  ctx - context for residual evaluation
177 
178    Level: advanced
179 
180    Note:
181    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
182    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
183    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
184 
185 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
186 @*/
187 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
188 {
189   PetscErrorCode ierr;
190   TSDM tsdm;
191 
192   PetscFunctionBegin;
193   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
194   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
195   if (func) tsdm->ifunction = func;
196   if (ctx)  tsdm->ifunctionctx = ctx;
197   PetscFunctionReturn(0);
198 }
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "DMTSGetIFunction"
202 /*@C
203    DMTSGetIFunction - get TS implicit residual evaluation function
204 
205    Not Collective
206 
207    Input Argument:
208 .  dm - DM to be used with TS
209 
210    Output Arguments:
211 +  func - function evaluation function, see TSSetIFunction() for calling sequence
212 -  ctx - context for residual evaluation
213 
214    Level: advanced
215 
216    Note:
217    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
218    associated with the DM.
219 
220 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
221 @*/
222 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
223 {
224   PetscErrorCode ierr;
225   TSDM tsdm;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
229   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
230   if (func) *func = tsdm->ifunction;
231   if (ctx)  *ctx = tsdm->ifunctionctx;
232   PetscFunctionReturn(0);
233 }
234 
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "DMTSSetRHSFunction"
238 /*@C
239    DMTSSetRHSFunction - set TS explicit residual evaluation function
240 
241    Not Collective
242 
243    Input Arguments:
244 +  dm - DM to be used with TS
245 .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
246 -  ctx - context for residual evaluation
247 
248    Level: advanced
249 
250    Note:
251    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
252    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
253    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
254 
255 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
256 @*/
257 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
258 {
259   PetscErrorCode ierr;
260   TSDM           tsdm;
261 
262   PetscFunctionBegin;
263   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
264   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
265   if (func) tsdm->rhsfunction = func;
266   if (ctx)  tsdm->rhsfunctionctx = ctx;
267   PetscFunctionReturn(0);
268 }
269 
270 #undef __FUNCT__
271 #define __FUNCT__ "DMTSGetSolutionFunction"
272 /*@C
273    DMTSGetSolutionFunction - gets the TS solution evaluation function
274 
275    Not Collective
276 
277    Input Arguments:
278 .  dm - DM to be used with TS
279 
280    Output Parameters:
281 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
282 -  ctx - context for solution evaluation
283 
284    Level: advanced
285 
286 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
287 @*/
288 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
289 {
290   PetscErrorCode ierr;
291   TSDM           tsdm;
292 
293   PetscFunctionBegin;
294   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
295   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
296   if (func) *func = tsdm->solution;
297   if (ctx)  *ctx  = tsdm->solutionctx;
298   PetscFunctionReturn(0);
299 }
300 
301 #undef __FUNCT__
302 #define __FUNCT__ "DMTSSetSolutionFunction"
303 /*@C
304    DMTSSetSolutionFunction - set TS solution evaluation function
305 
306    Not Collective
307 
308    Input Arguments:
309 +  dm - DM to be used with TS
310 .  func - solution function evaluation function, see TSSetSolution() for calling sequence
311 -  ctx - context for solution evaluation
312 
313    Level: advanced
314 
315    Note:
316    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
317    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
318    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
319 
320 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
321 @*/
322 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
323 {
324   PetscErrorCode ierr;
325   TSDM           tsdm;
326 
327   PetscFunctionBegin;
328   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
329   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
330   if (func) tsdm->solution    = func;
331   if (ctx)  tsdm->solutionctx = ctx;
332   PetscFunctionReturn(0);
333 }
334 
335 #undef __FUNCT__
336 #define __FUNCT__ "DMTSGetRHSFunction"
337 /*@C
338    DMTSGetRHSFunction - get TS explicit residual evaluation function
339 
340    Not Collective
341 
342    Input Argument:
343 .  dm - DM to be used with TS
344 
345    Output Arguments:
346 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
347 -  ctx - context for residual evaluation
348 
349    Level: advanced
350 
351    Note:
352    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
353    associated with the DM.
354 
355 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
356 @*/
357 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
358 {
359   PetscErrorCode ierr;
360   TSDM tsdm;
361 
362   PetscFunctionBegin;
363   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
364   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
365   if (func) *func = tsdm->rhsfunction;
366   if (ctx)  *ctx = tsdm->rhsfunctionctx;
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "DMTSSetIJacobian"
372 /*@C
373    DMTSSetIJacobian - set TS Jacobian evaluation function
374 
375    Not Collective
376 
377    Input Argument:
378 +  dm - DM to be used with TS
379 .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
380 -  ctx - context for residual evaluation
381 
382    Level: advanced
383 
384    Note:
385    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
386    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
387    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
388 
389 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
390 @*/
391 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
392 {
393   PetscErrorCode ierr;
394   TSDM sdm;
395 
396   PetscFunctionBegin;
397   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
398   ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr);
399   if (func) sdm->ijacobian = func;
400   if (ctx)  sdm->ijacobianctx = ctx;
401   PetscFunctionReturn(0);
402 }
403 
404 #undef __FUNCT__
405 #define __FUNCT__ "DMTSGetIJacobian"
406 /*@C
407    DMTSGetIJacobian - get TS Jacobian evaluation function
408 
409    Not Collective
410 
411    Input Argument:
412 .  dm - DM to be used with TS
413 
414    Output Arguments:
415 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
416 -  ctx - context for residual evaluation
417 
418    Level: advanced
419 
420    Note:
421    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
422    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
423    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
424 
425 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
426 @*/
427 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
428 {
429   PetscErrorCode ierr;
430   TSDM tsdm;
431 
432   PetscFunctionBegin;
433   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
434   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
435   if (func) *func = tsdm->ijacobian;
436   if (ctx)  *ctx = tsdm->ijacobianctx;
437   PetscFunctionReturn(0);
438 }
439 
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "DMTSSetRHSJacobian"
443 /*@C
444    DMTSSetRHSJacobian - set TS Jacobian evaluation function
445 
446    Not Collective
447 
448    Input Argument:
449 +  dm - DM to be used with TS
450 .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
451 -  ctx - context for residual evaluation
452 
453    Level: advanced
454 
455    Note:
456    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
457    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
458    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
459 
460 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
461 @*/
462 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
463 {
464   PetscErrorCode ierr;
465   TSDM tsdm;
466 
467   PetscFunctionBegin;
468   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
469   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
470   if (func) tsdm->rhsjacobian = func;
471   if (ctx)  tsdm->rhsjacobianctx = ctx;
472   PetscFunctionReturn(0);
473 }
474 
475 #undef __FUNCT__
476 #define __FUNCT__ "DMTSGetRHSJacobian"
477 /*@C
478    DMTSGetRHSJacobian - get TS Jacobian evaluation function
479 
480    Not Collective
481 
482    Input Argument:
483 .  dm - DM to be used with TS
484 
485    Output Arguments:
486 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
487 -  ctx - context for residual evaluation
488 
489    Level: advanced
490 
491    Note:
492    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
493    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
494    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
495 
496 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
497 @*/
498 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
499 {
500   PetscErrorCode ierr;
501   TSDM tsdm;
502 
503   PetscFunctionBegin;
504   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
505   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
506   if (func) *func = tsdm->rhsjacobian;
507   if (ctx)  *ctx = tsdm->rhsjacobianctx;
508   PetscFunctionReturn(0);
509 }
510