xref: /petsc/src/dm/interface/dm.c (revision c4721b0ef1399f9b4324b73a33db78f8d329fdf8)
147c6ae99SBarry Smith 
2c6db04a5SJed Brown #include <private/dmimpl.h>     /*I      "petscdm.h"     I*/
347c6ae99SBarry Smith 
4732e2eb9SMatthew G Knepley PetscClassId  DM_CLASSID;
567a56275SMatthew G Knepley PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal;
667a56275SMatthew G Knepley 
747c6ae99SBarry Smith #undef __FUNCT__
8a4121054SBarry Smith #define __FUNCT__ "DMCreate"
9a4121054SBarry Smith /*@
10de043629SMatthew G Knepley   DMCreate - Creates an empty DM object. The type can then be set with DMSetType().
11a4121054SBarry Smith 
12a4121054SBarry Smith    If you never  call DMSetType()  it will generate an
13a4121054SBarry Smith    error when you try to use the vector.
14a4121054SBarry Smith 
15a4121054SBarry Smith   Collective on MPI_Comm
16a4121054SBarry Smith 
17a4121054SBarry Smith   Input Parameter:
18a4121054SBarry Smith . comm - The communicator for the DM object
19a4121054SBarry Smith 
20a4121054SBarry Smith   Output Parameter:
21a4121054SBarry Smith . dm - The DM object
22a4121054SBarry Smith 
23a4121054SBarry Smith   Level: beginner
24a4121054SBarry Smith 
25a4121054SBarry Smith .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE
26a4121054SBarry Smith @*/
277087cfbeSBarry Smith PetscErrorCode  DMCreate(MPI_Comm comm,DM *dm)
28a4121054SBarry Smith {
29a4121054SBarry Smith   DM             v;
30a4121054SBarry Smith   PetscErrorCode ierr;
31a4121054SBarry Smith 
32a4121054SBarry Smith   PetscFunctionBegin;
331411c6eeSJed Brown   PetscValidPointer(dm,2);
341411c6eeSJed Brown   *dm = PETSC_NULL;
35a4121054SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
36b84caa0eSBarry Smith   ierr = VecInitializePackage(PETSC_NULL);CHKERRQ(ierr);
37b84caa0eSBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
38a4121054SBarry Smith   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
39a4121054SBarry Smith #endif
40a4121054SBarry Smith 
413194b578SJed Brown   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, -1, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
42a4121054SBarry Smith   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
431411c6eeSJed Brown 
44a89ea682SMatthew G Knepley   v->workSize     = 0;
45a89ea682SMatthew G Knepley   v->workArray    = PETSC_NULL;
461411c6eeSJed Brown   v->ltogmap      = PETSC_NULL;
471411c6eeSJed Brown   v->ltogmapb     = PETSC_NULL;
481411c6eeSJed Brown   v->bs           = 1;
49171400e9SBarry Smith   v->coloringtype = IS_COLORING_GLOBAL;
501411c6eeSJed Brown 
511411c6eeSJed Brown   *dm = v;
52a4121054SBarry Smith   PetscFunctionReturn(0);
53a4121054SBarry Smith }
54a4121054SBarry Smith 
55a4121054SBarry Smith 
56a4121054SBarry Smith #undef __FUNCT__
579a42bb27SBarry Smith #define __FUNCT__ "DMSetVecType"
589a42bb27SBarry Smith /*@C
59564755cdSBarry Smith        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
609a42bb27SBarry Smith 
61aa219208SBarry Smith    Logically Collective on DMDA
629a42bb27SBarry Smith 
639a42bb27SBarry Smith    Input Parameter:
649a42bb27SBarry Smith +  da - initial distributed array
658154be41SBarry Smith .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
669a42bb27SBarry Smith 
679a42bb27SBarry Smith    Options Database:
68dd85299cSBarry Smith .   -dm_vec_type ctype
699a42bb27SBarry Smith 
709a42bb27SBarry Smith    Level: intermediate
719a42bb27SBarry Smith 
72aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
739a42bb27SBarry Smith @*/
747087cfbeSBarry Smith PetscErrorCode  DMSetVecType(DM da,const VecType ctype)
759a42bb27SBarry Smith {
769a42bb27SBarry Smith   PetscErrorCode ierr;
779a42bb27SBarry Smith 
789a42bb27SBarry Smith   PetscFunctionBegin;
799a42bb27SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
809a42bb27SBarry Smith   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
819a42bb27SBarry Smith   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
829a42bb27SBarry Smith   PetscFunctionReturn(0);
839a42bb27SBarry Smith }
849a42bb27SBarry Smith 
859a42bb27SBarry Smith #undef __FUNCT__
869a42bb27SBarry Smith #define __FUNCT__ "DMSetOptionsPrefix"
879a42bb27SBarry Smith /*@C
889a42bb27SBarry Smith    DMSetOptionsPrefix - Sets the prefix used for searching for all
89aa219208SBarry Smith    DMDA options in the database.
909a42bb27SBarry Smith 
91aa219208SBarry Smith    Logically Collective on DMDA
929a42bb27SBarry Smith 
939a42bb27SBarry Smith    Input Parameter:
94aa219208SBarry Smith +  da - the DMDA context
959a42bb27SBarry Smith -  prefix - the prefix to prepend to all option names
969a42bb27SBarry Smith 
979a42bb27SBarry Smith    Notes:
989a42bb27SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
999a42bb27SBarry Smith    The first character of all runtime options is AUTOMATICALLY the hyphen.
1009a42bb27SBarry Smith 
1019a42bb27SBarry Smith    Level: advanced
1029a42bb27SBarry Smith 
103aa219208SBarry Smith .keywords: DMDA, set, options, prefix, database
1049a42bb27SBarry Smith 
1059a42bb27SBarry Smith .seealso: DMSetFromOptions()
1069a42bb27SBarry Smith @*/
1077087cfbeSBarry Smith PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
1089a42bb27SBarry Smith {
1099a42bb27SBarry Smith   PetscErrorCode ierr;
1109a42bb27SBarry Smith 
1119a42bb27SBarry Smith   PetscFunctionBegin;
1129a42bb27SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1139a42bb27SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
1149a42bb27SBarry Smith   PetscFunctionReturn(0);
1159a42bb27SBarry Smith }
1169a42bb27SBarry Smith 
1179a42bb27SBarry Smith #undef __FUNCT__
11847c6ae99SBarry Smith #define __FUNCT__ "DMDestroy"
11947c6ae99SBarry Smith /*@
120aa219208SBarry Smith     DMDestroy - Destroys a vector packer or DMDA.
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith     Collective on DM
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith     Input Parameter:
12547c6ae99SBarry Smith .   dm - the DM object to destroy
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith     Level: developer
12847c6ae99SBarry Smith 
129e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
13047c6ae99SBarry Smith 
13147c6ae99SBarry Smith @*/
132fcfd50ebSBarry Smith PetscErrorCode  DMDestroy(DM *dm)
13347c6ae99SBarry Smith {
134732e2eb9SMatthew G Knepley   PetscInt       i, cnt = 0;
13547c6ae99SBarry Smith   PetscErrorCode ierr;
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith   PetscFunctionBegin;
1386bf464f9SBarry Smith   if (!*dm) PetscFunctionReturn(0);
1396bf464f9SBarry Smith   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
14087e657c6SBarry Smith 
14187e657c6SBarry Smith   /* count all the circular references of DM and its contained Vecs */
142732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1436bf464f9SBarry Smith     if ((*dm)->localin[i])  {cnt++;}
1446bf464f9SBarry Smith     if ((*dm)->globalin[i]) {cnt++;}
145732e2eb9SMatthew G Knepley   }
14671cd77b2SBarry Smith   if ((*dm)->x) {
14771cd77b2SBarry Smith     PetscObject obj;
14871cd77b2SBarry Smith     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
14971cd77b2SBarry Smith     if (obj == (PetscObject)*dm) cnt++;
15071cd77b2SBarry Smith   }
151732e2eb9SMatthew G Knepley 
1526bf464f9SBarry Smith   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
153732e2eb9SMatthew G Knepley   /*
154732e2eb9SMatthew G Knepley      Need this test because the dm references the vectors that
155732e2eb9SMatthew G Knepley      reference the dm, so destroying the dm calls destroy on the
156732e2eb9SMatthew G Knepley      vectors that cause another destroy on the dm
157732e2eb9SMatthew G Knepley   */
1586bf464f9SBarry Smith   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
1596bf464f9SBarry Smith   ((PetscObject) (*dm))->refct = 0;
160732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1616bf464f9SBarry Smith     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
1626bf464f9SBarry Smith     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
163732e2eb9SMatthew G Knepley   }
1641a266240SBarry Smith 
1651a266240SBarry Smith   if ((*dm)->ctx && (*dm)->ctxdestroy) {
1661a266240SBarry Smith     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
1671a266240SBarry Smith   }
16887e657c6SBarry Smith   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
16971cd77b2SBarry Smith   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
1704dcab191SBarry Smith   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
1716bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
1726bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
1736bf464f9SBarry Smith   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
174073dac72SJed Brown   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
175a89ea682SMatthew G Knepley   ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr);
176732e2eb9SMatthew G Knepley   /* if memory was published with AMS then destroy it */
1776bf464f9SBarry Smith   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
178732e2eb9SMatthew G Knepley 
1796bf464f9SBarry Smith   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
1806bf464f9SBarry Smith   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
181732e2eb9SMatthew G Knepley   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
18247c6ae99SBarry Smith   PetscFunctionReturn(0);
18347c6ae99SBarry Smith }
18447c6ae99SBarry Smith 
18547c6ae99SBarry Smith #undef __FUNCT__
186d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp"
187d7bf68aeSBarry Smith /*@
188d7bf68aeSBarry Smith     DMSetUp - sets up the data structures inside a DM object
189d7bf68aeSBarry Smith 
190d7bf68aeSBarry Smith     Collective on DM
191d7bf68aeSBarry Smith 
192d7bf68aeSBarry Smith     Input Parameter:
193d7bf68aeSBarry Smith .   dm - the DM object to setup
194d7bf68aeSBarry Smith 
195d7bf68aeSBarry Smith     Level: developer
196d7bf68aeSBarry Smith 
197e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
198d7bf68aeSBarry Smith 
199d7bf68aeSBarry Smith @*/
2007087cfbeSBarry Smith PetscErrorCode  DMSetUp(DM dm)
201d7bf68aeSBarry Smith {
202d7bf68aeSBarry Smith   PetscErrorCode ierr;
203d7bf68aeSBarry Smith 
204d7bf68aeSBarry Smith   PetscFunctionBegin;
205171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2068387afaaSJed Brown   if (dm->setupcalled) PetscFunctionReturn(0);
207d7bf68aeSBarry Smith   if (dm->ops->setup) {
208d7bf68aeSBarry Smith     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
209d7bf68aeSBarry Smith   }
2108387afaaSJed Brown   dm->setupcalled = PETSC_TRUE;
211d7bf68aeSBarry Smith   PetscFunctionReturn(0);
212d7bf68aeSBarry Smith }
213d7bf68aeSBarry Smith 
214d7bf68aeSBarry Smith #undef __FUNCT__
215d7bf68aeSBarry Smith #define __FUNCT__ "DMSetFromOptions"
216d7bf68aeSBarry Smith /*@
217d7bf68aeSBarry Smith     DMSetFromOptions - sets parameters in a DM from the options database
218d7bf68aeSBarry Smith 
219d7bf68aeSBarry Smith     Collective on DM
220d7bf68aeSBarry Smith 
221d7bf68aeSBarry Smith     Input Parameter:
222d7bf68aeSBarry Smith .   dm - the DM object to set options for
223d7bf68aeSBarry Smith 
224732e2eb9SMatthew G Knepley     Options Database:
225dd85299cSBarry Smith +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
226dd85299cSBarry Smith .   -dm_vec_type <type>  type of vector to create inside DM
227171400e9SBarry Smith .   -dm_mat_type <type>  type of matrix to create inside DM
228171400e9SBarry Smith -   -dm_coloring_type <global or ghosted>
229732e2eb9SMatthew G Knepley 
230d7bf68aeSBarry Smith     Level: developer
231d7bf68aeSBarry Smith 
232e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
233d7bf68aeSBarry Smith 
234d7bf68aeSBarry Smith @*/
2357087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions(DM dm)
236d7bf68aeSBarry Smith {
237*c4721b0eSMatthew G Knepley   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg;
238d7bf68aeSBarry Smith   PetscErrorCode ierr;
239f9ba7244SBarry Smith   char           typeName[256] = MATAIJ;
240d7bf68aeSBarry Smith 
241d7bf68aeSBarry Smith   PetscFunctionBegin;
242171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2433194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
24482fcb398SMatthew G Knepley     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
245*c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr);
246*c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr);
247073dac72SJed Brown     ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,PETSC_NULL);CHKERRQ(ierr);
248f9ba7244SBarry Smith     ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
249f9ba7244SBarry Smith     if (flg) {
250f9ba7244SBarry Smith       ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
251f9ba7244SBarry Smith     }
252f9ba7244SBarry Smith     ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr);
253073dac72SJed Brown     if (flg) {
254073dac72SJed Brown       ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
255f9ba7244SBarry Smith       ierr = PetscStrallocpy(typeName,&dm->mattype);CHKERRQ(ierr);
256073dac72SJed Brown     }
2571b89239cSHong Zhang     ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr);
258f9ba7244SBarry Smith     if (dm->ops->setfromoptions) {
259f9ba7244SBarry Smith       ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
260f9ba7244SBarry Smith     }
261f9ba7244SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
262f9ba7244SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
26382fcb398SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
26482fcb398SMatthew G Knepley   if (flg1) {
26582fcb398SMatthew G Knepley     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
26682fcb398SMatthew G Knepley   }
267*c4721b0eSMatthew G Knepley   if (flg2) {
268*c4721b0eSMatthew G Knepley     PetscViewer viewer;
269*c4721b0eSMatthew G Knepley 
270*c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
271*c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
272*c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
273*c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
274*c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
275*c4721b0eSMatthew G Knepley   }
276*c4721b0eSMatthew G Knepley   if (flg3) {
277*c4721b0eSMatthew G Knepley     PetscViewer viewer;
278*c4721b0eSMatthew G Knepley 
279*c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
280*c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
281*c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
282*c4721b0eSMatthew G Knepley     ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr);
283*c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
284*c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
285*c4721b0eSMatthew G Knepley   }
286d7bf68aeSBarry Smith   PetscFunctionReturn(0);
287d7bf68aeSBarry Smith }
288d7bf68aeSBarry Smith 
289d7bf68aeSBarry Smith #undef __FUNCT__
29047c6ae99SBarry Smith #define __FUNCT__ "DMView"
291fc9bc008SSatish Balay /*@C
292aa219208SBarry Smith     DMView - Views a vector packer or DMDA.
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith     Collective on DM
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith     Input Parameter:
29747c6ae99SBarry Smith +   dm - the DM object to view
29847c6ae99SBarry Smith -   v - the viewer
29947c6ae99SBarry Smith 
30047c6ae99SBarry Smith     Level: developer
30147c6ae99SBarry Smith 
302e727c939SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
30347c6ae99SBarry Smith 
30447c6ae99SBarry Smith @*/
3057087cfbeSBarry Smith PetscErrorCode  DMView(DM dm,PetscViewer v)
30647c6ae99SBarry Smith {
30747c6ae99SBarry Smith   PetscErrorCode ierr;
30847c6ae99SBarry Smith 
30947c6ae99SBarry Smith   PetscFunctionBegin;
310171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3113014e516SBarry Smith  if (!v) {
3123014e516SBarry Smith     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
3133014e516SBarry Smith   }
3140c010503SBarry Smith   if (dm->ops->view) {
3150c010503SBarry Smith     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
31647c6ae99SBarry Smith   }
31747c6ae99SBarry Smith   PetscFunctionReturn(0);
31847c6ae99SBarry Smith }
31947c6ae99SBarry Smith 
32047c6ae99SBarry Smith #undef __FUNCT__
32147c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector"
32247c6ae99SBarry Smith /*@
323aa219208SBarry Smith     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith     Collective on DM
32647c6ae99SBarry Smith 
32747c6ae99SBarry Smith     Input Parameter:
32847c6ae99SBarry Smith .   dm - the DM object
32947c6ae99SBarry Smith 
33047c6ae99SBarry Smith     Output Parameter:
33147c6ae99SBarry Smith .   vec - the global vector
33247c6ae99SBarry Smith 
333073dac72SJed Brown     Level: beginner
33447c6ae99SBarry Smith 
335e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith @*/
3387087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
33947c6ae99SBarry Smith {
34047c6ae99SBarry Smith   PetscErrorCode ierr;
34147c6ae99SBarry Smith 
34247c6ae99SBarry Smith   PetscFunctionBegin;
343171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
34447c6ae99SBarry Smith   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
34547c6ae99SBarry Smith   PetscFunctionReturn(0);
34647c6ae99SBarry Smith }
34747c6ae99SBarry Smith 
34847c6ae99SBarry Smith #undef __FUNCT__
34947c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector"
35047c6ae99SBarry Smith /*@
351aa219208SBarry Smith     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
35247c6ae99SBarry Smith 
35347c6ae99SBarry Smith     Not Collective
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith     Input Parameter:
35647c6ae99SBarry Smith .   dm - the DM object
35747c6ae99SBarry Smith 
35847c6ae99SBarry Smith     Output Parameter:
35947c6ae99SBarry Smith .   vec - the local vector
36047c6ae99SBarry Smith 
361073dac72SJed Brown     Level: beginner
36247c6ae99SBarry Smith 
363e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
36447c6ae99SBarry Smith 
36547c6ae99SBarry Smith @*/
3667087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
36747c6ae99SBarry Smith {
36847c6ae99SBarry Smith   PetscErrorCode ierr;
36947c6ae99SBarry Smith 
37047c6ae99SBarry Smith   PetscFunctionBegin;
371171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37247c6ae99SBarry Smith   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
37347c6ae99SBarry Smith   PetscFunctionReturn(0);
37447c6ae99SBarry Smith }
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith #undef __FUNCT__
3771411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping"
3781411c6eeSJed Brown /*@
3791411c6eeSJed Brown    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
3801411c6eeSJed Brown 
3811411c6eeSJed Brown    Collective on DM
3821411c6eeSJed Brown 
3831411c6eeSJed Brown    Input Parameter:
3841411c6eeSJed Brown .  dm - the DM that provides the mapping
3851411c6eeSJed Brown 
3861411c6eeSJed Brown    Output Parameter:
3871411c6eeSJed Brown .  ltog - the mapping
3881411c6eeSJed Brown 
3891411c6eeSJed Brown    Level: intermediate
3901411c6eeSJed Brown 
3911411c6eeSJed Brown    Notes:
3921411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMapping() or
3931411c6eeSJed Brown    MatSetLocalToGlobalMapping().
3941411c6eeSJed Brown 
3951411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
3961411c6eeSJed Brown @*/
3977087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
3981411c6eeSJed Brown {
3991411c6eeSJed Brown   PetscErrorCode ierr;
4001411c6eeSJed Brown 
4011411c6eeSJed Brown   PetscFunctionBegin;
4021411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4031411c6eeSJed Brown   PetscValidPointer(ltog,2);
4041411c6eeSJed Brown   if (!dm->ltogmap) {
4051411c6eeSJed Brown     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
4061411c6eeSJed Brown     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
4071411c6eeSJed Brown   }
4081411c6eeSJed Brown   *ltog = dm->ltogmap;
4091411c6eeSJed Brown   PetscFunctionReturn(0);
4101411c6eeSJed Brown }
4111411c6eeSJed Brown 
4121411c6eeSJed Brown #undef __FUNCT__
4131411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
4141411c6eeSJed Brown /*@
4151411c6eeSJed Brown    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
4161411c6eeSJed Brown 
4171411c6eeSJed Brown    Collective on DM
4181411c6eeSJed Brown 
4191411c6eeSJed Brown    Input Parameter:
4201411c6eeSJed Brown .  da - the distributed array that provides the mapping
4211411c6eeSJed Brown 
4221411c6eeSJed Brown    Output Parameter:
4231411c6eeSJed Brown .  ltog - the block mapping
4241411c6eeSJed Brown 
4251411c6eeSJed Brown    Level: intermediate
4261411c6eeSJed Brown 
4271411c6eeSJed Brown    Notes:
4281411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
4291411c6eeSJed Brown    MatSetLocalToGlobalMappingBlock().
4301411c6eeSJed Brown 
4311411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
4321411c6eeSJed Brown @*/
4337087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
4341411c6eeSJed Brown {
4351411c6eeSJed Brown   PetscErrorCode ierr;
4361411c6eeSJed Brown 
4371411c6eeSJed Brown   PetscFunctionBegin;
4381411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4391411c6eeSJed Brown   PetscValidPointer(ltog,2);
4401411c6eeSJed Brown   if (!dm->ltogmapb) {
4411411c6eeSJed Brown     PetscInt bs;
4421411c6eeSJed Brown     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
4431411c6eeSJed Brown     if (bs > 1) {
4441411c6eeSJed Brown       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
4451411c6eeSJed Brown       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
4461411c6eeSJed Brown     } else {
4471411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
4481411c6eeSJed Brown       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
4491411c6eeSJed Brown     }
4501411c6eeSJed Brown   }
4511411c6eeSJed Brown   *ltog = dm->ltogmapb;
4521411c6eeSJed Brown   PetscFunctionReturn(0);
4531411c6eeSJed Brown }
4541411c6eeSJed Brown 
4551411c6eeSJed Brown #undef __FUNCT__
4561411c6eeSJed Brown #define __FUNCT__ "DMGetBlockSize"
4571411c6eeSJed Brown /*@
4581411c6eeSJed Brown    DMGetBlockSize - Gets the inherent block size associated with a DM
4591411c6eeSJed Brown 
4601411c6eeSJed Brown    Not Collective
4611411c6eeSJed Brown 
4621411c6eeSJed Brown    Input Parameter:
4631411c6eeSJed Brown .  dm - the DM with block structure
4641411c6eeSJed Brown 
4651411c6eeSJed Brown    Output Parameter:
4661411c6eeSJed Brown .  bs - the block size, 1 implies no exploitable block structure
4671411c6eeSJed Brown 
4681411c6eeSJed Brown    Level: intermediate
4691411c6eeSJed Brown 
4701411c6eeSJed Brown .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
4711411c6eeSJed Brown @*/
4727087cfbeSBarry Smith PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
4731411c6eeSJed Brown {
4741411c6eeSJed Brown   PetscFunctionBegin;
4751411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4761411c6eeSJed Brown   PetscValidPointer(bs,2);
4771411c6eeSJed Brown   if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet");
4781411c6eeSJed Brown   *bs = dm->bs;
4791411c6eeSJed Brown   PetscFunctionReturn(0);
4801411c6eeSJed Brown }
4811411c6eeSJed Brown 
4821411c6eeSJed Brown #undef __FUNCT__
483e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation"
48447c6ae99SBarry Smith /*@
485e727c939SJed Brown     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
48647c6ae99SBarry Smith 
48747c6ae99SBarry Smith     Collective on DM
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith     Input Parameter:
49047c6ae99SBarry Smith +   dm1 - the DM object
49147c6ae99SBarry Smith -   dm2 - the second, finer DM object
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith     Output Parameter:
49447c6ae99SBarry Smith +  mat - the interpolation
49547c6ae99SBarry Smith -  vec - the scaling (optional)
49647c6ae99SBarry Smith 
49747c6ae99SBarry Smith     Level: developer
49847c6ae99SBarry Smith 
49985afcc9aSBarry Smith     Notes:  For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by
50085afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
501d52bd9f3SBarry Smith 
502d52bd9f3SBarry Smith         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
503d52bd9f3SBarry Smith         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
50485afcc9aSBarry Smith 
50585afcc9aSBarry Smith 
506e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith @*/
509e727c939SJed Brown PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
51047c6ae99SBarry Smith {
51147c6ae99SBarry Smith   PetscErrorCode ierr;
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith   PetscFunctionBegin;
514171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
515171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
51625296bd5SBarry Smith   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
51747c6ae99SBarry Smith   PetscFunctionReturn(0);
51847c6ae99SBarry Smith }
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith #undef __FUNCT__
521e727c939SJed Brown #define __FUNCT__ "DMCreateInjection"
52247c6ae99SBarry Smith /*@
523e727c939SJed Brown     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
52447c6ae99SBarry Smith 
52547c6ae99SBarry Smith     Collective on DM
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith     Input Parameter:
52847c6ae99SBarry Smith +   dm1 - the DM object
52947c6ae99SBarry Smith -   dm2 - the second, finer DM object
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith     Output Parameter:
53247c6ae99SBarry Smith .   ctx - the injection
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith     Level: developer
53547c6ae99SBarry Smith 
53685afcc9aSBarry Smith    Notes:  For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by
53785afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
53885afcc9aSBarry Smith 
539e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith @*/
542e727c939SJed Brown PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
54347c6ae99SBarry Smith {
54447c6ae99SBarry Smith   PetscErrorCode ierr;
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith   PetscFunctionBegin;
547171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
548171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
54947c6ae99SBarry Smith   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
55047c6ae99SBarry Smith   PetscFunctionReturn(0);
55147c6ae99SBarry Smith }
55247c6ae99SBarry Smith 
55347c6ae99SBarry Smith #undef __FUNCT__
554e727c939SJed Brown #define __FUNCT__ "DMCreateColoring"
555d1e2c406SBarry Smith /*@C
556e727c939SJed Brown     DMCreateColoring - Gets coloring for a DMDA or DMComposite
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith     Collective on DM
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith     Input Parameter:
56147c6ae99SBarry Smith +   dm - the DM object
56247c6ae99SBarry Smith .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
56347c6ae99SBarry Smith -   matype - either MATAIJ or MATBAIJ
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith     Output Parameter:
56647c6ae99SBarry Smith .   coloring - the coloring
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith     Level: developer
56947c6ae99SBarry Smith 
570e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
57147c6ae99SBarry Smith 
572aab9d709SJed Brown @*/
573e727c939SJed Brown PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
57447c6ae99SBarry Smith {
57547c6ae99SBarry Smith   PetscErrorCode ierr;
57647c6ae99SBarry Smith 
57747c6ae99SBarry Smith   PetscFunctionBegin;
578171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
57947c6ae99SBarry Smith   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
58047c6ae99SBarry Smith   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
58147c6ae99SBarry Smith   PetscFunctionReturn(0);
58247c6ae99SBarry Smith }
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith #undef __FUNCT__
585950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix"
58647c6ae99SBarry Smith /*@C
587950540a4SJed Brown     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
58847c6ae99SBarry Smith 
58947c6ae99SBarry Smith     Collective on DM
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith     Input Parameter:
59247c6ae99SBarry Smith +   dm - the DM object
59347c6ae99SBarry Smith -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
59494013140SBarry Smith             any type which inherits from one of these (such as MATAIJ)
59547c6ae99SBarry Smith 
59647c6ae99SBarry Smith     Output Parameter:
59747c6ae99SBarry Smith .   mat - the empty Jacobian
59847c6ae99SBarry Smith 
599073dac72SJed Brown     Level: beginner
60047c6ae99SBarry Smith 
60194013140SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
60294013140SBarry Smith        do not need to do it yourself.
60394013140SBarry Smith 
60494013140SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
605aa219208SBarry Smith        the nonzero pattern call DMDASetMatPreallocateOnly()
60694013140SBarry Smith 
60794013140SBarry Smith        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
60894013140SBarry Smith        internally by PETSc.
60994013140SBarry Smith 
61094013140SBarry Smith        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
611aa219208SBarry Smith        the indices for the global numbering for DMDAs which is complicated.
61294013140SBarry Smith 
613e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
61447c6ae99SBarry Smith 
615aab9d709SJed Brown @*/
616950540a4SJed Brown PetscErrorCode  DMCreateMatrix(DM dm,const MatType mtype,Mat *mat)
61747c6ae99SBarry Smith {
61847c6ae99SBarry Smith   PetscErrorCode ierr;
61947c6ae99SBarry Smith 
62047c6ae99SBarry Smith   PetscFunctionBegin;
621171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
622235683edSBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
623235683edSBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
624235683edSBarry Smith #endif
625c7b7c8a4SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
626c7b7c8a4SJed Brown   PetscValidPointer(mat,3);
627073dac72SJed Brown   if (dm->mattype) {
62825296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
629073dac72SJed Brown   } else {
63025296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
631c7b7c8a4SJed Brown   }
63247c6ae99SBarry Smith   PetscFunctionReturn(0);
63347c6ae99SBarry Smith }
63447c6ae99SBarry Smith 
63547c6ae99SBarry Smith #undef __FUNCT__
636732e2eb9SMatthew G Knepley #define __FUNCT__ "DMSetMatrixPreallocateOnly"
637732e2eb9SMatthew G Knepley /*@
638950540a4SJed Brown   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
639732e2eb9SMatthew G Knepley     preallocated but the nonzero structure and zero values will not be set.
640732e2eb9SMatthew G Knepley 
641732e2eb9SMatthew G Knepley   Logically Collective on DMDA
642732e2eb9SMatthew G Knepley 
643732e2eb9SMatthew G Knepley   Input Parameter:
644732e2eb9SMatthew G Knepley + dm - the DM
645732e2eb9SMatthew G Knepley - only - PETSC_TRUE if only want preallocation
646732e2eb9SMatthew G Knepley 
647732e2eb9SMatthew G Knepley   Level: developer
648950540a4SJed Brown .seealso DMCreateMatrix()
649732e2eb9SMatthew G Knepley @*/
650732e2eb9SMatthew G Knepley PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
651732e2eb9SMatthew G Knepley {
652732e2eb9SMatthew G Knepley   PetscFunctionBegin;
653732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
654732e2eb9SMatthew G Knepley   dm->prealloc_only = only;
655732e2eb9SMatthew G Knepley   PetscFunctionReturn(0);
656732e2eb9SMatthew G Knepley }
657732e2eb9SMatthew G Knepley 
658732e2eb9SMatthew G Knepley #undef __FUNCT__
659a89ea682SMatthew G Knepley #define __FUNCT__ "DMGetWorkArray"
660a89ea682SMatthew G Knepley /*@C
661a89ea682SMatthew G Knepley   DMGetWorkArray - Gets a work array guaranteed to be at least the input size
662a89ea682SMatthew G Knepley 
663a89ea682SMatthew G Knepley   Not Collective
664a89ea682SMatthew G Knepley 
665a89ea682SMatthew G Knepley   Input Parameters:
666a89ea682SMatthew G Knepley + dm - the DM object
667a89ea682SMatthew G Knepley - size - The minium size
668a89ea682SMatthew G Knepley 
669a89ea682SMatthew G Knepley   Output Parameter:
670a89ea682SMatthew G Knepley . array - the work array
671a89ea682SMatthew G Knepley 
672a89ea682SMatthew G Knepley   Level: developer
673a89ea682SMatthew G Knepley 
674a89ea682SMatthew G Knepley .seealso DMDestroy(), DMCreate()
675a89ea682SMatthew G Knepley @*/
676a89ea682SMatthew G Knepley PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array)
677a89ea682SMatthew G Knepley {
678a89ea682SMatthew G Knepley   PetscErrorCode ierr;
679a89ea682SMatthew G Knepley 
680a89ea682SMatthew G Knepley   PetscFunctionBegin;
681a89ea682SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
682a89ea682SMatthew G Knepley   PetscValidPointer(array,3);
683a89ea682SMatthew G Knepley   if (size > dm->workSize) {
684a89ea682SMatthew G Knepley     dm->workSize = size;
685a89ea682SMatthew G Knepley     ierr = PetscFree(dm->workArray);CHKERRQ(ierr);
686a89ea682SMatthew G Knepley     ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr);
687a89ea682SMatthew G Knepley   }
688a89ea682SMatthew G Knepley   *array = dm->workArray;
689a89ea682SMatthew G Knepley   PetscFunctionReturn(0);
690a89ea682SMatthew G Knepley }
691a89ea682SMatthew G Knepley 
692a89ea682SMatthew G Knepley #undef __FUNCT__
69347c6ae99SBarry Smith #define __FUNCT__ "DMRefine"
69447c6ae99SBarry Smith /*@
69547c6ae99SBarry Smith     DMRefine - Refines a DM object
69647c6ae99SBarry Smith 
69747c6ae99SBarry Smith     Collective on DM
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith     Input Parameter:
70047c6ae99SBarry Smith +   dm - the DM object
70147c6ae99SBarry Smith -   comm - the communicator to contain the new DM object (or PETSC_NULL)
70247c6ae99SBarry Smith 
70347c6ae99SBarry Smith     Output Parameter:
70447c6ae99SBarry Smith .   dmf - the refined DM
70547c6ae99SBarry Smith 
70647c6ae99SBarry Smith     Level: developer
70747c6ae99SBarry Smith 
708e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith @*/
7117087cfbeSBarry Smith PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
71247c6ae99SBarry Smith {
71347c6ae99SBarry Smith   PetscErrorCode ierr;
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith   PetscFunctionBegin;
716732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
71747c6ae99SBarry Smith   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
718644e2e5bSBarry Smith   (*dmf)->ops->initialguess = dm->ops->initialguess;
719644e2e5bSBarry Smith   (*dmf)->ops->function     = dm->ops->function;
720644e2e5bSBarry Smith   (*dmf)->ops->functionj    = dm->ops->functionj;
721644e2e5bSBarry Smith   if (dm->ops->jacobian != DMComputeJacobianDefault) {
722644e2e5bSBarry Smith     (*dmf)->ops->jacobian     = dm->ops->jacobian;
723644e2e5bSBarry Smith   }
7248cd211a4SJed Brown   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
725644e2e5bSBarry Smith   (*dmf)->ctx     = dm->ctx;
726656b349aSBarry Smith   (*dmf)->levelup = dm->levelup + 1;
72747c6ae99SBarry Smith   PetscFunctionReturn(0);
72847c6ae99SBarry Smith }
72947c6ae99SBarry Smith 
73047c6ae99SBarry Smith #undef __FUNCT__
731eb3f98d2SBarry Smith #define __FUNCT__ "DMGetRefineLevel"
732eb3f98d2SBarry Smith /*@
733eb3f98d2SBarry Smith     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
734eb3f98d2SBarry Smith 
735eb3f98d2SBarry Smith     Not Collective
736eb3f98d2SBarry Smith 
737eb3f98d2SBarry Smith     Input Parameter:
738eb3f98d2SBarry Smith .   dm - the DM object
739eb3f98d2SBarry Smith 
740eb3f98d2SBarry Smith     Output Parameter:
741eb3f98d2SBarry Smith .   level - number of refinements
742eb3f98d2SBarry Smith 
743eb3f98d2SBarry Smith     Level: developer
744eb3f98d2SBarry Smith 
745e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
746eb3f98d2SBarry Smith 
747eb3f98d2SBarry Smith @*/
748eb3f98d2SBarry Smith PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
749eb3f98d2SBarry Smith {
750eb3f98d2SBarry Smith   PetscFunctionBegin;
751eb3f98d2SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
752eb3f98d2SBarry Smith   *level = dm->levelup;
753eb3f98d2SBarry Smith   PetscFunctionReturn(0);
754eb3f98d2SBarry Smith }
755eb3f98d2SBarry Smith 
756eb3f98d2SBarry Smith #undef __FUNCT__
75747c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin"
75847c6ae99SBarry Smith /*@
75947c6ae99SBarry Smith     DMGlobalToLocalBegin - Begins updating local vectors from global vector
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith     Neighbor-wise Collective on DM
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith     Input Parameters:
76447c6ae99SBarry Smith +   dm - the DM object
76547c6ae99SBarry Smith .   g - the global vector
76647c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
76747c6ae99SBarry Smith -   l - the local vector
76847c6ae99SBarry Smith 
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith     Level: beginner
77147c6ae99SBarry Smith 
772e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith @*/
7757087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
77647c6ae99SBarry Smith {
77747c6ae99SBarry Smith   PetscErrorCode ierr;
77847c6ae99SBarry Smith 
77947c6ae99SBarry Smith   PetscFunctionBegin;
780171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
781843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
78247c6ae99SBarry Smith   PetscFunctionReturn(0);
78347c6ae99SBarry Smith }
78447c6ae99SBarry Smith 
78547c6ae99SBarry Smith #undef __FUNCT__
78647c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd"
78747c6ae99SBarry Smith /*@
78847c6ae99SBarry Smith     DMGlobalToLocalEnd - Ends updating local vectors from global vector
78947c6ae99SBarry Smith 
79047c6ae99SBarry Smith     Neighbor-wise Collective on DM
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith     Input Parameters:
79347c6ae99SBarry Smith +   dm - the DM object
79447c6ae99SBarry Smith .   g - the global vector
79547c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
79647c6ae99SBarry Smith -   l - the local vector
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith     Level: beginner
80047c6ae99SBarry Smith 
801e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
80247c6ae99SBarry Smith 
80347c6ae99SBarry Smith @*/
8047087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
80547c6ae99SBarry Smith {
80647c6ae99SBarry Smith   PetscErrorCode ierr;
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith   PetscFunctionBegin;
809171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
810843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
81147c6ae99SBarry Smith   PetscFunctionReturn(0);
81247c6ae99SBarry Smith }
81347c6ae99SBarry Smith 
81447c6ae99SBarry Smith #undef __FUNCT__
8159a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalBegin"
81647c6ae99SBarry Smith /*@
8179a42bb27SBarry Smith     DMLocalToGlobalBegin - updates global vectors from local vectors
8189a42bb27SBarry Smith 
8199a42bb27SBarry Smith     Neighbor-wise Collective on DM
8209a42bb27SBarry Smith 
8219a42bb27SBarry Smith     Input Parameters:
8229a42bb27SBarry Smith +   dm - the DM object
823f6813fd5SJed Brown .   l - the local vector
8249a42bb27SBarry Smith .   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
8259a42bb27SBarry Smith            base point.
826f6813fd5SJed Brown - - the global vector
8279a42bb27SBarry Smith 
8289a42bb27SBarry Smith     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
8299a42bb27SBarry Smith            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
8309a42bb27SBarry Smith            global array to the final global array with VecAXPY().
8319a42bb27SBarry Smith 
8329a42bb27SBarry Smith     Level: beginner
8339a42bb27SBarry Smith 
834e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
8359a42bb27SBarry Smith 
8369a42bb27SBarry Smith @*/
8377087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
8389a42bb27SBarry Smith {
8399a42bb27SBarry Smith   PetscErrorCode ierr;
8409a42bb27SBarry Smith 
8419a42bb27SBarry Smith   PetscFunctionBegin;
842171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
843843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
8449a42bb27SBarry Smith   PetscFunctionReturn(0);
8459a42bb27SBarry Smith }
8469a42bb27SBarry Smith 
8479a42bb27SBarry Smith #undef __FUNCT__
8489a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalEnd"
8499a42bb27SBarry Smith /*@
8509a42bb27SBarry Smith     DMLocalToGlobalEnd - updates global vectors from local vectors
85147c6ae99SBarry Smith 
85247c6ae99SBarry Smith     Neighbor-wise Collective on DM
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith     Input Parameters:
85547c6ae99SBarry Smith +   dm - the DM object
856f6813fd5SJed Brown .   l - the local vector
85747c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
858f6813fd5SJed Brown -   g - the global vector
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith     Level: beginner
86247c6ae99SBarry Smith 
863e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
86447c6ae99SBarry Smith 
86547c6ae99SBarry Smith @*/
8667087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
86747c6ae99SBarry Smith {
86847c6ae99SBarry Smith   PetscErrorCode ierr;
86947c6ae99SBarry Smith 
87047c6ae99SBarry Smith   PetscFunctionBegin;
871171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
872843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
87347c6ae99SBarry Smith   PetscFunctionReturn(0);
87447c6ae99SBarry Smith }
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith #undef __FUNCT__
87747c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobianDefault"
87847c6ae99SBarry Smith /*@
87947c6ae99SBarry Smith     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
88047c6ae99SBarry Smith 
88147c6ae99SBarry Smith     Collective on DM
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith     Input Parameter:
88447c6ae99SBarry Smith +   dm - the DM object
88547c6ae99SBarry Smith .   x - location to compute Jacobian at; may be ignored for linear problems
88647c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
88747c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith     Level: developer
89047c6ae99SBarry Smith 
891e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
89247c6ae99SBarry Smith          DMSetFunction()
89347c6ae99SBarry Smith 
89447c6ae99SBarry Smith @*/
8957087cfbeSBarry Smith PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
89647c6ae99SBarry Smith {
89747c6ae99SBarry Smith   PetscErrorCode ierr;
898171400e9SBarry Smith 
89947c6ae99SBarry Smith   PetscFunctionBegin;
900171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
90147c6ae99SBarry Smith   *stflag = SAME_NONZERO_PATTERN;
90247c6ae99SBarry Smith   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
90347c6ae99SBarry Smith   if (A != B) {
90447c6ae99SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90547c6ae99SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90647c6ae99SBarry Smith   }
90747c6ae99SBarry Smith   PetscFunctionReturn(0);
90847c6ae99SBarry Smith }
90947c6ae99SBarry Smith 
91047c6ae99SBarry Smith #undef __FUNCT__
91147c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen"
91247c6ae99SBarry Smith /*@
91347c6ae99SBarry Smith     DMCoarsen - Coarsens a DM object
91447c6ae99SBarry Smith 
91547c6ae99SBarry Smith     Collective on DM
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith     Input Parameter:
91847c6ae99SBarry Smith +   dm - the DM object
91947c6ae99SBarry Smith -   comm - the communicator to contain the new DM object (or PETSC_NULL)
92047c6ae99SBarry Smith 
92147c6ae99SBarry Smith     Output Parameter:
92247c6ae99SBarry Smith .   dmc - the coarsened DM
92347c6ae99SBarry Smith 
92447c6ae99SBarry Smith     Level: developer
92547c6ae99SBarry Smith 
926e727c939SJed Brown .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
92747c6ae99SBarry Smith 
92847c6ae99SBarry Smith @*/
9297087cfbeSBarry Smith PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
93047c6ae99SBarry Smith {
93147c6ae99SBarry Smith   PetscErrorCode ierr;
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith   PetscFunctionBegin;
934171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
93547c6ae99SBarry Smith   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
93647c6ae99SBarry Smith   (*dmc)->ops->initialguess = dm->ops->initialguess;
93747c6ae99SBarry Smith   (*dmc)->ops->function     = dm->ops->function;
93847c6ae99SBarry Smith   (*dmc)->ops->functionj    = dm->ops->functionj;
93947c6ae99SBarry Smith   if (dm->ops->jacobian != DMComputeJacobianDefault) {
94047c6ae99SBarry Smith     (*dmc)->ops->jacobian     = dm->ops->jacobian;
94147c6ae99SBarry Smith   }
9428cd211a4SJed Brown   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
943644e2e5bSBarry Smith   (*dmc)->ctx       = dm->ctx;
944656b349aSBarry Smith   (*dmc)->leveldown = dm->leveldown + 1;
94547c6ae99SBarry Smith   PetscFunctionReturn(0);
94647c6ae99SBarry Smith }
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith #undef __FUNCT__
94947c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy"
95047c6ae99SBarry Smith /*@C
95147c6ae99SBarry Smith     DMRefineHierarchy - Refines a DM object, all levels at once
95247c6ae99SBarry Smith 
95347c6ae99SBarry Smith     Collective on DM
95447c6ae99SBarry Smith 
95547c6ae99SBarry Smith     Input Parameter:
95647c6ae99SBarry Smith +   dm - the DM object
95747c6ae99SBarry Smith -   nlevels - the number of levels of refinement
95847c6ae99SBarry Smith 
95947c6ae99SBarry Smith     Output Parameter:
96047c6ae99SBarry Smith .   dmf - the refined DM hierarchy
96147c6ae99SBarry Smith 
96247c6ae99SBarry Smith     Level: developer
96347c6ae99SBarry Smith 
964e727c939SJed Brown .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
96547c6ae99SBarry Smith 
96647c6ae99SBarry Smith @*/
9677087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
96847c6ae99SBarry Smith {
96947c6ae99SBarry Smith   PetscErrorCode ierr;
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith   PetscFunctionBegin;
972171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
97347c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
97447c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
97547c6ae99SBarry Smith   if (dm->ops->refinehierarchy) {
97647c6ae99SBarry Smith     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
97747c6ae99SBarry Smith   } else if (dm->ops->refine) {
97847c6ae99SBarry Smith     PetscInt i;
97947c6ae99SBarry Smith 
98047c6ae99SBarry Smith     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
98147c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
98247c6ae99SBarry Smith       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
98347c6ae99SBarry Smith     }
98447c6ae99SBarry Smith   } else {
98547c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
98647c6ae99SBarry Smith   }
98747c6ae99SBarry Smith   PetscFunctionReturn(0);
98847c6ae99SBarry Smith }
98947c6ae99SBarry Smith 
99047c6ae99SBarry Smith #undef __FUNCT__
99147c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy"
99247c6ae99SBarry Smith /*@C
99347c6ae99SBarry Smith     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
99447c6ae99SBarry Smith 
99547c6ae99SBarry Smith     Collective on DM
99647c6ae99SBarry Smith 
99747c6ae99SBarry Smith     Input Parameter:
99847c6ae99SBarry Smith +   dm - the DM object
99947c6ae99SBarry Smith -   nlevels - the number of levels of coarsening
100047c6ae99SBarry Smith 
100147c6ae99SBarry Smith     Output Parameter:
100247c6ae99SBarry Smith .   dmc - the coarsened DM hierarchy
100347c6ae99SBarry Smith 
100447c6ae99SBarry Smith     Level: developer
100547c6ae99SBarry Smith 
1006e727c939SJed Brown .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith @*/
10097087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
101047c6ae99SBarry Smith {
101147c6ae99SBarry Smith   PetscErrorCode ierr;
101247c6ae99SBarry Smith 
101347c6ae99SBarry Smith   PetscFunctionBegin;
1014171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
101547c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
101647c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
101747c6ae99SBarry Smith   PetscValidPointer(dmc,3);
101847c6ae99SBarry Smith   if (dm->ops->coarsenhierarchy) {
101947c6ae99SBarry Smith     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
102047c6ae99SBarry Smith   } else if (dm->ops->coarsen) {
102147c6ae99SBarry Smith     PetscInt i;
102247c6ae99SBarry Smith 
102347c6ae99SBarry Smith     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
102447c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
102547c6ae99SBarry Smith       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
102647c6ae99SBarry Smith     }
102747c6ae99SBarry Smith   } else {
102847c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
102947c6ae99SBarry Smith   }
103047c6ae99SBarry Smith   PetscFunctionReturn(0);
103147c6ae99SBarry Smith }
103247c6ae99SBarry Smith 
103347c6ae99SBarry Smith #undef __FUNCT__
1034e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates"
103547c6ae99SBarry Smith /*@
1036e727c939SJed Brown    DMCreateAggregates - Gets the aggregates that map between
103747c6ae99SBarry Smith    grids associated with two DMs.
103847c6ae99SBarry Smith 
103947c6ae99SBarry Smith    Collective on DM
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith    Input Parameters:
104247c6ae99SBarry Smith +  dmc - the coarse grid DM
104347c6ae99SBarry Smith -  dmf - the fine grid DM
104447c6ae99SBarry Smith 
104547c6ae99SBarry Smith    Output Parameters:
104647c6ae99SBarry Smith .  rest - the restriction matrix (transpose of the projection matrix)
104747c6ae99SBarry Smith 
104847c6ae99SBarry Smith    Level: intermediate
104947c6ae99SBarry Smith 
105047c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid
105147c6ae99SBarry Smith 
1052e727c939SJed Brown .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
105347c6ae99SBarry Smith @*/
1054e727c939SJed Brown PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
105547c6ae99SBarry Smith {
105647c6ae99SBarry Smith   PetscErrorCode ierr;
105747c6ae99SBarry Smith 
105847c6ae99SBarry Smith   PetscFunctionBegin;
1059171400e9SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1060171400e9SBarry Smith   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
106147c6ae99SBarry Smith   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
106247c6ae99SBarry Smith   PetscFunctionReturn(0);
106347c6ae99SBarry Smith }
106447c6ae99SBarry Smith 
106547c6ae99SBarry Smith #undef __FUNCT__
10661a266240SBarry Smith #define __FUNCT__ "DMSetApplicationContextDestroy"
10671a266240SBarry Smith /*@C
10681a266240SBarry Smith     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
10691a266240SBarry Smith 
10701a266240SBarry Smith     Not Collective
10711a266240SBarry Smith 
10721a266240SBarry Smith     Input Parameters:
10731a266240SBarry Smith +   dm - the DM object
10741a266240SBarry Smith -   destroy - the destroy function
10751a266240SBarry Smith 
10761a266240SBarry Smith     Level: intermediate
10771a266240SBarry Smith 
1078e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
10791a266240SBarry Smith 
1080f07f9ceaSJed Brown @*/
10811a266240SBarry Smith PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
10821a266240SBarry Smith {
10831a266240SBarry Smith   PetscFunctionBegin;
1084171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
10851a266240SBarry Smith   dm->ctxdestroy = destroy;
10861a266240SBarry Smith   PetscFunctionReturn(0);
10871a266240SBarry Smith }
10881a266240SBarry Smith 
10891a266240SBarry Smith #undef __FUNCT__
10901b2093e4SBarry Smith #define __FUNCT__ "DMSetApplicationContext"
1091b07ff414SBarry Smith /*@
10921b2093e4SBarry Smith     DMSetApplicationContext - Set a user context into a DM object
109347c6ae99SBarry Smith 
109447c6ae99SBarry Smith     Not Collective
109547c6ae99SBarry Smith 
109647c6ae99SBarry Smith     Input Parameters:
109747c6ae99SBarry Smith +   dm - the DM object
109847c6ae99SBarry Smith -   ctx - the user context
109947c6ae99SBarry Smith 
110047c6ae99SBarry Smith     Level: intermediate
110147c6ae99SBarry Smith 
1102e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
110347c6ae99SBarry Smith 
110447c6ae99SBarry Smith @*/
11051b2093e4SBarry Smith PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
110647c6ae99SBarry Smith {
110747c6ae99SBarry Smith   PetscFunctionBegin;
1108171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
110947c6ae99SBarry Smith   dm->ctx = ctx;
111047c6ae99SBarry Smith   PetscFunctionReturn(0);
111147c6ae99SBarry Smith }
111247c6ae99SBarry Smith 
111347c6ae99SBarry Smith #undef __FUNCT__
11141b2093e4SBarry Smith #define __FUNCT__ "DMGetApplicationContext"
111547c6ae99SBarry Smith /*@
11161b2093e4SBarry Smith     DMGetApplicationContext - Gets a user context from a DM object
111747c6ae99SBarry Smith 
111847c6ae99SBarry Smith     Not Collective
111947c6ae99SBarry Smith 
112047c6ae99SBarry Smith     Input Parameter:
112147c6ae99SBarry Smith .   dm - the DM object
112247c6ae99SBarry Smith 
112347c6ae99SBarry Smith     Output Parameter:
112447c6ae99SBarry Smith .   ctx - the user context
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith     Level: intermediate
112747c6ae99SBarry Smith 
1128e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
112947c6ae99SBarry Smith 
113047c6ae99SBarry Smith @*/
11311b2093e4SBarry Smith PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
113247c6ae99SBarry Smith {
113347c6ae99SBarry Smith   PetscFunctionBegin;
1134171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
11351b2093e4SBarry Smith   *(void**)ctx = dm->ctx;
113647c6ae99SBarry Smith   PetscFunctionReturn(0);
113747c6ae99SBarry Smith }
113847c6ae99SBarry Smith 
113947c6ae99SBarry Smith #undef __FUNCT__
114047c6ae99SBarry Smith #define __FUNCT__ "DMSetInitialGuess"
11417e833e3aSBarry Smith /*@C
114247c6ae99SBarry Smith     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
114347c6ae99SBarry Smith 
114447c6ae99SBarry Smith     Logically Collective on DM
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith     Input Parameter:
114747c6ae99SBarry Smith +   dm - the DM object to destroy
114847c6ae99SBarry Smith -   f - the function to compute the initial guess
114947c6ae99SBarry Smith 
115047c6ae99SBarry Smith     Level: intermediate
115147c6ae99SBarry Smith 
1152e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
115347c6ae99SBarry Smith 
1154f07f9ceaSJed Brown @*/
11557087cfbeSBarry Smith PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
115647c6ae99SBarry Smith {
115747c6ae99SBarry Smith   PetscFunctionBegin;
1158171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
115947c6ae99SBarry Smith   dm->ops->initialguess = f;
116047c6ae99SBarry Smith   PetscFunctionReturn(0);
116147c6ae99SBarry Smith }
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith #undef __FUNCT__
116447c6ae99SBarry Smith #define __FUNCT__ "DMSetFunction"
11657e833e3aSBarry Smith /*@C
116647c6ae99SBarry Smith     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
116747c6ae99SBarry Smith 
116847c6ae99SBarry Smith     Logically Collective on DM
116947c6ae99SBarry Smith 
117047c6ae99SBarry Smith     Input Parameter:
117147c6ae99SBarry Smith +   dm - the DM object
117247c6ae99SBarry Smith -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
117347c6ae99SBarry Smith 
117447c6ae99SBarry Smith     Level: intermediate
117547c6ae99SBarry Smith 
117647c6ae99SBarry Smith     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
117747c6ae99SBarry Smith            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
117847c6ae99SBarry Smith 
1179e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
118047c6ae99SBarry Smith          DMSetJacobian()
118147c6ae99SBarry Smith 
1182f07f9ceaSJed Brown @*/
11837087cfbeSBarry Smith PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
118447c6ae99SBarry Smith {
118547c6ae99SBarry Smith   PetscFunctionBegin;
1186171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
118747c6ae99SBarry Smith   dm->ops->function = f;
118847c6ae99SBarry Smith   if (f) {
118947c6ae99SBarry Smith     dm->ops->functionj = f;
119047c6ae99SBarry Smith   }
119147c6ae99SBarry Smith   PetscFunctionReturn(0);
119247c6ae99SBarry Smith }
119347c6ae99SBarry Smith 
119447c6ae99SBarry Smith #undef __FUNCT__
119547c6ae99SBarry Smith #define __FUNCT__ "DMSetJacobian"
11967e833e3aSBarry Smith /*@C
119747c6ae99SBarry Smith     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
119847c6ae99SBarry Smith 
119947c6ae99SBarry Smith     Logically Collective on DM
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith     Input Parameter:
120247c6ae99SBarry Smith +   dm - the DM object to destroy
120347c6ae99SBarry Smith -   f - the function to compute the matrix entries
120447c6ae99SBarry Smith 
120547c6ae99SBarry Smith     Level: intermediate
120647c6ae99SBarry Smith 
1207e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
120847c6ae99SBarry Smith          DMSetFunction()
120947c6ae99SBarry Smith 
1210f07f9ceaSJed Brown @*/
12117087cfbeSBarry Smith PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
121247c6ae99SBarry Smith {
121347c6ae99SBarry Smith   PetscFunctionBegin;
1214171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
121547c6ae99SBarry Smith   dm->ops->jacobian = f;
121647c6ae99SBarry Smith   PetscFunctionReturn(0);
121747c6ae99SBarry Smith }
121847c6ae99SBarry Smith 
121947c6ae99SBarry Smith #undef __FUNCT__
122047c6ae99SBarry Smith #define __FUNCT__ "DMComputeInitialGuess"
122147c6ae99SBarry Smith /*@
122247c6ae99SBarry Smith     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
122347c6ae99SBarry Smith 
122447c6ae99SBarry Smith     Collective on DM
122547c6ae99SBarry Smith 
122647c6ae99SBarry Smith     Input Parameter:
122747c6ae99SBarry Smith +   dm - the DM object to destroy
122847c6ae99SBarry Smith -   x - the vector to hold the initial guess values
122947c6ae99SBarry Smith 
123047c6ae99SBarry Smith     Level: developer
123147c6ae99SBarry Smith 
1232e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
123347c6ae99SBarry Smith 
123447c6ae99SBarry Smith @*/
12357087cfbeSBarry Smith PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
123647c6ae99SBarry Smith {
123747c6ae99SBarry Smith   PetscErrorCode ierr;
123847c6ae99SBarry Smith   PetscFunctionBegin;
1239171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
124047c6ae99SBarry Smith   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
124147c6ae99SBarry Smith   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
124247c6ae99SBarry Smith   PetscFunctionReturn(0);
124347c6ae99SBarry Smith }
124447c6ae99SBarry Smith 
124547c6ae99SBarry Smith #undef __FUNCT__
124647c6ae99SBarry Smith #define __FUNCT__ "DMHasInitialGuess"
124747c6ae99SBarry Smith /*@
124847c6ae99SBarry Smith     DMHasInitialGuess - does the DM object have an initial guess function
124947c6ae99SBarry Smith 
125047c6ae99SBarry Smith     Not Collective
125147c6ae99SBarry Smith 
125247c6ae99SBarry Smith     Input Parameter:
125347c6ae99SBarry Smith .   dm - the DM object to destroy
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith     Output Parameter:
125647c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
125747c6ae99SBarry Smith 
125847c6ae99SBarry Smith     Level: developer
125947c6ae99SBarry Smith 
1260e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
126147c6ae99SBarry Smith 
126247c6ae99SBarry Smith @*/
12637087cfbeSBarry Smith PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
126447c6ae99SBarry Smith {
126547c6ae99SBarry Smith   PetscFunctionBegin;
1266171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
126747c6ae99SBarry Smith   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
126847c6ae99SBarry Smith   PetscFunctionReturn(0);
126947c6ae99SBarry Smith }
127047c6ae99SBarry Smith 
127147c6ae99SBarry Smith #undef __FUNCT__
127247c6ae99SBarry Smith #define __FUNCT__ "DMHasFunction"
127347c6ae99SBarry Smith /*@
127447c6ae99SBarry Smith     DMHasFunction - does the DM object have a function
127547c6ae99SBarry Smith 
127647c6ae99SBarry Smith     Not Collective
127747c6ae99SBarry Smith 
127847c6ae99SBarry Smith     Input Parameter:
127947c6ae99SBarry Smith .   dm - the DM object to destroy
128047c6ae99SBarry Smith 
128147c6ae99SBarry Smith     Output Parameter:
128247c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
128347c6ae99SBarry Smith 
128447c6ae99SBarry Smith     Level: developer
128547c6ae99SBarry Smith 
1286e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
128747c6ae99SBarry Smith 
128847c6ae99SBarry Smith @*/
12897087cfbeSBarry Smith PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
129047c6ae99SBarry Smith {
129147c6ae99SBarry Smith   PetscFunctionBegin;
1292171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
129347c6ae99SBarry Smith   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
129447c6ae99SBarry Smith   PetscFunctionReturn(0);
129547c6ae99SBarry Smith }
129647c6ae99SBarry Smith 
129747c6ae99SBarry Smith #undef __FUNCT__
129847c6ae99SBarry Smith #define __FUNCT__ "DMHasJacobian"
129947c6ae99SBarry Smith /*@
130047c6ae99SBarry Smith     DMHasJacobian - does the DM object have a matrix function
130147c6ae99SBarry Smith 
130247c6ae99SBarry Smith     Not Collective
130347c6ae99SBarry Smith 
130447c6ae99SBarry Smith     Input Parameter:
130547c6ae99SBarry Smith .   dm - the DM object to destroy
130647c6ae99SBarry Smith 
130747c6ae99SBarry Smith     Output Parameter:
130847c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
130947c6ae99SBarry Smith 
131047c6ae99SBarry Smith     Level: developer
131147c6ae99SBarry Smith 
1312e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
131347c6ae99SBarry Smith 
131447c6ae99SBarry Smith @*/
13157087cfbeSBarry Smith PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
131647c6ae99SBarry Smith {
131747c6ae99SBarry Smith   PetscFunctionBegin;
1318171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
131947c6ae99SBarry Smith   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
132047c6ae99SBarry Smith   PetscFunctionReturn(0);
132147c6ae99SBarry Smith }
132247c6ae99SBarry Smith 
132347c6ae99SBarry Smith #undef __FUNCT__
132447c6ae99SBarry Smith #define __FUNCT__ "DMComputeFunction"
132547c6ae99SBarry Smith /*@
132647c6ae99SBarry Smith     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
132747c6ae99SBarry Smith 
132847c6ae99SBarry Smith     Collective on DM
132947c6ae99SBarry Smith 
133047c6ae99SBarry Smith     Input Parameter:
133147c6ae99SBarry Smith +   dm - the DM object to destroy
133247c6ae99SBarry Smith .   x - the location where the function is evaluationed, may be ignored for linear problems
133347c6ae99SBarry Smith -   b - the vector to hold the right hand side entries
133447c6ae99SBarry Smith 
133547c6ae99SBarry Smith     Level: developer
133647c6ae99SBarry Smith 
1337e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
133847c6ae99SBarry Smith          DMSetJacobian()
133947c6ae99SBarry Smith 
134047c6ae99SBarry Smith @*/
13417087cfbeSBarry Smith PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
134247c6ae99SBarry Smith {
134347c6ae99SBarry Smith   PetscErrorCode ierr;
134447c6ae99SBarry Smith   PetscFunctionBegin;
1345171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
134647c6ae99SBarry Smith   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1347644e2e5bSBarry Smith   PetscStackPush("DM user function");
134847c6ae99SBarry Smith   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1349644e2e5bSBarry Smith   PetscStackPop;
135047c6ae99SBarry Smith   PetscFunctionReturn(0);
135147c6ae99SBarry Smith }
135247c6ae99SBarry Smith 
135347c6ae99SBarry Smith 
135447c6ae99SBarry Smith #undef __FUNCT__
135547c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobian"
135647c6ae99SBarry Smith /*@
135747c6ae99SBarry Smith     DMComputeJacobian - compute the matrix entries for the solver
135847c6ae99SBarry Smith 
135947c6ae99SBarry Smith     Collective on DM
136047c6ae99SBarry Smith 
136147c6ae99SBarry Smith     Input Parameter:
136247c6ae99SBarry Smith +   dm - the DM object
1363cab2e9ccSBarry Smith .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
136447c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
136547c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
136647c6ae99SBarry Smith 
136747c6ae99SBarry Smith     Level: developer
136847c6ae99SBarry Smith 
1369e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
137047c6ae99SBarry Smith          DMSetFunction()
137147c6ae99SBarry Smith 
137247c6ae99SBarry Smith @*/
13737087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
137447c6ae99SBarry Smith {
137547c6ae99SBarry Smith   PetscErrorCode ierr;
137647c6ae99SBarry Smith 
137747c6ae99SBarry Smith   PetscFunctionBegin;
1378171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
137947c6ae99SBarry Smith   if (!dm->ops->jacobian) {
138047c6ae99SBarry Smith     ISColoring     coloring;
138147c6ae99SBarry Smith     MatFDColoring  fd;
138247c6ae99SBarry Smith 
1383171400e9SBarry Smith     ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
138447c6ae99SBarry Smith     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1385fcfd50ebSBarry Smith     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
138647c6ae99SBarry Smith     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
13870bdded8aSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
13880bdded8aSJed Brown     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
138971cd77b2SBarry Smith 
139047c6ae99SBarry Smith     dm->fd = fd;
139147c6ae99SBarry Smith     dm->ops->jacobian = DMComputeJacobianDefault;
13922533e041SBarry Smith 
139371cd77b2SBarry Smith     /* don't know why this is needed */
139471cd77b2SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
139547c6ae99SBarry Smith   }
139647c6ae99SBarry Smith   if (!x) x = dm->x;
139747c6ae99SBarry Smith   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1398cab2e9ccSBarry Smith 
139971cd77b2SBarry Smith   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1400649052a6SBarry Smith   if (x) {
1401cab2e9ccSBarry Smith     if (!dm->x) {
140271cd77b2SBarry Smith       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1403cab2e9ccSBarry Smith     }
1404cab2e9ccSBarry Smith     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1405649052a6SBarry Smith   }
1406a8248277SBarry Smith   if (A != B) {
1407a8248277SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1408a8248277SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1409a8248277SBarry Smith   }
141047c6ae99SBarry Smith   PetscFunctionReturn(0);
141147c6ae99SBarry Smith }
1412264ace61SBarry Smith 
1413cab2e9ccSBarry Smith 
1414264ace61SBarry Smith PetscFList DMList                       = PETSC_NULL;
1415264ace61SBarry Smith PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1416264ace61SBarry Smith 
1417264ace61SBarry Smith #undef __FUNCT__
1418264ace61SBarry Smith #define __FUNCT__ "DMSetType"
1419264ace61SBarry Smith /*@C
1420264ace61SBarry Smith   DMSetType - Builds a DM, for a particular DM implementation.
1421264ace61SBarry Smith 
1422264ace61SBarry Smith   Collective on DM
1423264ace61SBarry Smith 
1424264ace61SBarry Smith   Input Parameters:
1425264ace61SBarry Smith + dm     - The DM object
1426264ace61SBarry Smith - method - The name of the DM type
1427264ace61SBarry Smith 
1428264ace61SBarry Smith   Options Database Key:
1429264ace61SBarry Smith . -dm_type <type> - Sets the DM type; use -help for a list of available types
1430264ace61SBarry Smith 
1431264ace61SBarry Smith   Notes:
1432e1589f56SBarry Smith   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1433264ace61SBarry Smith 
1434264ace61SBarry Smith   Level: intermediate
1435264ace61SBarry Smith 
1436264ace61SBarry Smith .keywords: DM, set, type
1437264ace61SBarry Smith .seealso: DMGetType(), DMCreate()
1438264ace61SBarry Smith @*/
14397087cfbeSBarry Smith PetscErrorCode  DMSetType(DM dm, const DMType method)
1440264ace61SBarry Smith {
1441264ace61SBarry Smith   PetscErrorCode (*r)(DM);
1442264ace61SBarry Smith   PetscBool      match;
1443264ace61SBarry Smith   PetscErrorCode ierr;
1444264ace61SBarry Smith 
1445264ace61SBarry Smith   PetscFunctionBegin;
1446264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1447264ace61SBarry Smith   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1448264ace61SBarry Smith   if (match) PetscFunctionReturn(0);
1449264ace61SBarry Smith 
1450264ace61SBarry Smith   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
14514b91b6eaSBarry Smith   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1452264ace61SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1453264ace61SBarry Smith 
1454264ace61SBarry Smith   if (dm->ops->destroy) {
1455264ace61SBarry Smith     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1456264ace61SBarry Smith   }
1457264ace61SBarry Smith   ierr = (*r)(dm);CHKERRQ(ierr);
1458264ace61SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1459264ace61SBarry Smith   PetscFunctionReturn(0);
1460264ace61SBarry Smith }
1461264ace61SBarry Smith 
1462264ace61SBarry Smith #undef __FUNCT__
1463264ace61SBarry Smith #define __FUNCT__ "DMGetType"
1464264ace61SBarry Smith /*@C
1465264ace61SBarry Smith   DMGetType - Gets the DM type name (as a string) from the DM.
1466264ace61SBarry Smith 
1467264ace61SBarry Smith   Not Collective
1468264ace61SBarry Smith 
1469264ace61SBarry Smith   Input Parameter:
1470264ace61SBarry Smith . dm  - The DM
1471264ace61SBarry Smith 
1472264ace61SBarry Smith   Output Parameter:
1473264ace61SBarry Smith . type - The DM type name
1474264ace61SBarry Smith 
1475264ace61SBarry Smith   Level: intermediate
1476264ace61SBarry Smith 
1477264ace61SBarry Smith .keywords: DM, get, type, name
1478264ace61SBarry Smith .seealso: DMSetType(), DMCreate()
1479264ace61SBarry Smith @*/
14807087cfbeSBarry Smith PetscErrorCode  DMGetType(DM dm, const DMType *type)
1481264ace61SBarry Smith {
1482264ace61SBarry Smith   PetscErrorCode ierr;
1483264ace61SBarry Smith 
1484264ace61SBarry Smith   PetscFunctionBegin;
1485264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1486264ace61SBarry Smith   PetscValidCharPointer(type,2);
1487264ace61SBarry Smith   if (!DMRegisterAllCalled) {
1488264ace61SBarry Smith     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1489264ace61SBarry Smith   }
1490264ace61SBarry Smith   *type = ((PetscObject)dm)->type_name;
1491264ace61SBarry Smith   PetscFunctionReturn(0);
1492264ace61SBarry Smith }
1493264ace61SBarry Smith 
149467a56275SMatthew G Knepley #undef __FUNCT__
149567a56275SMatthew G Knepley #define __FUNCT__ "DMConvert"
149667a56275SMatthew G Knepley /*@C
149767a56275SMatthew G Knepley   DMConvert - Converts a DM to another DM, either of the same or different type.
149867a56275SMatthew G Knepley 
149967a56275SMatthew G Knepley   Collective on DM
150067a56275SMatthew G Knepley 
150167a56275SMatthew G Knepley   Input Parameters:
150267a56275SMatthew G Knepley + dm - the DM
150367a56275SMatthew G Knepley - newtype - new DM type (use "same" for the same type)
150467a56275SMatthew G Knepley 
150567a56275SMatthew G Knepley   Output Parameter:
150667a56275SMatthew G Knepley . M - pointer to new DM
150767a56275SMatthew G Knepley 
150867a56275SMatthew G Knepley   Notes:
150967a56275SMatthew G Knepley   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
151067a56275SMatthew G Knepley   the MPI communicator of the generated DM is always the same as the communicator
151167a56275SMatthew G Knepley   of the input DM.
151267a56275SMatthew G Knepley 
151367a56275SMatthew G Knepley   Level: intermediate
151467a56275SMatthew G Knepley 
151567a56275SMatthew G Knepley .seealso: DMCreate()
151667a56275SMatthew G Knepley @*/
151767a56275SMatthew G Knepley PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
151867a56275SMatthew G Knepley {
151967a56275SMatthew G Knepley   DM             B;
152067a56275SMatthew G Knepley   char           convname[256];
152167a56275SMatthew G Knepley   PetscBool      sametype, issame;
152267a56275SMatthew G Knepley   PetscErrorCode ierr;
152367a56275SMatthew G Knepley 
152467a56275SMatthew G Knepley   PetscFunctionBegin;
152567a56275SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
152667a56275SMatthew G Knepley   PetscValidType(dm,1);
152767a56275SMatthew G Knepley   PetscValidPointer(M,3);
152867a56275SMatthew G Knepley   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
152967a56275SMatthew G Knepley   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
153067a56275SMatthew G Knepley   {
153167a56275SMatthew G Knepley     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
153267a56275SMatthew G Knepley 
153367a56275SMatthew G Knepley     /*
153467a56275SMatthew G Knepley        Order of precedence:
153567a56275SMatthew G Knepley        1) See if a specialized converter is known to the current DM.
153667a56275SMatthew G Knepley        2) See if a specialized converter is known to the desired DM class.
153767a56275SMatthew G Knepley        3) See if a good general converter is registered for the desired class
153867a56275SMatthew G Knepley        4) See if a good general converter is known for the current matrix.
153967a56275SMatthew G Knepley        5) Use a really basic converter.
154067a56275SMatthew G Knepley     */
154167a56275SMatthew G Knepley 
154267a56275SMatthew G Knepley     /* 1) See if a specialized converter is known to the current DM and the desired class */
154367a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
154467a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
154567a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
154667a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
154767a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
154867a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
154967a56275SMatthew G Knepley     if (conv) goto foundconv;
155067a56275SMatthew G Knepley 
155167a56275SMatthew G Knepley     /* 2)  See if a specialized converter is known to the desired DM class. */
155267a56275SMatthew G Knepley     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
155367a56275SMatthew G Knepley     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
155467a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
155567a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
155667a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
155767a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
155867a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
155967a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
156067a56275SMatthew G Knepley     if (conv) {
1561fcfd50ebSBarry Smith       ierr = DMDestroy(&B);CHKERRQ(ierr);
156267a56275SMatthew G Knepley       goto foundconv;
156367a56275SMatthew G Knepley     }
156467a56275SMatthew G Knepley 
156567a56275SMatthew G Knepley #if 0
156667a56275SMatthew G Knepley     /* 3) See if a good general converter is registered for the desired class */
156767a56275SMatthew G Knepley     conv = B->ops->convertfrom;
1568fcfd50ebSBarry Smith     ierr = DMDestroy(&B);CHKERRQ(ierr);
156967a56275SMatthew G Knepley     if (conv) goto foundconv;
157067a56275SMatthew G Knepley 
157167a56275SMatthew G Knepley     /* 4) See if a good general converter is known for the current matrix */
157267a56275SMatthew G Knepley     if (dm->ops->convert) {
157367a56275SMatthew G Knepley       conv = dm->ops->convert;
157467a56275SMatthew G Knepley     }
157567a56275SMatthew G Knepley     if (conv) goto foundconv;
157667a56275SMatthew G Knepley #endif
157767a56275SMatthew G Knepley 
157867a56275SMatthew G Knepley     /* 5) Use a really basic converter. */
157967a56275SMatthew G Knepley     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
158067a56275SMatthew G Knepley 
158167a56275SMatthew G Knepley     foundconv:
158267a56275SMatthew G Knepley     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
158367a56275SMatthew G Knepley     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
158467a56275SMatthew G Knepley     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
158567a56275SMatthew G Knepley   }
158667a56275SMatthew G Knepley   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
158767a56275SMatthew G Knepley   PetscFunctionReturn(0);
158867a56275SMatthew G Knepley }
1589264ace61SBarry Smith 
1590264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
1591264ace61SBarry Smith 
1592264ace61SBarry Smith #undef __FUNCT__
1593264ace61SBarry Smith #define __FUNCT__ "DMRegister"
1594264ace61SBarry Smith /*@C
1595264ace61SBarry Smith   DMRegister - See DMRegisterDynamic()
1596264ace61SBarry Smith 
1597264ace61SBarry Smith   Level: advanced
1598264ace61SBarry Smith @*/
15997087cfbeSBarry Smith PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1600264ace61SBarry Smith {
1601264ace61SBarry Smith   char fullname[PETSC_MAX_PATH_LEN];
1602264ace61SBarry Smith   PetscErrorCode ierr;
1603264ace61SBarry Smith 
1604264ace61SBarry Smith   PetscFunctionBegin;
1605264ace61SBarry Smith   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1606264ace61SBarry Smith   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1607264ace61SBarry Smith   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1608264ace61SBarry Smith   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1609264ace61SBarry Smith   PetscFunctionReturn(0);
1610264ace61SBarry Smith }
1611264ace61SBarry Smith 
1612264ace61SBarry Smith 
1613264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
1614264ace61SBarry Smith #undef __FUNCT__
1615264ace61SBarry Smith #define __FUNCT__ "DMRegisterDestroy"
1616264ace61SBarry Smith /*@C
1617264ace61SBarry Smith    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1618264ace61SBarry Smith 
1619264ace61SBarry Smith    Not Collective
1620264ace61SBarry Smith 
1621264ace61SBarry Smith    Level: advanced
1622264ace61SBarry Smith 
1623264ace61SBarry Smith .keywords: DM, register, destroy
1624264ace61SBarry Smith .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1625264ace61SBarry Smith @*/
16267087cfbeSBarry Smith PetscErrorCode  DMRegisterDestroy(void)
1627264ace61SBarry Smith {
1628264ace61SBarry Smith   PetscErrorCode ierr;
1629264ace61SBarry Smith 
1630264ace61SBarry Smith   PetscFunctionBegin;
1631264ace61SBarry Smith   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1632264ace61SBarry Smith   DMRegisterAllCalled = PETSC_FALSE;
1633264ace61SBarry Smith   PetscFunctionReturn(0);
1634264ace61SBarry Smith }
163523f975d1SBarry Smith 
163623f975d1SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1637c6db04a5SJed Brown #include <mex.h>
163823f975d1SBarry Smith 
16393014e516SBarry Smith typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
164023f975d1SBarry Smith 
164123f975d1SBarry Smith #undef __FUNCT__
164223f975d1SBarry Smith #define __FUNCT__ "DMComputeFunction_Matlab"
164323f975d1SBarry Smith /*
164423f975d1SBarry Smith    DMComputeFunction_Matlab - Calls the function that has been set with
164523f975d1SBarry Smith                          DMSetFunctionMatlab().
164623f975d1SBarry Smith 
164723f975d1SBarry Smith    For linear problems x is null
164823f975d1SBarry Smith 
164923f975d1SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
165023f975d1SBarry Smith */
16517087cfbeSBarry Smith PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
165223f975d1SBarry Smith {
165323f975d1SBarry Smith   PetscErrorCode    ierr;
165423f975d1SBarry Smith   DMMatlabContext   *sctx;
165523f975d1SBarry Smith   int               nlhs = 1,nrhs = 4;
165623f975d1SBarry Smith   mxArray	    *plhs[1],*prhs[4];
165723f975d1SBarry Smith   long long int     lx = 0,ly = 0,ls = 0;
165823f975d1SBarry Smith 
165923f975d1SBarry Smith   PetscFunctionBegin;
166023f975d1SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
166123f975d1SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
166223f975d1SBarry Smith   PetscCheckSameComm(dm,1,y,3);
166323f975d1SBarry Smith 
166423f975d1SBarry Smith   /* call Matlab function in ctx with arguments x and y */
16651b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
166623f975d1SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
166723f975d1SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
16683014e516SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
166923f975d1SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
167023f975d1SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
167123f975d1SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)ly);
167223f975d1SBarry Smith   prhs[3] =  mxCreateString(sctx->funcname);
1673b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
167423f975d1SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
167523f975d1SBarry Smith   mxDestroyArray(prhs[0]);
167623f975d1SBarry Smith   mxDestroyArray(prhs[1]);
167723f975d1SBarry Smith   mxDestroyArray(prhs[2]);
167823f975d1SBarry Smith   mxDestroyArray(prhs[3]);
167923f975d1SBarry Smith   mxDestroyArray(plhs[0]);
168023f975d1SBarry Smith   PetscFunctionReturn(0);
168123f975d1SBarry Smith }
168223f975d1SBarry Smith 
168323f975d1SBarry Smith 
168423f975d1SBarry Smith #undef __FUNCT__
168523f975d1SBarry Smith #define __FUNCT__ "DMSetFunctionMatlab"
168623f975d1SBarry Smith /*
168723f975d1SBarry Smith    DMSetFunctionMatlab - Sets the function evaluation routine
168823f975d1SBarry Smith 
168923f975d1SBarry Smith */
16907087cfbeSBarry Smith PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
169123f975d1SBarry Smith {
169223f975d1SBarry Smith   PetscErrorCode    ierr;
169323f975d1SBarry Smith   DMMatlabContext   *sctx;
169423f975d1SBarry Smith 
169523f975d1SBarry Smith   PetscFunctionBegin;
1696171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
169723f975d1SBarry Smith   /* currently sctx is memory bleed */
16981b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
16993014e516SBarry Smith   if (!sctx) {
170023f975d1SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
17013014e516SBarry Smith   }
170223f975d1SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
17031b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
170423f975d1SBarry Smith   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
170523f975d1SBarry Smith   PetscFunctionReturn(0);
170623f975d1SBarry Smith }
17073014e516SBarry Smith 
17083014e516SBarry Smith #undef __FUNCT__
17093014e516SBarry Smith #define __FUNCT__ "DMComputeJacobian_Matlab"
17103014e516SBarry Smith /*
17113014e516SBarry Smith    DMComputeJacobian_Matlab - Calls the function that has been set with
17123014e516SBarry Smith                          DMSetJacobianMatlab().
17133014e516SBarry Smith 
17143014e516SBarry Smith    For linear problems x is null
17153014e516SBarry Smith 
17163014e516SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
17173014e516SBarry Smith */
17187087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
17193014e516SBarry Smith {
17203014e516SBarry Smith   PetscErrorCode    ierr;
17213014e516SBarry Smith   DMMatlabContext   *sctx;
17223014e516SBarry Smith   int               nlhs = 2,nrhs = 5;
17233014e516SBarry Smith   mxArray	    *plhs[2],*prhs[5];
17243014e516SBarry Smith   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
17253014e516SBarry Smith 
17263014e516SBarry Smith   PetscFunctionBegin;
17273014e516SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
17283014e516SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
17293014e516SBarry Smith 
1730e3c5b3baSBarry Smith   /* call MATLAB function in ctx with arguments x, A, and B */
17311b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
17323014e516SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
17333014e516SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
17343014e516SBarry Smith   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
17353014e516SBarry Smith   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
17363014e516SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
17373014e516SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
17383014e516SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lA);
17393014e516SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lB);
17403014e516SBarry Smith   prhs[4] =  mxCreateString(sctx->jacname);
1741b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
1742c980e822SBarry Smith   *str    =  (MatStructure) mxGetScalar(plhs[0]);
1743c088a8dcSBarry Smith   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
17443014e516SBarry Smith   mxDestroyArray(prhs[0]);
17453014e516SBarry Smith   mxDestroyArray(prhs[1]);
17463014e516SBarry Smith   mxDestroyArray(prhs[2]);
17473014e516SBarry Smith   mxDestroyArray(prhs[3]);
17483014e516SBarry Smith   mxDestroyArray(prhs[4]);
17493014e516SBarry Smith   mxDestroyArray(plhs[0]);
17503014e516SBarry Smith   mxDestroyArray(plhs[1]);
17513014e516SBarry Smith   PetscFunctionReturn(0);
17523014e516SBarry Smith }
17533014e516SBarry Smith 
17543014e516SBarry Smith 
17553014e516SBarry Smith #undef __FUNCT__
17563014e516SBarry Smith #define __FUNCT__ "DMSetJacobianMatlab"
17573014e516SBarry Smith /*
17583014e516SBarry Smith    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
17593014e516SBarry Smith 
17603014e516SBarry Smith */
17617087cfbeSBarry Smith PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
17623014e516SBarry Smith {
17633014e516SBarry Smith   PetscErrorCode    ierr;
17643014e516SBarry Smith   DMMatlabContext   *sctx;
17653014e516SBarry Smith 
17663014e516SBarry Smith   PetscFunctionBegin;
1767171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
17683014e516SBarry Smith   /* currently sctx is memory bleed */
17691b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
17703014e516SBarry Smith   if (!sctx) {
17713014e516SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
17723014e516SBarry Smith   }
17733014e516SBarry Smith   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
17741b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
17753014e516SBarry Smith   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
17763014e516SBarry Smith   PetscFunctionReturn(0);
17773014e516SBarry Smith }
177823f975d1SBarry Smith #endif
1779b859378eSBarry Smith 
1780b859378eSBarry Smith #undef __FUNCT__
1781b859378eSBarry Smith #define __FUNCT__ "DMLoad"
1782b859378eSBarry Smith /*@C
1783b859378eSBarry Smith   DMLoad - Loads a DM that has been stored in binary or HDF5 format
1784b859378eSBarry Smith   with DMView().
1785b859378eSBarry Smith 
1786b859378eSBarry Smith   Collective on PetscViewer
1787b859378eSBarry Smith 
1788b859378eSBarry Smith   Input Parameters:
1789b859378eSBarry Smith + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
1790b859378eSBarry Smith            some related function before a call to DMLoad().
1791b859378eSBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
1792b859378eSBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
1793b859378eSBarry Smith 
1794b859378eSBarry Smith    Level: intermediate
1795b859378eSBarry Smith 
1796b859378eSBarry Smith   Notes:
1797b859378eSBarry Smith   Defaults to the DM DA.
1798b859378eSBarry Smith 
1799b859378eSBarry Smith   Notes for advanced users:
1800b859378eSBarry Smith   Most users should not need to know the details of the binary storage
1801b859378eSBarry Smith   format, since DMLoad() and DMView() completely hide these details.
1802b859378eSBarry Smith   But for anyone who's interested, the standard binary matrix storage
1803b859378eSBarry Smith   format is
1804b859378eSBarry Smith .vb
1805b859378eSBarry Smith      has not yet been determined
1806b859378eSBarry Smith .ve
1807b859378eSBarry Smith 
1808b859378eSBarry Smith    In addition, PETSc automatically does the byte swapping for
1809b859378eSBarry Smith machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
1810b859378eSBarry Smith linux, Windows and the paragon; thus if you write your own binary
1811b859378eSBarry Smith read/write routines you have to swap the bytes; see PetscBinaryRead()
1812b859378eSBarry Smith and PetscBinaryWrite() to see how this may be done.
1813b859378eSBarry Smith 
1814b859378eSBarry Smith   Concepts: vector^loading from file
1815b859378eSBarry Smith 
1816b859378eSBarry Smith .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
1817b859378eSBarry Smith @*/
1818b859378eSBarry Smith PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
1819b859378eSBarry Smith {
1820b859378eSBarry Smith   PetscErrorCode ierr;
1821b859378eSBarry Smith 
1822b859378eSBarry Smith   PetscFunctionBegin;
1823b859378eSBarry Smith   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
1824b859378eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1825b859378eSBarry Smith 
1826b859378eSBarry Smith   if (!((PetscObject)newdm)->type_name) {
1827b859378eSBarry Smith     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
1828b859378eSBarry Smith   }
1829b859378eSBarry Smith   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1830b859378eSBarry Smith   PetscFunctionReturn(0);
1831b859378eSBarry Smith }
1832b859378eSBarry Smith 
1833