xref: /petsc/src/dm/interface/dm.c (revision 47c6ae997ffd1b2afd66b6474dff5950ae8613d1)
1*47c6ae99SBarry Smith #define PETSCDM_DLL
2*47c6ae99SBarry Smith 
3*47c6ae99SBarry Smith #include "private/dmimpl.h"     /*I      "petscda.h"     I*/
4*47c6ae99SBarry Smith 
5*47c6ae99SBarry Smith /*
6*47c6ae99SBarry Smith    Provides an interface for functionality needed by the DAMG routines.
7*47c6ae99SBarry Smith    Currently this interface is supported by the DA and DMComposite objects
8*47c6ae99SBarry Smith 
9*47c6ae99SBarry Smith    Note: this is actually no such thing as a DM object, rather it is
10*47c6ae99SBarry Smith    the common set of functions shared by DA and DMComposite.
11*47c6ae99SBarry Smith 
12*47c6ae99SBarry Smith */
13*47c6ae99SBarry Smith 
14*47c6ae99SBarry Smith #undef __FUNCT__
15*47c6ae99SBarry Smith #define __FUNCT__ "DMDestroy"
16*47c6ae99SBarry Smith /*@
17*47c6ae99SBarry Smith     DMDestroy - Destroys a vector packer or DA.
18*47c6ae99SBarry Smith 
19*47c6ae99SBarry Smith     Collective on DM
20*47c6ae99SBarry Smith 
21*47c6ae99SBarry Smith     Input Parameter:
22*47c6ae99SBarry Smith .   dm - the DM object to destroy
23*47c6ae99SBarry Smith 
24*47c6ae99SBarry Smith     Level: developer
25*47c6ae99SBarry Smith 
26*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
27*47c6ae99SBarry Smith 
28*47c6ae99SBarry Smith @*/
29*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDestroy(DM dm)
30*47c6ae99SBarry Smith {
31*47c6ae99SBarry Smith   PetscErrorCode ierr;
32*47c6ae99SBarry Smith 
33*47c6ae99SBarry Smith   PetscFunctionBegin;
34*47c6ae99SBarry Smith   ierr = (*((PetscObject)dm)->bops->destroy)((PetscObject)dm);CHKERRQ(ierr);
35*47c6ae99SBarry Smith   PetscFunctionReturn(0);
36*47c6ae99SBarry Smith }
37*47c6ae99SBarry Smith 
38*47c6ae99SBarry Smith #undef __FUNCT__
39*47c6ae99SBarry Smith #define __FUNCT__ "DMView"
40*47c6ae99SBarry Smith /*@
41*47c6ae99SBarry Smith     DMView - Views a vector packer or DA.
42*47c6ae99SBarry Smith 
43*47c6ae99SBarry Smith     Collective on DM
44*47c6ae99SBarry Smith 
45*47c6ae99SBarry Smith     Input Parameter:
46*47c6ae99SBarry Smith +   dm - the DM object to view
47*47c6ae99SBarry Smith -   v - the viewer
48*47c6ae99SBarry Smith 
49*47c6ae99SBarry Smith     Level: developer
50*47c6ae99SBarry Smith 
51*47c6ae99SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
52*47c6ae99SBarry Smith 
53*47c6ae99SBarry Smith @*/
54*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMView(DM dm,PetscViewer v)
55*47c6ae99SBarry Smith {
56*47c6ae99SBarry Smith   PetscErrorCode ierr;
57*47c6ae99SBarry Smith 
58*47c6ae99SBarry Smith   PetscFunctionBegin;
59*47c6ae99SBarry Smith   if (((PetscObject)dm)->bops->view) {
60*47c6ae99SBarry Smith     ierr = (*((PetscObject)dm)->bops->view)((PetscObject)dm,v);CHKERRQ(ierr);
61*47c6ae99SBarry Smith   }
62*47c6ae99SBarry Smith   PetscFunctionReturn(0);
63*47c6ae99SBarry Smith }
64*47c6ae99SBarry Smith 
65*47c6ae99SBarry Smith #undef __FUNCT__
66*47c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector"
67*47c6ae99SBarry Smith /*@
68*47c6ae99SBarry Smith     DMCreateGlobalVector - Creates a global vector from a DA or DMComposite object
69*47c6ae99SBarry Smith 
70*47c6ae99SBarry Smith     Collective on DM
71*47c6ae99SBarry Smith 
72*47c6ae99SBarry Smith     Input Parameter:
73*47c6ae99SBarry Smith .   dm - the DM object
74*47c6ae99SBarry Smith 
75*47c6ae99SBarry Smith     Output Parameter:
76*47c6ae99SBarry Smith .   vec - the global vector
77*47c6ae99SBarry Smith 
78*47c6ae99SBarry Smith     Level: developer
79*47c6ae99SBarry Smith 
80*47c6ae99SBarry Smith .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
81*47c6ae99SBarry Smith 
82*47c6ae99SBarry Smith @*/
83*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector(DM dm,Vec *vec)
84*47c6ae99SBarry Smith {
85*47c6ae99SBarry Smith   PetscErrorCode ierr;
86*47c6ae99SBarry Smith 
87*47c6ae99SBarry Smith   PetscFunctionBegin;
88*47c6ae99SBarry Smith   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
89*47c6ae99SBarry Smith   PetscFunctionReturn(0);
90*47c6ae99SBarry Smith }
91*47c6ae99SBarry Smith 
92*47c6ae99SBarry Smith #undef __FUNCT__
93*47c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector"
94*47c6ae99SBarry Smith /*@
95*47c6ae99SBarry Smith     DMCreateLocalVector - Creates a local vector from a DA or DMComposite object
96*47c6ae99SBarry Smith 
97*47c6ae99SBarry Smith     Not Collective
98*47c6ae99SBarry Smith 
99*47c6ae99SBarry Smith     Input Parameter:
100*47c6ae99SBarry Smith .   dm - the DM object
101*47c6ae99SBarry Smith 
102*47c6ae99SBarry Smith     Output Parameter:
103*47c6ae99SBarry Smith .   vec - the local vector
104*47c6ae99SBarry Smith 
105*47c6ae99SBarry Smith     Level: developer
106*47c6ae99SBarry Smith 
107*47c6ae99SBarry Smith .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
108*47c6ae99SBarry Smith 
109*47c6ae99SBarry Smith @*/
110*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector(DM dm,Vec *vec)
111*47c6ae99SBarry Smith {
112*47c6ae99SBarry Smith   PetscErrorCode ierr;
113*47c6ae99SBarry Smith 
114*47c6ae99SBarry Smith   PetscFunctionBegin;
115*47c6ae99SBarry Smith   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
116*47c6ae99SBarry Smith   PetscFunctionReturn(0);
117*47c6ae99SBarry Smith }
118*47c6ae99SBarry Smith 
119*47c6ae99SBarry Smith #undef __FUNCT__
120*47c6ae99SBarry Smith #define __FUNCT__ "DMGetInterpolation"
121*47c6ae99SBarry Smith /*@
122*47c6ae99SBarry Smith     DMGetInterpolation - Gets interpolation matrix between two DA or DMComposite objects
123*47c6ae99SBarry Smith 
124*47c6ae99SBarry Smith     Collective on DM
125*47c6ae99SBarry Smith 
126*47c6ae99SBarry Smith     Input Parameter:
127*47c6ae99SBarry Smith +   dm1 - the DM object
128*47c6ae99SBarry Smith -   dm2 - the second, finer DM object
129*47c6ae99SBarry Smith 
130*47c6ae99SBarry Smith     Output Parameter:
131*47c6ae99SBarry Smith +  mat - the interpolation
132*47c6ae99SBarry Smith -  vec - the scaling (optional)
133*47c6ae99SBarry Smith 
134*47c6ae99SBarry Smith     Level: developer
135*47c6ae99SBarry Smith 
136*47c6ae99SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix()
137*47c6ae99SBarry Smith 
138*47c6ae99SBarry Smith @*/
139*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
140*47c6ae99SBarry Smith {
141*47c6ae99SBarry Smith   PetscErrorCode ierr;
142*47c6ae99SBarry Smith 
143*47c6ae99SBarry Smith   PetscFunctionBegin;
144*47c6ae99SBarry Smith   ierr = (*dm1->ops->getinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
145*47c6ae99SBarry Smith   PetscFunctionReturn(0);
146*47c6ae99SBarry Smith }
147*47c6ae99SBarry Smith 
148*47c6ae99SBarry Smith #undef __FUNCT__
149*47c6ae99SBarry Smith #define __FUNCT__ "DMGetInjection"
150*47c6ae99SBarry Smith /*@
151*47c6ae99SBarry Smith     DMGetInjection - Gets injection matrix between two DA or DMComposite objects
152*47c6ae99SBarry Smith 
153*47c6ae99SBarry Smith     Collective on DM
154*47c6ae99SBarry Smith 
155*47c6ae99SBarry Smith     Input Parameter:
156*47c6ae99SBarry Smith +   dm1 - the DM object
157*47c6ae99SBarry Smith -   dm2 - the second, finer DM object
158*47c6ae99SBarry Smith 
159*47c6ae99SBarry Smith     Output Parameter:
160*47c6ae99SBarry Smith .   ctx - the injection
161*47c6ae99SBarry Smith 
162*47c6ae99SBarry Smith     Level: developer
163*47c6ae99SBarry Smith 
164*47c6ae99SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMGetInterpolation()
165*47c6ae99SBarry Smith 
166*47c6ae99SBarry Smith @*/
167*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInjection(DM dm1,DM dm2,VecScatter *ctx)
168*47c6ae99SBarry Smith {
169*47c6ae99SBarry Smith   PetscErrorCode ierr;
170*47c6ae99SBarry Smith 
171*47c6ae99SBarry Smith   PetscFunctionBegin;
172*47c6ae99SBarry Smith   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
173*47c6ae99SBarry Smith   PetscFunctionReturn(0);
174*47c6ae99SBarry Smith }
175*47c6ae99SBarry Smith 
176*47c6ae99SBarry Smith #undef __FUNCT__
177*47c6ae99SBarry Smith #define __FUNCT__ "DMGetColoring"
178*47c6ae99SBarry Smith /*@
179*47c6ae99SBarry Smith     DMGetColoring - Gets coloring for a DA or DMComposite
180*47c6ae99SBarry Smith 
181*47c6ae99SBarry Smith     Collective on DM
182*47c6ae99SBarry Smith 
183*47c6ae99SBarry Smith     Input Parameter:
184*47c6ae99SBarry Smith +   dm - the DM object
185*47c6ae99SBarry Smith .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
186*47c6ae99SBarry Smith -   matype - either MATAIJ or MATBAIJ
187*47c6ae99SBarry Smith 
188*47c6ae99SBarry Smith     Output Parameter:
189*47c6ae99SBarry Smith .   coloring - the coloring
190*47c6ae99SBarry Smith 
191*47c6ae99SBarry Smith     Level: developer
192*47c6ae99SBarry Smith 
193*47c6ae99SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
194*47c6ae99SBarry Smith 
195*47c6ae99SBarry Smith @*/
196*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
197*47c6ae99SBarry Smith {
198*47c6ae99SBarry Smith   PetscErrorCode ierr;
199*47c6ae99SBarry Smith 
200*47c6ae99SBarry Smith   PetscFunctionBegin;
201*47c6ae99SBarry Smith   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
202*47c6ae99SBarry Smith   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
203*47c6ae99SBarry Smith   PetscFunctionReturn(0);
204*47c6ae99SBarry Smith }
205*47c6ae99SBarry Smith 
206*47c6ae99SBarry Smith #undef __FUNCT__
207*47c6ae99SBarry Smith #define __FUNCT__ "DMGetMatrix"
208*47c6ae99SBarry Smith /*@C
209*47c6ae99SBarry Smith     DMGetMatrix - Gets empty Jacobian for a DA or DMComposite
210*47c6ae99SBarry Smith 
211*47c6ae99SBarry Smith     Collective on DM
212*47c6ae99SBarry Smith 
213*47c6ae99SBarry Smith     Input Parameter:
214*47c6ae99SBarry Smith +   dm - the DM object
215*47c6ae99SBarry Smith -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
216*47c6ae99SBarry Smith             any type which inherits from one of these (such as MATAIJ, MATLUSOL, etc.).
217*47c6ae99SBarry Smith 
218*47c6ae99SBarry Smith     Output Parameter:
219*47c6ae99SBarry Smith .   mat - the empty Jacobian
220*47c6ae99SBarry Smith 
221*47c6ae99SBarry Smith     Level: developer
222*47c6ae99SBarry Smith 
223*47c6ae99SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
224*47c6ae99SBarry Smith 
225*47c6ae99SBarry Smith @*/
226*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix(DM dm, const MatType mtype,Mat *mat)
227*47c6ae99SBarry Smith {
228*47c6ae99SBarry Smith   PetscErrorCode ierr;
229*47c6ae99SBarry Smith 
230*47c6ae99SBarry Smith   PetscFunctionBegin;
231*47c6ae99SBarry Smith   ierr = (*dm->ops->getmatrix)(dm,mtype,mat);CHKERRQ(ierr);
232*47c6ae99SBarry Smith   PetscFunctionReturn(0);
233*47c6ae99SBarry Smith }
234*47c6ae99SBarry Smith 
235*47c6ae99SBarry Smith #undef __FUNCT__
236*47c6ae99SBarry Smith #define __FUNCT__ "DMRefine"
237*47c6ae99SBarry Smith /*@
238*47c6ae99SBarry Smith     DMRefine - Refines a DM object
239*47c6ae99SBarry Smith 
240*47c6ae99SBarry Smith     Collective on DM
241*47c6ae99SBarry Smith 
242*47c6ae99SBarry Smith     Input Parameter:
243*47c6ae99SBarry Smith +   dm - the DM object
244*47c6ae99SBarry Smith -   comm - the communicator to contain the new DM object (or PETSC_NULL)
245*47c6ae99SBarry Smith 
246*47c6ae99SBarry Smith     Output Parameter:
247*47c6ae99SBarry Smith .   dmf - the refined DM
248*47c6ae99SBarry Smith 
249*47c6ae99SBarry Smith     Level: developer
250*47c6ae99SBarry Smith 
251*47c6ae99SBarry Smith .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
252*47c6ae99SBarry Smith 
253*47c6ae99SBarry Smith @*/
254*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefine(DM dm,MPI_Comm comm,DM *dmf)
255*47c6ae99SBarry Smith {
256*47c6ae99SBarry Smith   PetscErrorCode ierr;
257*47c6ae99SBarry Smith 
258*47c6ae99SBarry Smith   PetscFunctionBegin;
259*47c6ae99SBarry Smith   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
260*47c6ae99SBarry Smith   PetscFunctionReturn(0);
261*47c6ae99SBarry Smith }
262*47c6ae99SBarry Smith 
263*47c6ae99SBarry Smith #undef __FUNCT__
264*47c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin"
265*47c6ae99SBarry Smith /*@
266*47c6ae99SBarry Smith     DMGlobalToLocalBegin - Begins updating local vectors from global vector
267*47c6ae99SBarry Smith 
268*47c6ae99SBarry Smith     Neighbor-wise Collective on DM
269*47c6ae99SBarry Smith 
270*47c6ae99SBarry Smith     Input Parameters:
271*47c6ae99SBarry Smith +   dm - the DM object
272*47c6ae99SBarry Smith .   g - the global vector
273*47c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
274*47c6ae99SBarry Smith -   l - the local vector
275*47c6ae99SBarry Smith 
276*47c6ae99SBarry Smith 
277*47c6ae99SBarry Smith     Level: beginner
278*47c6ae99SBarry Smith 
279*47c6ae99SBarry Smith .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobal()
280*47c6ae99SBarry Smith 
281*47c6ae99SBarry Smith @*/
282*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
283*47c6ae99SBarry Smith {
284*47c6ae99SBarry Smith   PetscErrorCode ierr;
285*47c6ae99SBarry Smith 
286*47c6ae99SBarry Smith   PetscFunctionBegin;
287*47c6ae99SBarry Smith   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr);
288*47c6ae99SBarry Smith   PetscFunctionReturn(0);
289*47c6ae99SBarry Smith }
290*47c6ae99SBarry Smith 
291*47c6ae99SBarry Smith #undef __FUNCT__
292*47c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd"
293*47c6ae99SBarry Smith /*@
294*47c6ae99SBarry Smith     DMGlobalToLocalEnd - Ends updating local vectors from global vector
295*47c6ae99SBarry Smith 
296*47c6ae99SBarry Smith     Neighbor-wise Collective on DM
297*47c6ae99SBarry Smith 
298*47c6ae99SBarry Smith     Input Parameters:
299*47c6ae99SBarry Smith +   dm - the DM object
300*47c6ae99SBarry Smith .   g - the global vector
301*47c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
302*47c6ae99SBarry Smith -   l - the local vector
303*47c6ae99SBarry Smith 
304*47c6ae99SBarry Smith 
305*47c6ae99SBarry Smith     Level: beginner
306*47c6ae99SBarry Smith 
307*47c6ae99SBarry Smith .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobal()
308*47c6ae99SBarry Smith 
309*47c6ae99SBarry Smith @*/
310*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
311*47c6ae99SBarry Smith {
312*47c6ae99SBarry Smith   PetscErrorCode ierr;
313*47c6ae99SBarry Smith 
314*47c6ae99SBarry Smith   PetscFunctionBegin;
315*47c6ae99SBarry Smith   ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr);
316*47c6ae99SBarry Smith   PetscFunctionReturn(0);
317*47c6ae99SBarry Smith }
318*47c6ae99SBarry Smith 
319*47c6ae99SBarry Smith #undef __FUNCT__
320*47c6ae99SBarry Smith #define __FUNCT__ "DMLocalToGlobal"
321*47c6ae99SBarry Smith /*@
322*47c6ae99SBarry Smith     DMLocalToGlobal - updates global vectors from local vectors
323*47c6ae99SBarry Smith 
324*47c6ae99SBarry Smith     Neighbor-wise Collective on DM
325*47c6ae99SBarry Smith 
326*47c6ae99SBarry Smith     Input Parameters:
327*47c6ae99SBarry Smith +   dm - the DM object
328*47c6ae99SBarry Smith .   g - the global vector
329*47c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
330*47c6ae99SBarry Smith -   l - the local vector
331*47c6ae99SBarry Smith 
332*47c6ae99SBarry Smith 
333*47c6ae99SBarry Smith     Level: beginner
334*47c6ae99SBarry Smith 
335*47c6ae99SBarry Smith .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
336*47c6ae99SBarry Smith 
337*47c6ae99SBarry Smith @*/
338*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMLocalToGlobal(DM dm,Vec g,InsertMode mode,Vec l)
339*47c6ae99SBarry Smith {
340*47c6ae99SBarry Smith   PetscErrorCode ierr;
341*47c6ae99SBarry Smith 
342*47c6ae99SBarry Smith   PetscFunctionBegin;
343*47c6ae99SBarry Smith   ierr = (*dm->ops->localtoglobal)(dm,g,mode,l);CHKERRQ(ierr);
344*47c6ae99SBarry Smith   PetscFunctionReturn(0);
345*47c6ae99SBarry Smith }
346*47c6ae99SBarry Smith 
347*47c6ae99SBarry Smith #undef __FUNCT__
348*47c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobianDefault"
349*47c6ae99SBarry Smith /*@
350*47c6ae99SBarry Smith     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
351*47c6ae99SBarry Smith 
352*47c6ae99SBarry Smith     Collective on DM
353*47c6ae99SBarry Smith 
354*47c6ae99SBarry Smith     Input Parameter:
355*47c6ae99SBarry Smith +   dm - the DM object
356*47c6ae99SBarry Smith .   x - location to compute Jacobian at; may be ignored for linear problems
357*47c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
358*47c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
359*47c6ae99SBarry Smith 
360*47c6ae99SBarry Smith     Level: developer
361*47c6ae99SBarry Smith 
362*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
363*47c6ae99SBarry Smith          DMSetFunction()
364*47c6ae99SBarry Smith 
365*47c6ae99SBarry Smith @*/
366*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
367*47c6ae99SBarry Smith {
368*47c6ae99SBarry Smith   PetscErrorCode ierr;
369*47c6ae99SBarry Smith   PetscFunctionBegin;
370*47c6ae99SBarry Smith   *stflag = SAME_NONZERO_PATTERN;
371*47c6ae99SBarry Smith   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
372*47c6ae99SBarry Smith   if (A != B) {
373*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
374*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
375*47c6ae99SBarry Smith   }
376*47c6ae99SBarry Smith   PetscFunctionReturn(0);
377*47c6ae99SBarry Smith }
378*47c6ae99SBarry Smith 
379*47c6ae99SBarry Smith #undef __FUNCT__
380*47c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen"
381*47c6ae99SBarry Smith /*@
382*47c6ae99SBarry Smith     DMCoarsen - Coarsens a DM object
383*47c6ae99SBarry Smith 
384*47c6ae99SBarry Smith     Collective on DM
385*47c6ae99SBarry Smith 
386*47c6ae99SBarry Smith     Input Parameter:
387*47c6ae99SBarry Smith +   dm - the DM object
388*47c6ae99SBarry Smith -   comm - the communicator to contain the new DM object (or PETSC_NULL)
389*47c6ae99SBarry Smith 
390*47c6ae99SBarry Smith     Output Parameter:
391*47c6ae99SBarry Smith .   dmc - the coarsened DM
392*47c6ae99SBarry Smith 
393*47c6ae99SBarry Smith     Level: developer
394*47c6ae99SBarry Smith 
395*47c6ae99SBarry Smith .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
396*47c6ae99SBarry Smith 
397*47c6ae99SBarry Smith @*/
398*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
399*47c6ae99SBarry Smith {
400*47c6ae99SBarry Smith   PetscErrorCode ierr;
401*47c6ae99SBarry Smith 
402*47c6ae99SBarry Smith   PetscFunctionBegin;
403*47c6ae99SBarry Smith   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
404*47c6ae99SBarry Smith   (*dmc)->ops->initialguess = dm->ops->initialguess;
405*47c6ae99SBarry Smith   (*dmc)->ops->function     = dm->ops->function;
406*47c6ae99SBarry Smith   (*dmc)->ops->functionj    = dm->ops->functionj;
407*47c6ae99SBarry Smith   if (dm->ops->jacobian != DMComputeJacobianDefault) {
408*47c6ae99SBarry Smith     (*dmc)->ops->jacobian     = dm->ops->jacobian;
409*47c6ae99SBarry Smith   }
410*47c6ae99SBarry Smith   PetscFunctionReturn(0);
411*47c6ae99SBarry Smith }
412*47c6ae99SBarry Smith 
413*47c6ae99SBarry Smith #undef __FUNCT__
414*47c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy"
415*47c6ae99SBarry Smith /*@C
416*47c6ae99SBarry Smith     DMRefineHierarchy - Refines a DM object, all levels at once
417*47c6ae99SBarry Smith 
418*47c6ae99SBarry Smith     Collective on DM
419*47c6ae99SBarry Smith 
420*47c6ae99SBarry Smith     Input Parameter:
421*47c6ae99SBarry Smith +   dm - the DM object
422*47c6ae99SBarry Smith -   nlevels - the number of levels of refinement
423*47c6ae99SBarry Smith 
424*47c6ae99SBarry Smith     Output Parameter:
425*47c6ae99SBarry Smith .   dmf - the refined DM hierarchy
426*47c6ae99SBarry Smith 
427*47c6ae99SBarry Smith     Level: developer
428*47c6ae99SBarry Smith 
429*47c6ae99SBarry Smith .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
430*47c6ae99SBarry Smith 
431*47c6ae99SBarry Smith @*/
432*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
433*47c6ae99SBarry Smith {
434*47c6ae99SBarry Smith   PetscErrorCode ierr;
435*47c6ae99SBarry Smith 
436*47c6ae99SBarry Smith   PetscFunctionBegin;
437*47c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
438*47c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
439*47c6ae99SBarry Smith   if (dm->ops->refinehierarchy) {
440*47c6ae99SBarry Smith     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
441*47c6ae99SBarry Smith   } else if (dm->ops->refine) {
442*47c6ae99SBarry Smith     PetscInt i;
443*47c6ae99SBarry Smith 
444*47c6ae99SBarry Smith     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
445*47c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
446*47c6ae99SBarry Smith       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
447*47c6ae99SBarry Smith     }
448*47c6ae99SBarry Smith   } else {
449*47c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
450*47c6ae99SBarry Smith   }
451*47c6ae99SBarry Smith   PetscFunctionReturn(0);
452*47c6ae99SBarry Smith }
453*47c6ae99SBarry Smith 
454*47c6ae99SBarry Smith #undef __FUNCT__
455*47c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy"
456*47c6ae99SBarry Smith /*@C
457*47c6ae99SBarry Smith     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
458*47c6ae99SBarry Smith 
459*47c6ae99SBarry Smith     Collective on DM
460*47c6ae99SBarry Smith 
461*47c6ae99SBarry Smith     Input Parameter:
462*47c6ae99SBarry Smith +   dm - the DM object
463*47c6ae99SBarry Smith -   nlevels - the number of levels of coarsening
464*47c6ae99SBarry Smith 
465*47c6ae99SBarry Smith     Output Parameter:
466*47c6ae99SBarry Smith .   dmc - the coarsened DM hierarchy
467*47c6ae99SBarry Smith 
468*47c6ae99SBarry Smith     Level: developer
469*47c6ae99SBarry Smith 
470*47c6ae99SBarry Smith .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
471*47c6ae99SBarry Smith 
472*47c6ae99SBarry Smith @*/
473*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
474*47c6ae99SBarry Smith {
475*47c6ae99SBarry Smith   PetscErrorCode ierr;
476*47c6ae99SBarry Smith 
477*47c6ae99SBarry Smith   PetscFunctionBegin;
478*47c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
479*47c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
480*47c6ae99SBarry Smith   PetscValidPointer(dmc,3);
481*47c6ae99SBarry Smith   if (dm->ops->coarsenhierarchy) {
482*47c6ae99SBarry Smith     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
483*47c6ae99SBarry Smith   } else if (dm->ops->coarsen) {
484*47c6ae99SBarry Smith     PetscInt i;
485*47c6ae99SBarry Smith 
486*47c6ae99SBarry Smith     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
487*47c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
488*47c6ae99SBarry Smith       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
489*47c6ae99SBarry Smith     }
490*47c6ae99SBarry Smith   } else {
491*47c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
492*47c6ae99SBarry Smith   }
493*47c6ae99SBarry Smith   PetscFunctionReturn(0);
494*47c6ae99SBarry Smith }
495*47c6ae99SBarry Smith 
496*47c6ae99SBarry Smith #undef __FUNCT__
497*47c6ae99SBarry Smith #define __FUNCT__ "DMGetAggregates"
498*47c6ae99SBarry Smith /*@
499*47c6ae99SBarry Smith    DMGetAggregates - Gets the aggregates that map between
500*47c6ae99SBarry Smith    grids associated with two DMs.
501*47c6ae99SBarry Smith 
502*47c6ae99SBarry Smith    Collective on DM
503*47c6ae99SBarry Smith 
504*47c6ae99SBarry Smith    Input Parameters:
505*47c6ae99SBarry Smith +  dmc - the coarse grid DM
506*47c6ae99SBarry Smith -  dmf - the fine grid DM
507*47c6ae99SBarry Smith 
508*47c6ae99SBarry Smith    Output Parameters:
509*47c6ae99SBarry Smith .  rest - the restriction matrix (transpose of the projection matrix)
510*47c6ae99SBarry Smith 
511*47c6ae99SBarry Smith    Level: intermediate
512*47c6ae99SBarry Smith 
513*47c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid
514*47c6ae99SBarry Smith 
515*47c6ae99SBarry Smith .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation()
516*47c6ae99SBarry Smith @*/
517*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetAggregates(DM dmc, DM dmf, Mat *rest)
518*47c6ae99SBarry Smith {
519*47c6ae99SBarry Smith   PetscErrorCode ierr;
520*47c6ae99SBarry Smith 
521*47c6ae99SBarry Smith   PetscFunctionBegin;
522*47c6ae99SBarry Smith   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
523*47c6ae99SBarry Smith   PetscFunctionReturn(0);
524*47c6ae99SBarry Smith }
525*47c6ae99SBarry Smith 
526*47c6ae99SBarry Smith #undef __FUNCT__
527*47c6ae99SBarry Smith #define __FUNCT__ "DMSetContext"
528*47c6ae99SBarry Smith /*@
529*47c6ae99SBarry Smith     DMSetContext - Set a user context into a DM object
530*47c6ae99SBarry Smith 
531*47c6ae99SBarry Smith     Not Collective
532*47c6ae99SBarry Smith 
533*47c6ae99SBarry Smith     Input Parameters:
534*47c6ae99SBarry Smith +   dm - the DM object
535*47c6ae99SBarry Smith -   ctx - the user context
536*47c6ae99SBarry Smith 
537*47c6ae99SBarry Smith     Level: intermediate
538*47c6ae99SBarry Smith 
539*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
540*47c6ae99SBarry Smith 
541*47c6ae99SBarry Smith @*/
542*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetContext(DM dm,void *ctx)
543*47c6ae99SBarry Smith {
544*47c6ae99SBarry Smith   PetscFunctionBegin;
545*47c6ae99SBarry Smith   dm->ctx = ctx;
546*47c6ae99SBarry Smith   PetscFunctionReturn(0);
547*47c6ae99SBarry Smith }
548*47c6ae99SBarry Smith 
549*47c6ae99SBarry Smith #undef __FUNCT__
550*47c6ae99SBarry Smith #define __FUNCT__ "DMGetContext"
551*47c6ae99SBarry Smith /*@
552*47c6ae99SBarry Smith     DMGetContext - Gets a user context from a DM object
553*47c6ae99SBarry Smith 
554*47c6ae99SBarry Smith     Not Collective
555*47c6ae99SBarry Smith 
556*47c6ae99SBarry Smith     Input Parameter:
557*47c6ae99SBarry Smith .   dm - the DM object
558*47c6ae99SBarry Smith 
559*47c6ae99SBarry Smith     Output Parameter:
560*47c6ae99SBarry Smith .   ctx - the user context
561*47c6ae99SBarry Smith 
562*47c6ae99SBarry Smith     Level: intermediate
563*47c6ae99SBarry Smith 
564*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
565*47c6ae99SBarry Smith 
566*47c6ae99SBarry Smith @*/
567*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetContext(DM dm,void **ctx)
568*47c6ae99SBarry Smith {
569*47c6ae99SBarry Smith   PetscFunctionBegin;
570*47c6ae99SBarry Smith   *ctx = dm->ctx;
571*47c6ae99SBarry Smith   PetscFunctionReturn(0);
572*47c6ae99SBarry Smith }
573*47c6ae99SBarry Smith 
574*47c6ae99SBarry Smith #undef __FUNCT__
575*47c6ae99SBarry Smith #define __FUNCT__ "DMSetInitialGuess"
576*47c6ae99SBarry Smith /*@
577*47c6ae99SBarry Smith     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
578*47c6ae99SBarry Smith 
579*47c6ae99SBarry Smith     Logically Collective on DM
580*47c6ae99SBarry Smith 
581*47c6ae99SBarry Smith     Input Parameter:
582*47c6ae99SBarry Smith +   dm - the DM object to destroy
583*47c6ae99SBarry Smith -   f - the function to compute the initial guess
584*47c6ae99SBarry Smith 
585*47c6ae99SBarry Smith     Level: intermediate
586*47c6ae99SBarry Smith 
587*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
588*47c6ae99SBarry Smith 
589*47c6ae99SBarry Smith @*/
590*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
591*47c6ae99SBarry Smith {
592*47c6ae99SBarry Smith   PetscFunctionBegin;
593*47c6ae99SBarry Smith   dm->ops->initialguess = f;
594*47c6ae99SBarry Smith   PetscFunctionReturn(0);
595*47c6ae99SBarry Smith }
596*47c6ae99SBarry Smith 
597*47c6ae99SBarry Smith #undef __FUNCT__
598*47c6ae99SBarry Smith #define __FUNCT__ "DMSetFunction"
599*47c6ae99SBarry Smith /*@
600*47c6ae99SBarry Smith     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
601*47c6ae99SBarry Smith 
602*47c6ae99SBarry Smith     Logically Collective on DM
603*47c6ae99SBarry Smith 
604*47c6ae99SBarry Smith     Input Parameter:
605*47c6ae99SBarry Smith +   dm - the DM object
606*47c6ae99SBarry Smith -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
607*47c6ae99SBarry Smith 
608*47c6ae99SBarry Smith     Level: intermediate
609*47c6ae99SBarry Smith 
610*47c6ae99SBarry Smith     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
611*47c6ae99SBarry Smith            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
612*47c6ae99SBarry Smith 
613*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
614*47c6ae99SBarry Smith          DMSetJacobian()
615*47c6ae99SBarry Smith 
616*47c6ae99SBarry Smith @*/
617*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
618*47c6ae99SBarry Smith {
619*47c6ae99SBarry Smith   PetscFunctionBegin;
620*47c6ae99SBarry Smith   dm->ops->function = f;
621*47c6ae99SBarry Smith   if (f) {
622*47c6ae99SBarry Smith     dm->ops->functionj = f;
623*47c6ae99SBarry Smith   }
624*47c6ae99SBarry Smith   PetscFunctionReturn(0);
625*47c6ae99SBarry Smith }
626*47c6ae99SBarry Smith 
627*47c6ae99SBarry Smith #undef __FUNCT__
628*47c6ae99SBarry Smith #define __FUNCT__ "DMSetJacobian"
629*47c6ae99SBarry Smith /*@
630*47c6ae99SBarry Smith     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
631*47c6ae99SBarry Smith 
632*47c6ae99SBarry Smith     Logically Collective on DM
633*47c6ae99SBarry Smith 
634*47c6ae99SBarry Smith     Input Parameter:
635*47c6ae99SBarry Smith +   dm - the DM object to destroy
636*47c6ae99SBarry Smith -   f - the function to compute the matrix entries
637*47c6ae99SBarry Smith 
638*47c6ae99SBarry Smith     Level: intermediate
639*47c6ae99SBarry Smith 
640*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
641*47c6ae99SBarry Smith          DMSetFunction()
642*47c6ae99SBarry Smith 
643*47c6ae99SBarry Smith @*/
644*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
645*47c6ae99SBarry Smith {
646*47c6ae99SBarry Smith   PetscFunctionBegin;
647*47c6ae99SBarry Smith   dm->ops->jacobian = f;
648*47c6ae99SBarry Smith   PetscFunctionReturn(0);
649*47c6ae99SBarry Smith }
650*47c6ae99SBarry Smith 
651*47c6ae99SBarry Smith #undef __FUNCT__
652*47c6ae99SBarry Smith #define __FUNCT__ "DMComputeInitialGuess"
653*47c6ae99SBarry Smith /*@
654*47c6ae99SBarry Smith     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
655*47c6ae99SBarry Smith 
656*47c6ae99SBarry Smith     Collective on DM
657*47c6ae99SBarry Smith 
658*47c6ae99SBarry Smith     Input Parameter:
659*47c6ae99SBarry Smith +   dm - the DM object to destroy
660*47c6ae99SBarry Smith -   x - the vector to hold the initial guess values
661*47c6ae99SBarry Smith 
662*47c6ae99SBarry Smith     Level: developer
663*47c6ae99SBarry Smith 
664*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetRhs(), DMSetMat()
665*47c6ae99SBarry Smith 
666*47c6ae99SBarry Smith @*/
667*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMComputeInitialGuess(DM dm,Vec x)
668*47c6ae99SBarry Smith {
669*47c6ae99SBarry Smith   PetscErrorCode ierr;
670*47c6ae99SBarry Smith   PetscFunctionBegin;
671*47c6ae99SBarry Smith   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
672*47c6ae99SBarry Smith   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
673*47c6ae99SBarry Smith   PetscFunctionReturn(0);
674*47c6ae99SBarry Smith }
675*47c6ae99SBarry Smith 
676*47c6ae99SBarry Smith #undef __FUNCT__
677*47c6ae99SBarry Smith #define __FUNCT__ "DMHasInitialGuess"
678*47c6ae99SBarry Smith /*@
679*47c6ae99SBarry Smith     DMHasInitialGuess - does the DM object have an initial guess function
680*47c6ae99SBarry Smith 
681*47c6ae99SBarry Smith     Not Collective
682*47c6ae99SBarry Smith 
683*47c6ae99SBarry Smith     Input Parameter:
684*47c6ae99SBarry Smith .   dm - the DM object to destroy
685*47c6ae99SBarry Smith 
686*47c6ae99SBarry Smith     Output Parameter:
687*47c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
688*47c6ae99SBarry Smith 
689*47c6ae99SBarry Smith     Level: developer
690*47c6ae99SBarry Smith 
691*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
692*47c6ae99SBarry Smith 
693*47c6ae99SBarry Smith @*/
694*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMHasInitialGuess(DM dm,PetscBool  *flg)
695*47c6ae99SBarry Smith {
696*47c6ae99SBarry Smith   PetscFunctionBegin;
697*47c6ae99SBarry Smith   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
698*47c6ae99SBarry Smith   PetscFunctionReturn(0);
699*47c6ae99SBarry Smith }
700*47c6ae99SBarry Smith 
701*47c6ae99SBarry Smith #undef __FUNCT__
702*47c6ae99SBarry Smith #define __FUNCT__ "DMHasFunction"
703*47c6ae99SBarry Smith /*@
704*47c6ae99SBarry Smith     DMHasFunction - does the DM object have a function
705*47c6ae99SBarry Smith 
706*47c6ae99SBarry Smith     Not Collective
707*47c6ae99SBarry Smith 
708*47c6ae99SBarry Smith     Input Parameter:
709*47c6ae99SBarry Smith .   dm - the DM object to destroy
710*47c6ae99SBarry Smith 
711*47c6ae99SBarry Smith     Output Parameter:
712*47c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
713*47c6ae99SBarry Smith 
714*47c6ae99SBarry Smith     Level: developer
715*47c6ae99SBarry Smith 
716*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
717*47c6ae99SBarry Smith 
718*47c6ae99SBarry Smith @*/
719*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMHasFunction(DM dm,PetscBool  *flg)
720*47c6ae99SBarry Smith {
721*47c6ae99SBarry Smith   PetscFunctionBegin;
722*47c6ae99SBarry Smith   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
723*47c6ae99SBarry Smith   PetscFunctionReturn(0);
724*47c6ae99SBarry Smith }
725*47c6ae99SBarry Smith 
726*47c6ae99SBarry Smith #undef __FUNCT__
727*47c6ae99SBarry Smith #define __FUNCT__ "DMHasJacobian"
728*47c6ae99SBarry Smith /*@
729*47c6ae99SBarry Smith     DMHasJacobian - does the DM object have a matrix function
730*47c6ae99SBarry Smith 
731*47c6ae99SBarry Smith     Not Collective
732*47c6ae99SBarry Smith 
733*47c6ae99SBarry Smith     Input Parameter:
734*47c6ae99SBarry Smith .   dm - the DM object to destroy
735*47c6ae99SBarry Smith 
736*47c6ae99SBarry Smith     Output Parameter:
737*47c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
738*47c6ae99SBarry Smith 
739*47c6ae99SBarry Smith     Level: developer
740*47c6ae99SBarry Smith 
741*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
742*47c6ae99SBarry Smith 
743*47c6ae99SBarry Smith @*/
744*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMHasJacobian(DM dm,PetscBool  *flg)
745*47c6ae99SBarry Smith {
746*47c6ae99SBarry Smith   PetscFunctionBegin;
747*47c6ae99SBarry Smith   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
748*47c6ae99SBarry Smith   PetscFunctionReturn(0);
749*47c6ae99SBarry Smith }
750*47c6ae99SBarry Smith 
751*47c6ae99SBarry Smith #undef __FUNCT__
752*47c6ae99SBarry Smith #define __FUNCT__ "DMComputeFunction"
753*47c6ae99SBarry Smith /*@
754*47c6ae99SBarry Smith     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
755*47c6ae99SBarry Smith 
756*47c6ae99SBarry Smith     Collective on DM
757*47c6ae99SBarry Smith 
758*47c6ae99SBarry Smith     Input Parameter:
759*47c6ae99SBarry Smith +   dm - the DM object to destroy
760*47c6ae99SBarry Smith .   x - the location where the function is evaluationed, may be ignored for linear problems
761*47c6ae99SBarry Smith -   b - the vector to hold the right hand side entries
762*47c6ae99SBarry Smith 
763*47c6ae99SBarry Smith     Level: developer
764*47c6ae99SBarry Smith 
765*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
766*47c6ae99SBarry Smith          DMSetJacobian()
767*47c6ae99SBarry Smith 
768*47c6ae99SBarry Smith @*/
769*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMComputeFunction(DM dm,Vec x,Vec b)
770*47c6ae99SBarry Smith {
771*47c6ae99SBarry Smith   PetscErrorCode ierr;
772*47c6ae99SBarry Smith   PetscFunctionBegin;
773*47c6ae99SBarry Smith   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
774*47c6ae99SBarry Smith   if (!x) x = dm->x;
775*47c6ae99SBarry Smith   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
776*47c6ae99SBarry Smith   PetscFunctionReturn(0);
777*47c6ae99SBarry Smith }
778*47c6ae99SBarry Smith 
779*47c6ae99SBarry Smith 
780*47c6ae99SBarry Smith #undef __FUNCT__
781*47c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobian"
782*47c6ae99SBarry Smith /*@
783*47c6ae99SBarry Smith     DMComputeJacobian - compute the matrix entries for the solver
784*47c6ae99SBarry Smith 
785*47c6ae99SBarry Smith     Collective on DM
786*47c6ae99SBarry Smith 
787*47c6ae99SBarry Smith     Input Parameter:
788*47c6ae99SBarry Smith +   dm - the DM object
789*47c6ae99SBarry Smith .   x - location to compute Jacobian at; may be ignored for linear problems
790*47c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
791*47c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
792*47c6ae99SBarry Smith 
793*47c6ae99SBarry Smith     Level: developer
794*47c6ae99SBarry Smith 
795*47c6ae99SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
796*47c6ae99SBarry Smith          DMSetFunction()
797*47c6ae99SBarry Smith 
798*47c6ae99SBarry Smith @*/
799*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
800*47c6ae99SBarry Smith {
801*47c6ae99SBarry Smith   PetscErrorCode ierr;
802*47c6ae99SBarry Smith 
803*47c6ae99SBarry Smith   PetscFunctionBegin;
804*47c6ae99SBarry Smith   if (!dm->ops->jacobian) {
805*47c6ae99SBarry Smith     ISColoring    coloring;
806*47c6ae99SBarry Smith     MatFDColoring fd;
807*47c6ae99SBarry Smith 
808*47c6ae99SBarry Smith     ierr = DMGetColoring(dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr);
809*47c6ae99SBarry Smith     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
810*47c6ae99SBarry Smith     ierr = ISColoringDestroy(coloring);CHKERRQ(ierr);
811*47c6ae99SBarry Smith     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
812*47c6ae99SBarry Smith     dm->fd = fd;
813*47c6ae99SBarry Smith     dm->ops->jacobian = DMComputeJacobianDefault;
814*47c6ae99SBarry Smith 
815*47c6ae99SBarry Smith     if (!dm->x) {
816*47c6ae99SBarry Smith       ierr = MatGetVecs(B,&dm->x,PETSC_NULL);CHKERRQ(ierr);
817*47c6ae99SBarry Smith     }
818*47c6ae99SBarry Smith   }
819*47c6ae99SBarry Smith   if (!x) x = dm->x;
820*47c6ae99SBarry Smith   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
821*47c6ae99SBarry Smith   PetscFunctionReturn(0);
822*47c6ae99SBarry Smith }
823