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