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