xref: /petsc/src/ts/utils/dmts.c (revision 24989b8cf3e165ebaa37f4b2d30151131f137c58)
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_SNESDM"
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 = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr);
84     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
85   }
86   PetscFunctionReturn(0);
87 }
88 
89 #undef __FUNCT__
90 #define __FUNCT__ "DMTSGetContextWrite"
91 /*@C
92    DMTSGetContextWrite - get write access to private TSDM context from a DM
93 
94    Not Collective
95 
96    Input Argument:
97 .  dm - DM to be used with TS
98 
99    Output Argument:
100 .  tsdm - private TSDM context
101 
102    Level: developer
103 
104 .seealso: DMTSGetContext()
105 @*/
106 PetscErrorCode DMTSGetContextWrite(DM dm,TSDM *tsdm)
107 {
108   PetscErrorCode ierr;
109   TSDM           sdm;
110 
111   PetscFunctionBegin;
112   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
113   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
114   if (!sdm->originaldm) sdm->originaldm = dm;
115   if (sdm->originaldm != dm) {  /* Copy on write */
116     PetscContainer container;
117     TSDM         oldsdm = sdm;
118     ierr = PetscInfo(dm,"Copying TSDM due to write\n");CHKERRQ(ierr);
119     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
120     ierr = PetscNewLog(dm,struct _n_TSDM,&sdm);CHKERRQ(ierr);
121     ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr);
122     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
123     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_TSDM);CHKERRQ(ierr);
124     ierr = PetscObjectCompose((PetscObject)dm,"TSDM",(PetscObject)container);CHKERRQ(ierr);
125     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
126   }
127   *tsdm = sdm;
128   PetscFunctionReturn(0);
129 }
130 
131 #undef __FUNCT__
132 #define __FUNCT__ "DMTSCopyContext"
133 /*@C
134    DMTSCopyContext - copies a DM context to a new DM
135 
136    Logically Collective
137 
138    Input Arguments:
139 +  dmsrc - DM to obtain context from
140 -  dmdest - DM to add context to
141 
142    Level: developer
143 
144    Note:
145    The context is copied by reference. This function does not ensure that a context exists.
146 
147 .seealso: DMTSGetContext(), TSSetDM()
148 @*/
149 PetscErrorCode DMTSCopyContext(DM dmsrc,DM dmdest)
150 {
151   PetscErrorCode ierr;
152   PetscContainer container;
153 
154   PetscFunctionBegin;
155   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
156   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
157   ierr = PetscObjectQuery((PetscObject)dmsrc,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
158   if (container) {
159     ierr = PetscObjectCompose((PetscObject)dmdest,"TSDM",(PetscObject)container);CHKERRQ(ierr);
160     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_TSDM,DMRestrictHook_TSDM,PETSC_NULL);CHKERRQ(ierr);
161   }
162   PetscFunctionReturn(0);
163 }
164 
165 #undef __FUNCT__
166 #define __FUNCT__ "DMTSSetIFunction"
167 /*@C
168    DMTSSetIFunction - set TS implicit function evaluation function
169 
170    Not Collective
171 
172    Input Arguments:
173 +  dm - DM to be used with TS
174 .  func - function evaluation function, see TSSetIFunction() for calling sequence
175 -  ctx - context for residual evaluation
176 
177    Level: advanced
178 
179    Note:
180    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
181    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
182    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
183 
184 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
185 @*/
186 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
187 {
188   PetscErrorCode ierr;
189   TSDM tsdm;
190 
191   PetscFunctionBegin;
192   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
193   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
194   if (func) tsdm->ifunction = func;
195   if (ctx)  tsdm->ifunctionctx = ctx;
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "DMTSGetIFunction"
201 /*@C
202    DMTSGetIFunction - get TS implicit residual evaluation function
203 
204    Not Collective
205 
206    Input Argument:
207 .  dm - DM to be used with TS
208 
209    Output Arguments:
210 +  func - function evaluation function, see TSSetIFunction() for calling sequence
211 -  ctx - context for residual evaluation
212 
213    Level: advanced
214 
215    Note:
216    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
217    associated with the DM.
218 
219 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
220 @*/
221 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
222 {
223   PetscErrorCode ierr;
224   TSDM tsdm;
225 
226   PetscFunctionBegin;
227   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
228   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
229   if (func) *func = tsdm->ifunction;
230   if (ctx)  *ctx = tsdm->ifunctionctx;
231   PetscFunctionReturn(0);
232 }
233 
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "DMTSSetRHSFunction"
237 /*@C
238    DMTSSetRHSFunction - set TS explicit residual evaluation function
239 
240    Not Collective
241 
242    Input Arguments:
243 +  dm - DM to be used with TS
244 .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
245 -  ctx - context for residual evaluation
246 
247    Level: advanced
248 
249    Note:
250    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
251    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
252    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
253 
254 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
255 @*/
256 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
257 {
258   PetscErrorCode ierr;
259   TSDM tsdm;
260 
261   PetscFunctionBegin;
262   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
263   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
264   if (func) tsdm->rhsfunction = func;
265   if (ctx)  tsdm->rhsfunctionctx = ctx;
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "DMTSGetRHSFunction"
271 /*@C
272    DMTSGetRHSFunction - get TS explicit residual evaluation function
273 
274    Not Collective
275 
276    Input Argument:
277 .  dm - DM to be used with TS
278 
279    Output Arguments:
280 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
281 -  ctx - context for residual evaluation
282 
283    Level: advanced
284 
285    Note:
286    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
287    associated with the DM.
288 
289 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
290 @*/
291 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
292 {
293   PetscErrorCode ierr;
294   TSDM tsdm;
295 
296   PetscFunctionBegin;
297   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
298   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
299   if (func) *func = tsdm->rhsfunction;
300   if (ctx)  *ctx = tsdm->rhsfunctionctx;
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "DMTSSetIJacobian"
306 /*@C
307    DMTSSetIJacobian - set TS Jacobian evaluation function
308 
309    Not Collective
310 
311    Input Argument:
312 +  dm - DM to be used with TS
313 .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
314 -  ctx - context for residual evaluation
315 
316    Level: advanced
317 
318    Note:
319    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
320    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
321    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
322 
323 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
324 @*/
325 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
326 {
327   PetscErrorCode ierr;
328   TSDM sdm;
329 
330   PetscFunctionBegin;
331   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
332   ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr);
333   if (func) sdm->ijacobian = func;
334   if (ctx)  sdm->ijacobianctx = ctx;
335   PetscFunctionReturn(0);
336 }
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "DMTSGetIJacobian"
340 /*@C
341    DMTSGetIJacobian - get TS Jacobian evaluation function
342 
343    Not Collective
344 
345    Input Argument:
346 .  dm - DM to be used with TS
347 
348    Output Arguments:
349 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
350 -  ctx - context for residual evaluation
351 
352    Level: advanced
353 
354    Note:
355    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
356    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
357    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
358 
359 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
360 @*/
361 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
362 {
363   PetscErrorCode ierr;
364   TSDM tsdm;
365 
366   PetscFunctionBegin;
367   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
368   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
369   if (func) *func = tsdm->ijacobian;
370   if (ctx)  *ctx = tsdm->ijacobianctx;
371   PetscFunctionReturn(0);
372 }
373 
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "DMTSSetRHSJacobian"
377 /*@C
378    DMTSSetRHSJacobian - set TS Jacobian evaluation function
379 
380    Not Collective
381 
382    Input Argument:
383 +  dm - DM to be used with TS
384 .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
385 -  ctx - context for residual evaluation
386 
387    Level: advanced
388 
389    Note:
390    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
391    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
392    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
393 
394 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
395 @*/
396 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
397 {
398   PetscErrorCode ierr;
399   TSDM tsdm;
400 
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
403   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
404   if (func) tsdm->rhsjacobian = func;
405   if (ctx)  tsdm->rhsjacobianctx = ctx;
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "DMTSGetRHSJacobian"
411 /*@C
412    DMTSGetRHSJacobian - get TS Jacobian evaluation function
413 
414    Not Collective
415 
416    Input Argument:
417 .  dm - DM to be used with TS
418 
419    Output Arguments:
420 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
421 -  ctx - context for residual evaluation
422 
423    Level: advanced
424 
425    Note:
426    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
427    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
428    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
429 
430 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
431 @*/
432 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
433 {
434   PetscErrorCode ierr;
435   TSDM tsdm;
436 
437   PetscFunctionBegin;
438   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
439   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
440   if (func) *func = tsdm->rhsjacobian;
441   if (ctx)  *ctx = tsdm->rhsjacobianctx;
442   PetscFunctionReturn(0);
443 }
444 
445 #if 0
446 
447 #undef __FUNCT__
448 #define __FUNCT__ "TSDefaultComputeFunction_DMLegacy"
449 static PetscErrorCode TSDefaultComputeFunction_DMLegacy(TS ts,Vec X,Vec F,void *ctx)
450 {
451   PetscErrorCode ierr;
452   DM             dm;
453 
454   PetscFunctionBegin;
455   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
456   ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr);
457   PetscFunctionReturn(0);
458 }
459 
460 #undef __FUNCT__
461 #define __FUNCT__ "TSDefaultComputeJacobian_DMLegacy"
462 static PetscErrorCode TSDefaultComputeJacobian_DMLegacy(TS ts,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
463 {
464   PetscErrorCode ierr;
465   DM             dm;
466 
467   PetscFunctionBegin;
468   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
469   ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr);
470   PetscFunctionReturn(0);
471 }
472 
473 #undef __FUNCT__
474 #define __FUNCT__ "DMTSSetUpLegacy"
475 /* Sets up calling of legacy DM routines */
476 PetscErrorCode DMTSSetUpLegacy(DM dm)
477 {
478   PetscErrorCode ierr;
479   TSDM         sdm;
480 
481   PetscFunctionBegin;
482   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
483   if (!sdm->computefunction) {ierr = DMTSSetFunction(dm,TSDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
484   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
485   if (!sdm->computejacobian) {
486     if (dm->ops->functionj) {
487       ierr = DMTSSetJacobian(dm,TSDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);
488     } else {
489       ierr = DMTSSetJacobian(dm,TSDefaultComputeJacobianColor,PETSC_NULL);CHKERRQ(ierr);
490     }
491   }
492   PetscFunctionReturn(0);
493 }
494 
495 #endif
496