xref: /petsc/src/dm/interface/dm.c (revision efe48dd8cf8b935accbbb9f4bcb20bc83865fa4d)
1 #define PETSCDM_DLL
2 
3 #include "private/dmimpl.h"     /*I      "petscdm.h"     I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMCreate"
7 /*@
8   DMCreate - Creates an empty vector object. The type can then be set with DMetType().
9 
10    If you never  call DMSetType()  it will generate an
11    error when you try to use the vector.
12 
13   Collective on MPI_Comm
14 
15   Input Parameter:
16 . comm - The communicator for the DM object
17 
18   Output Parameter:
19 . dm - The DM object
20 
21   Level: beginner
22 
23 .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE
24 @*/
25 PetscErrorCode PETSCVEC_DLLEXPORT DMCreate(MPI_Comm comm, DM *vec)
26 {
27   DM             v;
28   PetscErrorCode ierr;
29 
30   PetscFunctionBegin;
31   PetscValidPointer(vec,2);
32   *vec = PETSC_NULL;
33 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
34   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
35 #endif
36 
37   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, -1, "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
38   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
39   *vec = v;
40   PetscFunctionReturn(0);
41 }
42 
43 
44 #undef __FUNCT__
45 #define __FUNCT__ "DMSetVecType"
46 /*@C
47        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
48 
49    Logically Collective on DMDA
50 
51    Input Parameter:
52 +  da - initial distributed array
53 .  ctype - the vector type, currently either VECSTANDARD or VECCUDA
54 
55    Options Database:
56 .   -da_vec_type ctype
57 
58    Level: intermediate
59 
60 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
61 @*/
62 PetscErrorCode PETSCDM_DLLEXPORT DMSetVecType(DM da,const VecType ctype)
63 {
64   PetscErrorCode ierr;
65 
66   PetscFunctionBegin;
67   PetscValidHeaderSpecific(da,DM_CLASSID,1);
68   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
69   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "DMSetOptionsPrefix"
75 /*@C
76    DMSetOptionsPrefix - Sets the prefix used for searching for all
77    DMDA options in the database.
78 
79    Logically Collective on DMDA
80 
81    Input Parameter:
82 +  da - the DMDA context
83 -  prefix - the prefix to prepend to all option names
84 
85    Notes:
86    A hyphen (-) must NOT be given at the beginning of the prefix name.
87    The first character of all runtime options is AUTOMATICALLY the hyphen.
88 
89    Level: advanced
90 
91 .keywords: DMDA, set, options, prefix, database
92 
93 .seealso: DMSetFromOptions()
94 @*/
95 PetscErrorCode PETSCDM_DLLEXPORT DMSetOptionsPrefix(DM dm,const char prefix[])
96 {
97   PetscErrorCode ierr;
98 
99   PetscFunctionBegin;
100   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
101   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
102   PetscFunctionReturn(0);
103 }
104 
105 #undef __FUNCT__
106 #define __FUNCT__ "DMDestroy"
107 /*@
108     DMDestroy - Destroys a vector packer or DMDA.
109 
110     Collective on DM
111 
112     Input Parameter:
113 .   dm - the DM object to destroy
114 
115     Level: developer
116 
117 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
118 
119 @*/
120 PetscErrorCode PETSCDM_DLLEXPORT DMDestroy(DM dm)
121 {
122   PetscErrorCode ierr;
123 
124   PetscFunctionBegin;
125   ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "DMSetUp"
131 /*@
132     DMSetUp - sets up the data structures inside a DM object
133 
134     Collective on DM
135 
136     Input Parameter:
137 .   dm - the DM object to setup
138 
139     Level: developer
140 
141 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
142 
143 @*/
144 PetscErrorCode PETSCDM_DLLEXPORT DMSetUp(DM dm)
145 {
146   PetscErrorCode ierr;
147 
148   PetscFunctionBegin;
149   if (dm->ops->setup) {
150     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
151   }
152   PetscFunctionReturn(0);
153 }
154 
155 #undef __FUNCT__
156 #define __FUNCT__ "DMSetFromOptions"
157 /*@
158     DMSetFromOptions - sets parameters in a DM from the options database
159 
160     Collective on DM
161 
162     Input Parameter:
163 .   dm - the DM object to set options for
164 
165     Level: developer
166 
167 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
168 
169 @*/
170 PetscErrorCode PETSCDM_DLLEXPORT DMSetFromOptions(DM dm)
171 {
172   PetscErrorCode ierr;
173 
174   PetscFunctionBegin;
175   if (dm->ops->setfromoptions) {
176     ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
177   }
178   PetscFunctionReturn(0);
179 }
180 
181 #undef __FUNCT__
182 #define __FUNCT__ "DMView"
183 /*@
184     DMView - Views a vector packer or DMDA.
185 
186     Collective on DM
187 
188     Input Parameter:
189 +   dm - the DM object to view
190 -   v - the viewer
191 
192     Level: developer
193 
194 .seealso DMDestroy(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
195 
196 @*/
197 PetscErrorCode PETSCDM_DLLEXPORT DMView(DM dm,PetscViewer v)
198 {
199   PetscErrorCode ierr;
200 
201   PetscFunctionBegin;
202  if (!v) {
203     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
204   }
205   if (dm->ops->view) {
206     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
207   }
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "DMCreateGlobalVector"
213 /*@
214     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
215 
216     Collective on DM
217 
218     Input Parameter:
219 .   dm - the DM object
220 
221     Output Parameter:
222 .   vec - the global vector
223 
224     Level: developer
225 
226 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
227 
228 @*/
229 PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector(DM dm,Vec *vec)
230 {
231   PetscErrorCode ierr;
232 
233   PetscFunctionBegin;
234   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "DMCreateLocalVector"
240 /*@
241     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
242 
243     Not Collective
244 
245     Input Parameter:
246 .   dm - the DM object
247 
248     Output Parameter:
249 .   vec - the local vector
250 
251     Level: developer
252 
253 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
254 
255 @*/
256 PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector(DM dm,Vec *vec)
257 {
258   PetscErrorCode ierr;
259 
260   PetscFunctionBegin;
261   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 
265 #undef __FUNCT__
266 #define __FUNCT__ "DMGetInterpolation"
267 /*@
268     DMGetInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
269 
270     Collective on DM
271 
272     Input Parameter:
273 +   dm1 - the DM object
274 -   dm2 - the second, finer DM object
275 
276     Output Parameter:
277 +  mat - the interpolation
278 -  vec - the scaling (optional)
279 
280     Level: developer
281 
282 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix()
283 
284 @*/
285 PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
286 {
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   ierr = (*dm1->ops->getinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "DMGetInjection"
296 /*@
297     DMGetInjection - Gets injection matrix between two DMDA or DMComposite objects
298 
299     Collective on DM
300 
301     Input Parameter:
302 +   dm1 - the DM object
303 -   dm2 - the second, finer DM object
304 
305     Output Parameter:
306 .   ctx - the injection
307 
308     Level: developer
309 
310 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMGetInterpolation()
311 
312 @*/
313 PetscErrorCode PETSCDM_DLLEXPORT DMGetInjection(DM dm1,DM dm2,VecScatter *ctx)
314 {
315   PetscErrorCode ierr;
316 
317   PetscFunctionBegin;
318   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
319   PetscFunctionReturn(0);
320 }
321 
322 #undef __FUNCT__
323 #define __FUNCT__ "DMGetColoring"
324 /*@
325     DMGetColoring - Gets coloring for a DMDA or DMComposite
326 
327     Collective on DM
328 
329     Input Parameter:
330 +   dm - the DM object
331 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
332 -   matype - either MATAIJ or MATBAIJ
333 
334     Output Parameter:
335 .   coloring - the coloring
336 
337     Level: developer
338 
339 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
340 
341 @*/
342 PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
343 {
344   PetscErrorCode ierr;
345 
346   PetscFunctionBegin;
347   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
348   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
349   PetscFunctionReturn(0);
350 }
351 
352 #undef __FUNCT__
353 #define __FUNCT__ "DMGetMatrix"
354 /*@C
355     DMGetMatrix - Gets empty Jacobian for a DMDA or DMComposite
356 
357     Collective on DM
358 
359     Input Parameter:
360 +   dm - the DM object
361 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
362             any type which inherits from one of these (such as MATAIJ)
363 
364     Output Parameter:
365 .   mat - the empty Jacobian
366 
367     Level: developer
368 
369     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
370        do not need to do it yourself.
371 
372        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
373        the nonzero pattern call DMDASetMatPreallocateOnly()
374 
375        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
376        internally by PETSc.
377 
378        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
379        the indices for the global numbering for DMDAs which is complicated.
380 
381 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
382 
383 @*/
384 PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix(DM dm, const MatType mtype,Mat *mat)
385 {
386   PetscErrorCode ierr;
387 
388   PetscFunctionBegin;
389   ierr = (*dm->ops->getmatrix)(dm,mtype,mat);CHKERRQ(ierr);
390   PetscFunctionReturn(0);
391 }
392 
393 #undef __FUNCT__
394 #define __FUNCT__ "DMRefine"
395 /*@
396     DMRefine - Refines a DM object
397 
398     Collective on DM
399 
400     Input Parameter:
401 +   dm - the DM object
402 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
403 
404     Output Parameter:
405 .   dmf - the refined DM
406 
407     Level: developer
408 
409 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
410 
411 @*/
412 PetscErrorCode PETSCDM_DLLEXPORT DMRefine(DM dm,MPI_Comm comm,DM *dmf)
413 {
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin;
417   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
418   PetscFunctionReturn(0);
419 }
420 
421 #undef __FUNCT__
422 #define __FUNCT__ "DMGlobalToLocalBegin"
423 /*@
424     DMGlobalToLocalBegin - Begins updating local vectors from global vector
425 
426     Neighbor-wise Collective on DM
427 
428     Input Parameters:
429 +   dm - the DM object
430 .   g - the global vector
431 .   mode - INSERT_VALUES or ADD_VALUES
432 -   l - the local vector
433 
434 
435     Level: beginner
436 
437 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
438 
439 @*/
440 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
441 {
442   PetscErrorCode ierr;
443 
444   PetscFunctionBegin;
445   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr);
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "DMGlobalToLocalEnd"
451 /*@
452     DMGlobalToLocalEnd - Ends updating local vectors from global vector
453 
454     Neighbor-wise Collective on DM
455 
456     Input Parameters:
457 +   dm - the DM object
458 .   g - the global vector
459 .   mode - INSERT_VALUES or ADD_VALUES
460 -   l - the local vector
461 
462 
463     Level: beginner
464 
465 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
466 
467 @*/
468 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
469 {
470   PetscErrorCode ierr;
471 
472   PetscFunctionBegin;
473   ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr);
474   PetscFunctionReturn(0);
475 }
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "DMLocalToGlobalBegin"
479 /*@
480     DMLocalToGlobalBegin - updates global vectors from local vectors
481 
482     Neighbor-wise Collective on DM
483 
484     Input Parameters:
485 +   dm - the DM object
486 .   l - the local vector
487 .   mode - if INSERT_VALUES then no parallel communication is used, if ADD_VALUES then all ghost points from the same base point accumulate into that
488            base point.
489 - - the global vector
490 
491     Notes: In the ADD_VALUES case you normally would zero the receiving vector before beginning this operation. If you would like to simply add the non-ghosted values in the local
492            array into the global array you need to either (1) zero the ghosted locations and use ADD_VALUES or (2) use INSERT_VALUES into a work global array and then add the work
493            global array to the final global array with VecAXPY().
494 
495     Level: beginner
496 
497 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
498 
499 @*/
500 PetscErrorCode PETSCDM_DLLEXPORT DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
501 {
502   PetscErrorCode ierr;
503 
504   PetscFunctionBegin;
505   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode,g);CHKERRQ(ierr);
506   PetscFunctionReturn(0);
507 }
508 
509 #undef __FUNCT__
510 #define __FUNCT__ "DMLocalToGlobalEnd"
511 /*@
512     DMLocalToGlobalEnd - updates global vectors from local vectors
513 
514     Neighbor-wise Collective on DM
515 
516     Input Parameters:
517 +   dm - the DM object
518 .   l - the local vector
519 .   mode - INSERT_VALUES or ADD_VALUES
520 -   g - the global vector
521 
522 
523     Level: beginner
524 
525 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
526 
527 @*/
528 PetscErrorCode PETSCDM_DLLEXPORT DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
529 {
530   PetscErrorCode ierr;
531 
532   PetscFunctionBegin;
533   ierr = (*dm->ops->localtoglobalend)(dm,l,mode,g);CHKERRQ(ierr);
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "DMComputeJacobianDefault"
539 /*@
540     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
541 
542     Collective on DM
543 
544     Input Parameter:
545 +   dm - the DM object
546 .   x - location to compute Jacobian at; may be ignored for linear problems
547 .   A - matrix that defines the operator for the linear solve
548 -   B - the matrix used to construct the preconditioner
549 
550     Level: developer
551 
552 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
553          DMSetFunction()
554 
555 @*/
556 PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
557 {
558   PetscErrorCode ierr;
559   PetscFunctionBegin;
560   *stflag = SAME_NONZERO_PATTERN;
561   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
562   if (A != B) {
563     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
564     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
565   }
566   PetscFunctionReturn(0);
567 }
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "DMCoarsen"
571 /*@
572     DMCoarsen - Coarsens a DM object
573 
574     Collective on DM
575 
576     Input Parameter:
577 +   dm - the DM object
578 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
579 
580     Output Parameter:
581 .   dmc - the coarsened DM
582 
583     Level: developer
584 
585 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
586 
587 @*/
588 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
589 {
590   PetscErrorCode ierr;
591 
592   PetscFunctionBegin;
593   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
594   (*dmc)->ops->initialguess = dm->ops->initialguess;
595   (*dmc)->ops->function     = dm->ops->function;
596   (*dmc)->ops->functionj    = dm->ops->functionj;
597   if (dm->ops->jacobian != DMComputeJacobianDefault) {
598     (*dmc)->ops->jacobian     = dm->ops->jacobian;
599   }
600   PetscFunctionReturn(0);
601 }
602 
603 #undef __FUNCT__
604 #define __FUNCT__ "DMRefineHierarchy"
605 /*@C
606     DMRefineHierarchy - Refines a DM object, all levels at once
607 
608     Collective on DM
609 
610     Input Parameter:
611 +   dm - the DM object
612 -   nlevels - the number of levels of refinement
613 
614     Output Parameter:
615 .   dmf - the refined DM hierarchy
616 
617     Level: developer
618 
619 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
620 
621 @*/
622 PetscErrorCode PETSCDM_DLLEXPORT DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
623 {
624   PetscErrorCode ierr;
625 
626   PetscFunctionBegin;
627   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
628   if (nlevels == 0) PetscFunctionReturn(0);
629   if (dm->ops->refinehierarchy) {
630     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
631   } else if (dm->ops->refine) {
632     PetscInt i;
633 
634     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
635     for (i=1; i<nlevels; i++) {
636       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
637     }
638   } else {
639     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
640   }
641   PetscFunctionReturn(0);
642 }
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "DMCoarsenHierarchy"
646 /*@C
647     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
648 
649     Collective on DM
650 
651     Input Parameter:
652 +   dm - the DM object
653 -   nlevels - the number of levels of coarsening
654 
655     Output Parameter:
656 .   dmc - the coarsened DM hierarchy
657 
658     Level: developer
659 
660 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
661 
662 @*/
663 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
664 {
665   PetscErrorCode ierr;
666 
667   PetscFunctionBegin;
668   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
669   if (nlevels == 0) PetscFunctionReturn(0);
670   PetscValidPointer(dmc,3);
671   if (dm->ops->coarsenhierarchy) {
672     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
673   } else if (dm->ops->coarsen) {
674     PetscInt i;
675 
676     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
677     for (i=1; i<nlevels; i++) {
678       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
679     }
680   } else {
681     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
682   }
683   PetscFunctionReturn(0);
684 }
685 
686 #undef __FUNCT__
687 #define __FUNCT__ "DMGetAggregates"
688 /*@
689    DMGetAggregates - Gets the aggregates that map between
690    grids associated with two DMs.
691 
692    Collective on DM
693 
694    Input Parameters:
695 +  dmc - the coarse grid DM
696 -  dmf - the fine grid DM
697 
698    Output Parameters:
699 .  rest - the restriction matrix (transpose of the projection matrix)
700 
701    Level: intermediate
702 
703 .keywords: interpolation, restriction, multigrid
704 
705 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation()
706 @*/
707 PetscErrorCode PETSCDM_DLLEXPORT DMGetAggregates(DM dmc, DM dmf, Mat *rest)
708 {
709   PetscErrorCode ierr;
710 
711   PetscFunctionBegin;
712   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
713   PetscFunctionReturn(0);
714 }
715 
716 #undef __FUNCT__
717 #define __FUNCT__ "DMSetContext"
718 /*@
719     DMSetContext - Set a user context into a DM object
720 
721     Not Collective
722 
723     Input Parameters:
724 +   dm - the DM object
725 -   ctx - the user context
726 
727     Level: intermediate
728 
729 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
730 
731 @*/
732 PetscErrorCode PETSCDM_DLLEXPORT DMSetContext(DM dm,void *ctx)
733 {
734   PetscFunctionBegin;
735   dm->ctx = ctx;
736   PetscFunctionReturn(0);
737 }
738 
739 #undef __FUNCT__
740 #define __FUNCT__ "DMGetContext"
741 /*@
742     DMGetContext - Gets a user context from a DM object
743 
744     Not Collective
745 
746     Input Parameter:
747 .   dm - the DM object
748 
749     Output Parameter:
750 .   ctx - the user context
751 
752     Level: intermediate
753 
754 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
755 
756 @*/
757 PetscErrorCode PETSCDM_DLLEXPORT DMGetContext(DM dm,void **ctx)
758 {
759   PetscFunctionBegin;
760   *ctx = dm->ctx;
761   PetscFunctionReturn(0);
762 }
763 
764 #undef __FUNCT__
765 #define __FUNCT__ "DMSetInitialGuess"
766 /*@
767     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
768 
769     Logically Collective on DM
770 
771     Input Parameter:
772 +   dm - the DM object to destroy
773 -   f - the function to compute the initial guess
774 
775     Level: intermediate
776 
777 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
778 
779 @*/
780 PetscErrorCode PETSCDM_DLLEXPORT DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
781 {
782   PetscFunctionBegin;
783   dm->ops->initialguess = f;
784   PetscFunctionReturn(0);
785 }
786 
787 #undef __FUNCT__
788 #define __FUNCT__ "DMSetFunction"
789 /*@
790     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
791 
792     Logically Collective on DM
793 
794     Input Parameter:
795 +   dm - the DM object
796 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
797 
798     Level: intermediate
799 
800     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
801            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
802 
803 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
804          DMSetJacobian()
805 
806 @*/
807 PetscErrorCode PETSCDM_DLLEXPORT DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
808 {
809   PetscFunctionBegin;
810   dm->ops->function = f;
811   if (f) {
812     dm->ops->functionj = f;
813   }
814   PetscFunctionReturn(0);
815 }
816 
817 #undef __FUNCT__
818 #define __FUNCT__ "DMSetJacobian"
819 /*@
820     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
821 
822     Logically Collective on DM
823 
824     Input Parameter:
825 +   dm - the DM object to destroy
826 -   f - the function to compute the matrix entries
827 
828     Level: intermediate
829 
830 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
831          DMSetFunction()
832 
833 @*/
834 PetscErrorCode PETSCDM_DLLEXPORT DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
835 {
836   PetscFunctionBegin;
837   dm->ops->jacobian = f;
838   PetscFunctionReturn(0);
839 }
840 
841 #undef __FUNCT__
842 #define __FUNCT__ "DMComputeInitialGuess"
843 /*@
844     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
845 
846     Collective on DM
847 
848     Input Parameter:
849 +   dm - the DM object to destroy
850 -   x - the vector to hold the initial guess values
851 
852     Level: developer
853 
854 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetRhs(), DMSetMat()
855 
856 @*/
857 PetscErrorCode PETSCDM_DLLEXPORT DMComputeInitialGuess(DM dm,Vec x)
858 {
859   PetscErrorCode ierr;
860   PetscFunctionBegin;
861   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
862   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
863   PetscFunctionReturn(0);
864 }
865 
866 #undef __FUNCT__
867 #define __FUNCT__ "DMHasInitialGuess"
868 /*@
869     DMHasInitialGuess - does the DM object have an initial guess function
870 
871     Not Collective
872 
873     Input Parameter:
874 .   dm - the DM object to destroy
875 
876     Output Parameter:
877 .   flg - PETSC_TRUE if function exists
878 
879     Level: developer
880 
881 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
882 
883 @*/
884 PetscErrorCode PETSCDM_DLLEXPORT DMHasInitialGuess(DM dm,PetscBool  *flg)
885 {
886   PetscFunctionBegin;
887   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
888   PetscFunctionReturn(0);
889 }
890 
891 #undef __FUNCT__
892 #define __FUNCT__ "DMHasFunction"
893 /*@
894     DMHasFunction - does the DM object have a function
895 
896     Not Collective
897 
898     Input Parameter:
899 .   dm - the DM object to destroy
900 
901     Output Parameter:
902 .   flg - PETSC_TRUE if function exists
903 
904     Level: developer
905 
906 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
907 
908 @*/
909 PetscErrorCode PETSCDM_DLLEXPORT DMHasFunction(DM dm,PetscBool  *flg)
910 {
911   PetscFunctionBegin;
912   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
913   PetscFunctionReturn(0);
914 }
915 
916 #undef __FUNCT__
917 #define __FUNCT__ "DMHasJacobian"
918 /*@
919     DMHasJacobian - does the DM object have a matrix function
920 
921     Not Collective
922 
923     Input Parameter:
924 .   dm - the DM object to destroy
925 
926     Output Parameter:
927 .   flg - PETSC_TRUE if function exists
928 
929     Level: developer
930 
931 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
932 
933 @*/
934 PetscErrorCode PETSCDM_DLLEXPORT DMHasJacobian(DM dm,PetscBool  *flg)
935 {
936   PetscFunctionBegin;
937   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
938   PetscFunctionReturn(0);
939 }
940 
941 #undef __FUNCT__
942 #define __FUNCT__ "DMComputeFunction"
943 /*@
944     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
945 
946     Collective on DM
947 
948     Input Parameter:
949 +   dm - the DM object to destroy
950 .   x - the location where the function is evaluationed, may be ignored for linear problems
951 -   b - the vector to hold the right hand side entries
952 
953     Level: developer
954 
955 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
956          DMSetJacobian()
957 
958 @*/
959 PetscErrorCode PETSCDM_DLLEXPORT DMComputeFunction(DM dm,Vec x,Vec b)
960 {
961   PetscErrorCode ierr;
962   PetscFunctionBegin;
963   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
964   if (!x) x = dm->x;
965   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
966   PetscFunctionReturn(0);
967 }
968 
969 
970 #undef __FUNCT__
971 #define __FUNCT__ "DMComputeJacobian"
972 /*@
973     DMComputeJacobian - compute the matrix entries for the solver
974 
975     Collective on DM
976 
977     Input Parameter:
978 +   dm - the DM object
979 .   x - location to compute Jacobian at; may be ignored for linear problems
980 .   A - matrix that defines the operator for the linear solve
981 -   B - the matrix used to construct the preconditioner
982 
983     Level: developer
984 
985 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
986          DMSetFunction()
987 
988 @*/
989 PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
990 {
991   PetscErrorCode ierr;
992 
993   PetscFunctionBegin;
994   if (!dm->ops->jacobian) {
995     ISColoring    coloring;
996     MatFDColoring fd;
997 
998     ierr = DMGetColoring(dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr);
999     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1000     ierr = ISColoringDestroy(coloring);CHKERRQ(ierr);
1001     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1002     dm->fd = fd;
1003     dm->ops->jacobian = DMComputeJacobianDefault;
1004 
1005     if (!dm->x) {
1006       ierr = MatGetVecs(B,&dm->x,PETSC_NULL);CHKERRQ(ierr);
1007     }
1008   }
1009   if (!x) x = dm->x;
1010   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1011   PetscFunctionReturn(0);
1012 }
1013 
1014 PetscFList DMList                       = PETSC_NULL;
1015 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1016 
1017 #undef __FUNCT__
1018 #define __FUNCT__ "DMSetType"
1019 /*@C
1020   DMSetType - Builds a DM, for a particular DM implementation.
1021 
1022   Collective on DM
1023 
1024   Input Parameters:
1025 + dm     - The DM object
1026 - method - The name of the DM type
1027 
1028   Options Database Key:
1029 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1030 
1031   Notes:
1032   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1033 
1034   Level: intermediate
1035 
1036 .keywords: DM, set, type
1037 .seealso: DMGetType(), DMCreate()
1038 @*/
1039 PetscErrorCode PETSCDM_DLLEXPORT DMSetType(DM dm, const DMType method)
1040 {
1041   PetscErrorCode (*r)(DM);
1042   PetscBool      match;
1043   PetscErrorCode ierr;
1044 
1045   PetscFunctionBegin;
1046   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1047   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1048   if (match) PetscFunctionReturn(0);
1049 
1050   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1051   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,(void (**)(void)) &r);CHKERRQ(ierr);
1052   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1053 
1054   if (dm->ops->destroy) {
1055     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1056   }
1057   ierr = (*r)(dm);CHKERRQ(ierr);
1058   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1059   PetscFunctionReturn(0);
1060 }
1061 
1062 #undef __FUNCT__
1063 #define __FUNCT__ "DMGetType"
1064 /*@C
1065   DMGetType - Gets the DM type name (as a string) from the DM.
1066 
1067   Not Collective
1068 
1069   Input Parameter:
1070 . dm  - The DM
1071 
1072   Output Parameter:
1073 . type - The DM type name
1074 
1075   Level: intermediate
1076 
1077 .keywords: DM, get, type, name
1078 .seealso: DMSetType(), DMCreate()
1079 @*/
1080 PetscErrorCode PETSCDM_DLLEXPORT DMGetType(DM dm, const DMType *type)
1081 {
1082   PetscErrorCode ierr;
1083 
1084   PetscFunctionBegin;
1085   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1086   PetscValidCharPointer(type,2);
1087   if (!DMRegisterAllCalled) {
1088     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1089   }
1090   *type = ((PetscObject)dm)->type_name;
1091   PetscFunctionReturn(0);
1092 }
1093 
1094 
1095 /*--------------------------------------------------------------------------------------------------------------------*/
1096 
1097 #undef __FUNCT__
1098 #define __FUNCT__ "DMRegister"
1099 /*@C
1100   DMRegister - See DMRegisterDynamic()
1101 
1102   Level: advanced
1103 @*/
1104 PetscErrorCode PETSCDM_DLLEXPORT DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1105 {
1106   char fullname[PETSC_MAX_PATH_LEN];
1107   PetscErrorCode ierr;
1108 
1109   PetscFunctionBegin;
1110   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1111   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1112   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1113   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1114   PetscFunctionReturn(0);
1115 }
1116 
1117 
1118 /*--------------------------------------------------------------------------------------------------------------------*/
1119 #undef __FUNCT__
1120 #define __FUNCT__ "DMRegisterDestroy"
1121 /*@C
1122    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1123 
1124    Not Collective
1125 
1126    Level: advanced
1127 
1128 .keywords: DM, register, destroy
1129 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1130 @*/
1131 PetscErrorCode PETSCDM_DLLEXPORT DMRegisterDestroy(void)
1132 {
1133   PetscErrorCode ierr;
1134 
1135   PetscFunctionBegin;
1136   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1137   DMRegisterAllCalled = PETSC_FALSE;
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1142 #include "mex.h"
1143 
1144 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
1145 
1146 #undef __FUNCT__
1147 #define __FUNCT__ "DMComputeFunction_Matlab"
1148 /*
1149    DMComputeFunction_Matlab - Calls the function that has been set with
1150                          DMSetFunctionMatlab().
1151 
1152    For linear problems x is null
1153 
1154 .seealso: DMSetFunction(), DMGetFunction()
1155 */
1156 PetscErrorCode PETSCDM_DLLEXPORT DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
1157 {
1158   PetscErrorCode    ierr;
1159   DMMatlabContext   *sctx;
1160   int               nlhs = 1,nrhs = 4;
1161   mxArray	    *plhs[1],*prhs[4];
1162   long long int     lx = 0,ly = 0,ls = 0;
1163 
1164   PetscFunctionBegin;
1165   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1166   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1167   PetscCheckSameComm(dm,1,y,3);
1168 
1169   /* call Matlab function in ctx with arguments x and y */
1170   ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr);
1171   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1172   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1173   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
1174   prhs[0] =  mxCreateDoubleScalar((double)ls);
1175   prhs[1] =  mxCreateDoubleScalar((double)lx);
1176   prhs[2] =  mxCreateDoubleScalar((double)ly);
1177   prhs[3] =  mxCreateString(sctx->funcname);
1178   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
1179   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
1180   mxDestroyArray(prhs[0]);
1181   mxDestroyArray(prhs[1]);
1182   mxDestroyArray(prhs[2]);
1183   mxDestroyArray(prhs[3]);
1184   mxDestroyArray(plhs[0]);
1185   PetscFunctionReturn(0);
1186 }
1187 
1188 
1189 #undef __FUNCT__
1190 #define __FUNCT__ "DMSetFunctionMatlab"
1191 /*
1192    DMSetFunctionMatlab - Sets the function evaluation routine
1193 
1194 */
1195 PetscErrorCode PETSCDM_DLLEXPORT DMSetFunctionMatlab(DM dm,const char *func)
1196 {
1197   PetscErrorCode    ierr;
1198   DMMatlabContext   *sctx;
1199 
1200   PetscFunctionBegin;
1201   /* currently sctx is memory bleed */
1202   ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr);
1203   if (!sctx) {
1204     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1205   }
1206   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
1207   ierr = DMSetContext(dm,sctx);CHKERRQ(ierr);
1208   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 #undef __FUNCT__
1213 #define __FUNCT__ "DMComputeJacobian_Matlab"
1214 /*
1215    DMComputeJacobian_Matlab - Calls the function that has been set with
1216                          DMSetJacobianMatlab().
1217 
1218    For linear problems x is null
1219 
1220 .seealso: DMSetFunction(), DMGetFunction()
1221 */
1222 PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
1223 {
1224   PetscErrorCode    ierr;
1225   DMMatlabContext   *sctx;
1226   int               nlhs = 2,nrhs = 5;
1227   mxArray	    *plhs[2],*prhs[5];
1228   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
1229 
1230   PetscFunctionBegin;
1231   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1232   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
1233 
1234   /* call Matlab function in ctx with arguments x, A, and B */
1235   ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr);
1236   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1237   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1238   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
1239   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
1240   prhs[0] =  mxCreateDoubleScalar((double)ls);
1241   prhs[1] =  mxCreateDoubleScalar((double)lx);
1242   prhs[2] =  mxCreateDoubleScalar((double)lA);
1243   prhs[3] =  mxCreateDoubleScalar((double)lB);
1244   prhs[4] =  mxCreateString(sctx->jacname);
1245   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
1246   *str    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
1247   ierr    =  mxGetScalar(plhs[1]);CHKERRQ(ierr);
1248   mxDestroyArray(prhs[0]);
1249   mxDestroyArray(prhs[1]);
1250   mxDestroyArray(prhs[2]);
1251   mxDestroyArray(prhs[3]);
1252   mxDestroyArray(prhs[4]);
1253   mxDestroyArray(plhs[0]);
1254   mxDestroyArray(plhs[1]);
1255   PetscFunctionReturn(0);
1256 }
1257 
1258 
1259 #undef __FUNCT__
1260 #define __FUNCT__ "DMSetJacobianMatlab"
1261 /*
1262    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
1263 
1264 */
1265 PetscErrorCode PETSCDM_DLLEXPORT DMSetJacobianMatlab(DM dm,const char *func)
1266 {
1267   PetscErrorCode    ierr;
1268   DMMatlabContext   *sctx;
1269 
1270   PetscFunctionBegin;
1271   /* currently sctx is memory bleed */
1272   ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr);
1273   if (!sctx) {
1274     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1275   }
1276   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
1277   ierr = DMSetContext(dm,sctx);CHKERRQ(ierr);
1278   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
1279   PetscFunctionReturn(0);
1280 }
1281 #endif
1282