xref: /petsc/src/ts/utils/dmts.c (revision 2a34c10ce821da49e32a0f08e4bd95dc44083eea)
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 DMTS 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   *tsdm = (DMTS) dm->dmts;
135   if (!*tsdm) {
136     ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr);
137     ierr = DMTSCreate(((PetscObject)dm)->comm,tsdm);CHKERRQ(ierr);
138     dm->dmts = (PetscObject) *tsdm;
139     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr);
140   }
141   PetscFunctionReturn(0);
142 }
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "DMGetDMTSWrite"
146 /*@C
147    DMGetDMTSWrite - get write access to private DMTS context from a DM
148 
149    Not Collective
150 
151    Input Argument:
152 .  dm - DM to be used with TS
153 
154    Output Argument:
155 .  tsdm - private DMTS context
156 
157    Level: developer
158 
159 .seealso: DMGetDMTS()
160 @*/
161 PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
162 {
163   PetscErrorCode ierr;
164   DMTS           sdm;
165 
166   PetscFunctionBegin;
167   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
168   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
169   if (!sdm->originaldm) sdm->originaldm = dm;
170   if (sdm->originaldm != dm) {  /* Copy on write */
171     DMTS          oldsdm = sdm;
172     ierr = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr);
173     ierr = DMTSCreate(((PetscObject)dm)->comm,&sdm);CHKERRQ(ierr);
174     ierr = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr);
175     ierr = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr);
176     dm->dmts = (PetscObject) sdm;
177   }
178   *tsdm = sdm;
179   PetscFunctionReturn(0);
180 }
181 
182 #undef __FUNCT__
183 #define __FUNCT__ "DMCopyDMTS"
184 /*@C
185    DMCopyDMTS - copies a DM context to a new DM
186 
187    Logically Collective
188 
189    Input Arguments:
190 +  dmsrc - DM to obtain context from
191 -  dmdest - DM to add context to
192 
193    Level: developer
194 
195    Note:
196    The context is copied by reference. This function does not ensure that a context exists.
197 
198 .seealso: DMGetDMTS(), TSSetDM()
199 @*/
200 PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
201 {
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
206   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
207   ierr = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr);
208   dmdest->dmts = dmsrc->dmts;
209   ierr = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr);
210   ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "DMTSSetIFunction"
216 /*@C
217    DMTSSetIFunction - set TS implicit function evaluation function
218 
219    Not Collective
220 
221    Input Arguments:
222 +  dm - DM to be used with TS
223 .  func - function evaluation function, see TSSetIFunction() for calling sequence
224 -  ctx - context for residual evaluation
225 
226    Level: advanced
227 
228    Note:
229    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
230    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
231    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
232 
233 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
234 @*/
235 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
236 {
237   PetscErrorCode ierr;
238   DMTS           tsdm;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
242   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
243   if (func) tsdm->ops->ifunction = func;
244   if (ctx)  tsdm->ifunctionctx = ctx;
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "DMTSGetIFunction"
250 /*@C
251    DMTSGetIFunction - get TS implicit residual evaluation function
252 
253    Not Collective
254 
255    Input Argument:
256 .  dm - DM to be used with TS
257 
258    Output Arguments:
259 +  func - function evaluation function, see TSSetIFunction() for calling sequence
260 -  ctx - context for residual evaluation
261 
262    Level: advanced
263 
264    Note:
265    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
266    associated with the DM.
267 
268 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
269 @*/
270 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
271 {
272   PetscErrorCode ierr;
273   DMTS           tsdm;
274 
275   PetscFunctionBegin;
276   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
277   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
278   if (func) *func = tsdm->ops->ifunction;
279   if (ctx)  *ctx = tsdm->ifunctionctx;
280   PetscFunctionReturn(0);
281 }
282 
283 
284 #undef __FUNCT__
285 #define __FUNCT__ "DMTSSetRHSFunction"
286 /*@C
287    DMTSSetRHSFunction - set TS explicit residual evaluation function
288 
289    Not Collective
290 
291    Input Arguments:
292 +  dm - DM to be used with TS
293 .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
294 -  ctx - context for residual evaluation
295 
296    Level: advanced
297 
298    Note:
299    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
300    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
301    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
302 
303 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
304 @*/
305 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
306 {
307   PetscErrorCode ierr;
308   DMTS           tsdm;
309 
310   PetscFunctionBegin;
311   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
312   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
313   if (func) tsdm->ops->rhsfunction = func;
314   if (ctx)  tsdm->rhsfunctionctx = ctx;
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "DMTSGetSolutionFunction"
320 /*@C
321    DMTSGetSolutionFunction - gets the TS solution evaluation function
322 
323    Not Collective
324 
325    Input Arguments:
326 .  dm - DM to be used with TS
327 
328    Output Parameters:
329 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
330 -  ctx - context for solution evaluation
331 
332    Level: advanced
333 
334 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
335 @*/
336 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
337 {
338   PetscErrorCode ierr;
339   DMTS           tsdm;
340 
341   PetscFunctionBegin;
342   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
343   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
344   if (func) *func = tsdm->ops->solution;
345   if (ctx)  *ctx  = tsdm->solutionctx;
346   PetscFunctionReturn(0);
347 }
348 
349 #undef __FUNCT__
350 #define __FUNCT__ "DMTSSetSolutionFunction"
351 /*@C
352    DMTSSetSolutionFunction - set TS solution evaluation function
353 
354    Not Collective
355 
356    Input Arguments:
357 +  dm - DM to be used with TS
358 .  func - solution function evaluation function, see TSSetSolution() for calling sequence
359 -  ctx - context for solution evaluation
360 
361    Level: advanced
362 
363    Note:
364    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
365    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
366    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
367 
368 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
369 @*/
370 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
371 {
372   PetscErrorCode ierr;
373   DMTS           tsdm;
374 
375   PetscFunctionBegin;
376   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
377   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
378   if (func) tsdm->ops->solution = func;
379   if (ctx)  tsdm->solutionctx   = ctx;
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "DMTSGetRHSFunction"
385 /*@C
386    DMTSGetRHSFunction - get TS explicit residual evaluation function
387 
388    Not Collective
389 
390    Input Argument:
391 .  dm - DM to be used with TS
392 
393    Output Arguments:
394 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
395 -  ctx - context for residual evaluation
396 
397    Level: advanced
398 
399    Note:
400    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
401    associated with the DM.
402 
403 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
404 @*/
405 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
406 {
407   PetscErrorCode ierr;
408   DMTS           tsdm;
409 
410   PetscFunctionBegin;
411   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
412   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
413   if (func) *func = tsdm->ops->rhsfunction;
414   if (ctx)  *ctx = tsdm->rhsfunctionctx;
415   PetscFunctionReturn(0);
416 }
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "DMTSSetIJacobian"
420 /*@C
421    DMTSSetIJacobian - set TS Jacobian evaluation function
422 
423    Not Collective
424 
425    Input Argument:
426 +  dm - DM to be used with TS
427 .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
428 -  ctx - context for residual evaluation
429 
430    Level: advanced
431 
432    Note:
433    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
434    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
435    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
436 
437 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
438 @*/
439 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
440 {
441   PetscErrorCode ierr;
442   DMTS           sdm;
443 
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
446   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
447   if (func) sdm->ops->ijacobian = func;
448   if (ctx)  sdm->ijacobianctx   = ctx;
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "DMTSGetIJacobian"
454 /*@C
455    DMTSGetIJacobian - get TS Jacobian evaluation function
456 
457    Not Collective
458 
459    Input Argument:
460 .  dm - DM to be used with TS
461 
462    Output Arguments:
463 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
464 -  ctx - context for residual evaluation
465 
466    Level: advanced
467 
468    Note:
469    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
470    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
471    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
472 
473 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
474 @*/
475 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
476 {
477   PetscErrorCode ierr;
478   DMTS           tsdm;
479 
480   PetscFunctionBegin;
481   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
482   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
483   if (func) *func = tsdm->ops->ijacobian;
484   if (ctx)  *ctx = tsdm->ijacobianctx;
485   PetscFunctionReturn(0);
486 }
487 
488 
489 #undef __FUNCT__
490 #define __FUNCT__ "DMTSSetRHSJacobian"
491 /*@C
492    DMTSSetRHSJacobian - set TS Jacobian evaluation function
493 
494    Not Collective
495 
496    Input Argument:
497 +  dm - DM to be used with TS
498 .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
499 -  ctx - context for residual evaluation
500 
501    Level: advanced
502 
503    Note:
504    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
505    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
506    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
507 
508 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
509 @*/
510 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
511 {
512   PetscErrorCode ierr;
513   DMTS           tsdm;
514 
515   PetscFunctionBegin;
516   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
517   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
518   if (func) tsdm->ops->rhsjacobian = func;
519   if (ctx)  tsdm->rhsjacobianctx = ctx;
520   PetscFunctionReturn(0);
521 }
522 
523 #undef __FUNCT__
524 #define __FUNCT__ "DMTSGetRHSJacobian"
525 /*@C
526    DMTSGetRHSJacobian - get TS Jacobian evaluation function
527 
528    Not Collective
529 
530    Input Argument:
531 .  dm - DM to be used with TS
532 
533    Output Arguments:
534 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
535 -  ctx - context for residual evaluation
536 
537    Level: advanced
538 
539    Note:
540    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
541    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
542    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
543 
544 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
545 @*/
546 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
547 {
548   PetscErrorCode ierr;
549   DMTS           tsdm;
550 
551   PetscFunctionBegin;
552   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
553   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
554   if (func) *func = tsdm->ops->rhsjacobian;
555   if (ctx)  *ctx = tsdm->rhsjacobianctx;
556   PetscFunctionReturn(0);
557 }
558