xref: /petsc/src/dm/interface/dm.c (revision 08da532b2310464c3bf489aef60f0ccd6d61141d)
1*08da532bSDmitry Karpeev #include <petscsnes.h>
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;
135b17ce1afSJed Brown   DMCoarsenHookLink link,next;
13647c6ae99SBarry Smith   PetscErrorCode ierr;
13747c6ae99SBarry Smith 
13847c6ae99SBarry Smith   PetscFunctionBegin;
1396bf464f9SBarry Smith   if (!*dm) PetscFunctionReturn(0);
1406bf464f9SBarry Smith   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
14187e657c6SBarry Smith 
14287e657c6SBarry Smith   /* count all the circular references of DM and its contained Vecs */
143732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1446bf464f9SBarry Smith     if ((*dm)->localin[i])  {cnt++;}
1456bf464f9SBarry Smith     if ((*dm)->globalin[i]) {cnt++;}
146732e2eb9SMatthew G Knepley   }
14771cd77b2SBarry Smith   if ((*dm)->x) {
14871cd77b2SBarry Smith     PetscObject obj;
14971cd77b2SBarry Smith     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
15071cd77b2SBarry Smith     if (obj == (PetscObject)*dm) cnt++;
15171cd77b2SBarry Smith   }
152732e2eb9SMatthew G Knepley 
1536bf464f9SBarry Smith   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
154732e2eb9SMatthew G Knepley   /*
155732e2eb9SMatthew G Knepley      Need this test because the dm references the vectors that
156732e2eb9SMatthew G Knepley      reference the dm, so destroying the dm calls destroy on the
157732e2eb9SMatthew G Knepley      vectors that cause another destroy on the dm
158732e2eb9SMatthew G Knepley   */
1596bf464f9SBarry Smith   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
1606bf464f9SBarry Smith   ((PetscObject) (*dm))->refct = 0;
161732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1626bf464f9SBarry Smith     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
1636bf464f9SBarry Smith     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
164732e2eb9SMatthew G Knepley   }
1651a266240SBarry Smith 
166b17ce1afSJed Brown   /* Destroy the list of hooks */
167b17ce1afSJed Brown   for (link=(*dm)->coarsenhook; link; link=next) {
168b17ce1afSJed Brown     next = link->next;
169b17ce1afSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
170b17ce1afSJed Brown   }
171b17ce1afSJed Brown   (*dm)->coarsenhook = PETSC_NULL;
172b17ce1afSJed Brown 
1731a266240SBarry Smith   if ((*dm)->ctx && (*dm)->ctxdestroy) {
1741a266240SBarry Smith     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
1751a266240SBarry Smith   }
17687e657c6SBarry Smith   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
17771cd77b2SBarry Smith   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
1784dcab191SBarry Smith   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
1796bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
1806bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
1816bf464f9SBarry Smith   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
182073dac72SJed Brown   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
183a89ea682SMatthew G Knepley   ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr);
184732e2eb9SMatthew G Knepley   /* if memory was published with AMS then destroy it */
1856bf464f9SBarry Smith   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
186732e2eb9SMatthew G Knepley 
1876bf464f9SBarry Smith   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
1886bf464f9SBarry Smith   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
189732e2eb9SMatthew G Knepley   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
19047c6ae99SBarry Smith   PetscFunctionReturn(0);
19147c6ae99SBarry Smith }
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith #undef __FUNCT__
194d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp"
195d7bf68aeSBarry Smith /*@
196d7bf68aeSBarry Smith     DMSetUp - sets up the data structures inside a DM object
197d7bf68aeSBarry Smith 
198d7bf68aeSBarry Smith     Collective on DM
199d7bf68aeSBarry Smith 
200d7bf68aeSBarry Smith     Input Parameter:
201d7bf68aeSBarry Smith .   dm - the DM object to setup
202d7bf68aeSBarry Smith 
203d7bf68aeSBarry Smith     Level: developer
204d7bf68aeSBarry Smith 
205e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
206d7bf68aeSBarry Smith 
207d7bf68aeSBarry Smith @*/
2087087cfbeSBarry Smith PetscErrorCode  DMSetUp(DM dm)
209d7bf68aeSBarry Smith {
210d7bf68aeSBarry Smith   PetscErrorCode ierr;
211d7bf68aeSBarry Smith 
212d7bf68aeSBarry Smith   PetscFunctionBegin;
213171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2148387afaaSJed Brown   if (dm->setupcalled) PetscFunctionReturn(0);
215d7bf68aeSBarry Smith   if (dm->ops->setup) {
216d7bf68aeSBarry Smith     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
217d7bf68aeSBarry Smith   }
2188387afaaSJed Brown   dm->setupcalled = PETSC_TRUE;
219d7bf68aeSBarry Smith   PetscFunctionReturn(0);
220d7bf68aeSBarry Smith }
221d7bf68aeSBarry Smith 
222d7bf68aeSBarry Smith #undef __FUNCT__
223d7bf68aeSBarry Smith #define __FUNCT__ "DMSetFromOptions"
224d7bf68aeSBarry Smith /*@
225d7bf68aeSBarry Smith     DMSetFromOptions - sets parameters in a DM from the options database
226d7bf68aeSBarry Smith 
227d7bf68aeSBarry Smith     Collective on DM
228d7bf68aeSBarry Smith 
229d7bf68aeSBarry Smith     Input Parameter:
230d7bf68aeSBarry Smith .   dm - the DM object to set options for
231d7bf68aeSBarry Smith 
232732e2eb9SMatthew G Knepley     Options Database:
233dd85299cSBarry Smith +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
234dd85299cSBarry Smith .   -dm_vec_type <type>  type of vector to create inside DM
235171400e9SBarry Smith .   -dm_mat_type <type>  type of matrix to create inside DM
236171400e9SBarry Smith -   -dm_coloring_type <global or ghosted>
237732e2eb9SMatthew G Knepley 
238d7bf68aeSBarry Smith     Level: developer
239d7bf68aeSBarry Smith 
240e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
241d7bf68aeSBarry Smith 
242d7bf68aeSBarry Smith @*/
2437087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions(DM dm)
244d7bf68aeSBarry Smith {
24567ad5babSMatthew G Knepley   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg;
246d7bf68aeSBarry Smith   PetscErrorCode ierr;
247f9ba7244SBarry Smith   char           typeName[256] = MATAIJ;
248d7bf68aeSBarry Smith 
249d7bf68aeSBarry Smith   PetscFunctionBegin;
250171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2513194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
25282fcb398SMatthew G Knepley     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
253c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr);
254c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr);
25567ad5babSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr);
256073dac72SJed 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);
257f9ba7244SBarry Smith     ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
258f9ba7244SBarry Smith     if (flg) {
259f9ba7244SBarry Smith       ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
260f9ba7244SBarry Smith     }
261f9ba7244SBarry Smith     ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr);
262073dac72SJed Brown     if (flg) {
263073dac72SJed Brown       ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
264f9ba7244SBarry Smith       ierr = PetscStrallocpy(typeName,&dm->mattype);CHKERRQ(ierr);
265073dac72SJed Brown     }
2661b89239cSHong 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);
267f9ba7244SBarry Smith     if (dm->ops->setfromoptions) {
268f9ba7244SBarry Smith       ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
269f9ba7244SBarry Smith     }
270f9ba7244SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
271f9ba7244SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
27282fcb398SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
27382fcb398SMatthew G Knepley   if (flg1) {
27482fcb398SMatthew G Knepley     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
27582fcb398SMatthew G Knepley   }
276c4721b0eSMatthew G Knepley   if (flg2) {
277c4721b0eSMatthew G Knepley     PetscViewer viewer;
278c4721b0eSMatthew G Knepley 
279c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
280c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
281c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
282c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
283c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
284c4721b0eSMatthew G Knepley   }
285c4721b0eSMatthew G Knepley   if (flg3) {
286c4721b0eSMatthew G Knepley     PetscViewer viewer;
287c4721b0eSMatthew G Knepley 
288c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
289c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
290c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
291c4721b0eSMatthew G Knepley     ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr);
292c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
293c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
294c4721b0eSMatthew G Knepley   }
29567ad5babSMatthew G Knepley   if (flg4) {
29667ad5babSMatthew G Knepley     PetscViewer viewer;
29767ad5babSMatthew G Knepley 
29867ad5babSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
29967ad5babSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
30067ad5babSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr);
30167ad5babSMatthew G Knepley     ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr);
30267ad5babSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
30367ad5babSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
30467ad5babSMatthew G Knepley   }
305d7bf68aeSBarry Smith   PetscFunctionReturn(0);
306d7bf68aeSBarry Smith }
307d7bf68aeSBarry Smith 
308d7bf68aeSBarry Smith #undef __FUNCT__
30947c6ae99SBarry Smith #define __FUNCT__ "DMView"
310fc9bc008SSatish Balay /*@C
311aa219208SBarry Smith     DMView - Views a vector packer or DMDA.
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith     Collective on DM
31447c6ae99SBarry Smith 
31547c6ae99SBarry Smith     Input Parameter:
31647c6ae99SBarry Smith +   dm - the DM object to view
31747c6ae99SBarry Smith -   v - the viewer
31847c6ae99SBarry Smith 
31947c6ae99SBarry Smith     Level: developer
32047c6ae99SBarry Smith 
321e727c939SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
32247c6ae99SBarry Smith 
32347c6ae99SBarry Smith @*/
3247087cfbeSBarry Smith PetscErrorCode  DMView(DM dm,PetscViewer v)
32547c6ae99SBarry Smith {
32647c6ae99SBarry Smith   PetscErrorCode ierr;
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith   PetscFunctionBegin;
329171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3303014e516SBarry Smith  if (!v) {
3313014e516SBarry Smith     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
3323014e516SBarry Smith   }
3330c010503SBarry Smith   if (dm->ops->view) {
3340c010503SBarry Smith     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
33547c6ae99SBarry Smith   }
33647c6ae99SBarry Smith   PetscFunctionReturn(0);
33747c6ae99SBarry Smith }
33847c6ae99SBarry Smith 
33947c6ae99SBarry Smith #undef __FUNCT__
34047c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector"
34147c6ae99SBarry Smith /*@
342aa219208SBarry Smith     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
34347c6ae99SBarry Smith 
34447c6ae99SBarry Smith     Collective on DM
34547c6ae99SBarry Smith 
34647c6ae99SBarry Smith     Input Parameter:
34747c6ae99SBarry Smith .   dm - the DM object
34847c6ae99SBarry Smith 
34947c6ae99SBarry Smith     Output Parameter:
35047c6ae99SBarry Smith .   vec - the global vector
35147c6ae99SBarry Smith 
352073dac72SJed Brown     Level: beginner
35347c6ae99SBarry Smith 
354e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith @*/
3577087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
35847c6ae99SBarry Smith {
35947c6ae99SBarry Smith   PetscErrorCode ierr;
36047c6ae99SBarry Smith 
36147c6ae99SBarry Smith   PetscFunctionBegin;
362171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36347c6ae99SBarry Smith   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
36447c6ae99SBarry Smith   PetscFunctionReturn(0);
36547c6ae99SBarry Smith }
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith #undef __FUNCT__
36847c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector"
36947c6ae99SBarry Smith /*@
370aa219208SBarry Smith     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
37147c6ae99SBarry Smith 
37247c6ae99SBarry Smith     Not Collective
37347c6ae99SBarry Smith 
37447c6ae99SBarry Smith     Input Parameter:
37547c6ae99SBarry Smith .   dm - the DM object
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith     Output Parameter:
37847c6ae99SBarry Smith .   vec - the local vector
37947c6ae99SBarry Smith 
380073dac72SJed Brown     Level: beginner
38147c6ae99SBarry Smith 
382e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
38347c6ae99SBarry Smith 
38447c6ae99SBarry Smith @*/
3857087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
38647c6ae99SBarry Smith {
38747c6ae99SBarry Smith   PetscErrorCode ierr;
38847c6ae99SBarry Smith 
38947c6ae99SBarry Smith   PetscFunctionBegin;
390171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
39147c6ae99SBarry Smith   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
39247c6ae99SBarry Smith   PetscFunctionReturn(0);
39347c6ae99SBarry Smith }
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith #undef __FUNCT__
3961411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping"
3971411c6eeSJed Brown /*@
3981411c6eeSJed Brown    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
3991411c6eeSJed Brown 
4001411c6eeSJed Brown    Collective on DM
4011411c6eeSJed Brown 
4021411c6eeSJed Brown    Input Parameter:
4031411c6eeSJed Brown .  dm - the DM that provides the mapping
4041411c6eeSJed Brown 
4051411c6eeSJed Brown    Output Parameter:
4061411c6eeSJed Brown .  ltog - the mapping
4071411c6eeSJed Brown 
4081411c6eeSJed Brown    Level: intermediate
4091411c6eeSJed Brown 
4101411c6eeSJed Brown    Notes:
4111411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMapping() or
4121411c6eeSJed Brown    MatSetLocalToGlobalMapping().
4131411c6eeSJed Brown 
4141411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
4151411c6eeSJed Brown @*/
4167087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
4171411c6eeSJed Brown {
4181411c6eeSJed Brown   PetscErrorCode ierr;
4191411c6eeSJed Brown 
4201411c6eeSJed Brown   PetscFunctionBegin;
4211411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4221411c6eeSJed Brown   PetscValidPointer(ltog,2);
4231411c6eeSJed Brown   if (!dm->ltogmap) {
4241411c6eeSJed Brown     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
4251411c6eeSJed Brown     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
4261411c6eeSJed Brown   }
4271411c6eeSJed Brown   *ltog = dm->ltogmap;
4281411c6eeSJed Brown   PetscFunctionReturn(0);
4291411c6eeSJed Brown }
4301411c6eeSJed Brown 
4311411c6eeSJed Brown #undef __FUNCT__
4321411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
4331411c6eeSJed Brown /*@
4341411c6eeSJed Brown    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
4351411c6eeSJed Brown 
4361411c6eeSJed Brown    Collective on DM
4371411c6eeSJed Brown 
4381411c6eeSJed Brown    Input Parameter:
4391411c6eeSJed Brown .  da - the distributed array that provides the mapping
4401411c6eeSJed Brown 
4411411c6eeSJed Brown    Output Parameter:
4421411c6eeSJed Brown .  ltog - the block mapping
4431411c6eeSJed Brown 
4441411c6eeSJed Brown    Level: intermediate
4451411c6eeSJed Brown 
4461411c6eeSJed Brown    Notes:
4471411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
4481411c6eeSJed Brown    MatSetLocalToGlobalMappingBlock().
4491411c6eeSJed Brown 
4501411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
4511411c6eeSJed Brown @*/
4527087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
4531411c6eeSJed Brown {
4541411c6eeSJed Brown   PetscErrorCode ierr;
4551411c6eeSJed Brown 
4561411c6eeSJed Brown   PetscFunctionBegin;
4571411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4581411c6eeSJed Brown   PetscValidPointer(ltog,2);
4591411c6eeSJed Brown   if (!dm->ltogmapb) {
4601411c6eeSJed Brown     PetscInt bs;
4611411c6eeSJed Brown     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
4621411c6eeSJed Brown     if (bs > 1) {
4631411c6eeSJed Brown       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
4641411c6eeSJed Brown       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
4651411c6eeSJed Brown     } else {
4661411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
4671411c6eeSJed Brown       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
4681411c6eeSJed Brown     }
4691411c6eeSJed Brown   }
4701411c6eeSJed Brown   *ltog = dm->ltogmapb;
4711411c6eeSJed Brown   PetscFunctionReturn(0);
4721411c6eeSJed Brown }
4731411c6eeSJed Brown 
4741411c6eeSJed Brown #undef __FUNCT__
4751411c6eeSJed Brown #define __FUNCT__ "DMGetBlockSize"
4761411c6eeSJed Brown /*@
4771411c6eeSJed Brown    DMGetBlockSize - Gets the inherent block size associated with a DM
4781411c6eeSJed Brown 
4791411c6eeSJed Brown    Not Collective
4801411c6eeSJed Brown 
4811411c6eeSJed Brown    Input Parameter:
4821411c6eeSJed Brown .  dm - the DM with block structure
4831411c6eeSJed Brown 
4841411c6eeSJed Brown    Output Parameter:
4851411c6eeSJed Brown .  bs - the block size, 1 implies no exploitable block structure
4861411c6eeSJed Brown 
4871411c6eeSJed Brown    Level: intermediate
4881411c6eeSJed Brown 
4891411c6eeSJed Brown .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
4901411c6eeSJed Brown @*/
4917087cfbeSBarry Smith PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
4921411c6eeSJed Brown {
4931411c6eeSJed Brown   PetscFunctionBegin;
4941411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4951411c6eeSJed Brown   PetscValidPointer(bs,2);
4961411c6eeSJed 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");
4971411c6eeSJed Brown   *bs = dm->bs;
4981411c6eeSJed Brown   PetscFunctionReturn(0);
4991411c6eeSJed Brown }
5001411c6eeSJed Brown 
5011411c6eeSJed Brown #undef __FUNCT__
502e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation"
50347c6ae99SBarry Smith /*@
504e727c939SJed Brown     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith     Collective on DM
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith     Input Parameter:
50947c6ae99SBarry Smith +   dm1 - the DM object
51047c6ae99SBarry Smith -   dm2 - the second, finer DM object
51147c6ae99SBarry Smith 
51247c6ae99SBarry Smith     Output Parameter:
51347c6ae99SBarry Smith +  mat - the interpolation
51447c6ae99SBarry Smith -  vec - the scaling (optional)
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith     Level: developer
51747c6ae99SBarry Smith 
51885afcc9aSBarry 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
51985afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
520d52bd9f3SBarry Smith 
521d52bd9f3SBarry Smith         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
522d52bd9f3SBarry Smith         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
52385afcc9aSBarry Smith 
52485afcc9aSBarry Smith 
525e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith @*/
528e727c939SJed Brown PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
52947c6ae99SBarry Smith {
53047c6ae99SBarry Smith   PetscErrorCode ierr;
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith   PetscFunctionBegin;
533171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
534171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
53525296bd5SBarry Smith   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
53647c6ae99SBarry Smith   PetscFunctionReturn(0);
53747c6ae99SBarry Smith }
53847c6ae99SBarry Smith 
53947c6ae99SBarry Smith #undef __FUNCT__
540e727c939SJed Brown #define __FUNCT__ "DMCreateInjection"
54147c6ae99SBarry Smith /*@
542e727c939SJed Brown     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith     Collective on DM
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith     Input Parameter:
54747c6ae99SBarry Smith +   dm1 - the DM object
54847c6ae99SBarry Smith -   dm2 - the second, finer DM object
54947c6ae99SBarry Smith 
55047c6ae99SBarry Smith     Output Parameter:
55147c6ae99SBarry Smith .   ctx - the injection
55247c6ae99SBarry Smith 
55347c6ae99SBarry Smith     Level: developer
55447c6ae99SBarry Smith 
55585afcc9aSBarry 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
55685afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
55785afcc9aSBarry Smith 
558e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith @*/
561e727c939SJed Brown PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
56247c6ae99SBarry Smith {
56347c6ae99SBarry Smith   PetscErrorCode ierr;
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith   PetscFunctionBegin;
566171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
567171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
56847c6ae99SBarry Smith   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
56947c6ae99SBarry Smith   PetscFunctionReturn(0);
57047c6ae99SBarry Smith }
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith #undef __FUNCT__
573e727c939SJed Brown #define __FUNCT__ "DMCreateColoring"
574d1e2c406SBarry Smith /*@C
575e727c939SJed Brown     DMCreateColoring - Gets coloring for a DMDA or DMComposite
57647c6ae99SBarry Smith 
57747c6ae99SBarry Smith     Collective on DM
57847c6ae99SBarry Smith 
57947c6ae99SBarry Smith     Input Parameter:
58047c6ae99SBarry Smith +   dm - the DM object
58147c6ae99SBarry Smith .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
58247c6ae99SBarry Smith -   matype - either MATAIJ or MATBAIJ
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith     Output Parameter:
58547c6ae99SBarry Smith .   coloring - the coloring
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith     Level: developer
58847c6ae99SBarry Smith 
589e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
59047c6ae99SBarry Smith 
591aab9d709SJed Brown @*/
592e727c939SJed Brown PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
59347c6ae99SBarry Smith {
59447c6ae99SBarry Smith   PetscErrorCode ierr;
59547c6ae99SBarry Smith 
59647c6ae99SBarry Smith   PetscFunctionBegin;
597171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
59847c6ae99SBarry Smith   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
59947c6ae99SBarry Smith   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
60047c6ae99SBarry Smith   PetscFunctionReturn(0);
60147c6ae99SBarry Smith }
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith #undef __FUNCT__
604950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix"
60547c6ae99SBarry Smith /*@C
606950540a4SJed Brown     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith     Collective on DM
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith     Input Parameter:
61147c6ae99SBarry Smith +   dm - the DM object
61247c6ae99SBarry Smith -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
61394013140SBarry Smith             any type which inherits from one of these (such as MATAIJ)
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith     Output Parameter:
61647c6ae99SBarry Smith .   mat - the empty Jacobian
61747c6ae99SBarry Smith 
618073dac72SJed Brown     Level: beginner
61947c6ae99SBarry Smith 
62094013140SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
62194013140SBarry Smith        do not need to do it yourself.
62294013140SBarry Smith 
62394013140SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
624aa219208SBarry Smith        the nonzero pattern call DMDASetMatPreallocateOnly()
62594013140SBarry Smith 
62694013140SBarry 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
62794013140SBarry Smith        internally by PETSc.
62894013140SBarry Smith 
62994013140SBarry Smith        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
630aa219208SBarry Smith        the indices for the global numbering for DMDAs which is complicated.
63194013140SBarry Smith 
632e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
63347c6ae99SBarry Smith 
634aab9d709SJed Brown @*/
635950540a4SJed Brown PetscErrorCode  DMCreateMatrix(DM dm,const MatType mtype,Mat *mat)
63647c6ae99SBarry Smith {
63747c6ae99SBarry Smith   PetscErrorCode ierr;
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith   PetscFunctionBegin;
640171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
641235683edSBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
642235683edSBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
643235683edSBarry Smith #endif
644c7b7c8a4SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
645c7b7c8a4SJed Brown   PetscValidPointer(mat,3);
646073dac72SJed Brown   if (dm->mattype) {
64725296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
648073dac72SJed Brown   } else {
64925296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
650c7b7c8a4SJed Brown   }
65147c6ae99SBarry Smith   PetscFunctionReturn(0);
65247c6ae99SBarry Smith }
65347c6ae99SBarry Smith 
65447c6ae99SBarry Smith #undef __FUNCT__
655732e2eb9SMatthew G Knepley #define __FUNCT__ "DMSetMatrixPreallocateOnly"
656732e2eb9SMatthew G Knepley /*@
657950540a4SJed Brown   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
658732e2eb9SMatthew G Knepley     preallocated but the nonzero structure and zero values will not be set.
659732e2eb9SMatthew G Knepley 
660732e2eb9SMatthew G Knepley   Logically Collective on DMDA
661732e2eb9SMatthew G Knepley 
662732e2eb9SMatthew G Knepley   Input Parameter:
663732e2eb9SMatthew G Knepley + dm - the DM
664732e2eb9SMatthew G Knepley - only - PETSC_TRUE if only want preallocation
665732e2eb9SMatthew G Knepley 
666732e2eb9SMatthew G Knepley   Level: developer
667950540a4SJed Brown .seealso DMCreateMatrix()
668732e2eb9SMatthew G Knepley @*/
669732e2eb9SMatthew G Knepley PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
670732e2eb9SMatthew G Knepley {
671732e2eb9SMatthew G Knepley   PetscFunctionBegin;
672732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
673732e2eb9SMatthew G Knepley   dm->prealloc_only = only;
674732e2eb9SMatthew G Knepley   PetscFunctionReturn(0);
675732e2eb9SMatthew G Knepley }
676732e2eb9SMatthew G Knepley 
677732e2eb9SMatthew G Knepley #undef __FUNCT__
678a89ea682SMatthew G Knepley #define __FUNCT__ "DMGetWorkArray"
679a89ea682SMatthew G Knepley /*@C
680a89ea682SMatthew G Knepley   DMGetWorkArray - Gets a work array guaranteed to be at least the input size
681a89ea682SMatthew G Knepley 
682a89ea682SMatthew G Knepley   Not Collective
683a89ea682SMatthew G Knepley 
684a89ea682SMatthew G Knepley   Input Parameters:
685a89ea682SMatthew G Knepley + dm - the DM object
686a89ea682SMatthew G Knepley - size - The minium size
687a89ea682SMatthew G Knepley 
688a89ea682SMatthew G Knepley   Output Parameter:
689a89ea682SMatthew G Knepley . array - the work array
690a89ea682SMatthew G Knepley 
691a89ea682SMatthew G Knepley   Level: developer
692a89ea682SMatthew G Knepley 
693a89ea682SMatthew G Knepley .seealso DMDestroy(), DMCreate()
694a89ea682SMatthew G Knepley @*/
695a89ea682SMatthew G Knepley PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array)
696a89ea682SMatthew G Knepley {
697a89ea682SMatthew G Knepley   PetscErrorCode ierr;
698a89ea682SMatthew G Knepley 
699a89ea682SMatthew G Knepley   PetscFunctionBegin;
700a89ea682SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
701a89ea682SMatthew G Knepley   PetscValidPointer(array,3);
702a89ea682SMatthew G Knepley   if (size > dm->workSize) {
703a89ea682SMatthew G Knepley     dm->workSize = size;
704a89ea682SMatthew G Knepley     ierr = PetscFree(dm->workArray);CHKERRQ(ierr);
705a89ea682SMatthew G Knepley     ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr);
706a89ea682SMatthew G Knepley   }
707a89ea682SMatthew G Knepley   *array = dm->workArray;
708a89ea682SMatthew G Knepley   PetscFunctionReturn(0);
709a89ea682SMatthew G Knepley }
710a89ea682SMatthew G Knepley 
7115fe1f584SPeter Brune 
712a89ea682SMatthew G Knepley #undef __FUNCT__
71347c6ae99SBarry Smith #define __FUNCT__ "DMRefine"
71447c6ae99SBarry Smith /*@
71547c6ae99SBarry Smith     DMRefine - Refines a DM object
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith     Collective on DM
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith     Input Parameter:
72047c6ae99SBarry Smith +   dm - the DM object
72147c6ae99SBarry Smith -   comm - the communicator to contain the new DM object (or PETSC_NULL)
72247c6ae99SBarry Smith 
72347c6ae99SBarry Smith     Output Parameter:
72447c6ae99SBarry Smith .   dmf - the refined DM
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith     Level: developer
72747c6ae99SBarry Smith 
728e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
72947c6ae99SBarry Smith 
73047c6ae99SBarry Smith @*/
7317087cfbeSBarry Smith PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
73247c6ae99SBarry Smith {
73347c6ae99SBarry Smith   PetscErrorCode ierr;
73447c6ae99SBarry Smith 
73547c6ae99SBarry Smith   PetscFunctionBegin;
736732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
73747c6ae99SBarry Smith   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
7384057135bSMatthew G Knepley   if (*dmf) {
739644e2e5bSBarry Smith     (*dmf)->ops->initialguess = dm->ops->initialguess;
740644e2e5bSBarry Smith     (*dmf)->ops->function     = dm->ops->function;
741644e2e5bSBarry Smith     (*dmf)->ops->functionj    = dm->ops->functionj;
742644e2e5bSBarry Smith     if (dm->ops->jacobian != DMComputeJacobianDefault) {
743644e2e5bSBarry Smith       (*dmf)->ops->jacobian     = dm->ops->jacobian;
744644e2e5bSBarry Smith     }
7458cd211a4SJed Brown     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
746644e2e5bSBarry Smith     (*dmf)->ctx     = dm->ctx;
747656b349aSBarry Smith     (*dmf)->levelup = dm->levelup + 1;
7484057135bSMatthew G Knepley   }
74947c6ae99SBarry Smith   PetscFunctionReturn(0);
75047c6ae99SBarry Smith }
75147c6ae99SBarry Smith 
75247c6ae99SBarry Smith #undef __FUNCT__
753eb3f98d2SBarry Smith #define __FUNCT__ "DMGetRefineLevel"
754eb3f98d2SBarry Smith /*@
755eb3f98d2SBarry Smith     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
756eb3f98d2SBarry Smith 
757eb3f98d2SBarry Smith     Not Collective
758eb3f98d2SBarry Smith 
759eb3f98d2SBarry Smith     Input Parameter:
760eb3f98d2SBarry Smith .   dm - the DM object
761eb3f98d2SBarry Smith 
762eb3f98d2SBarry Smith     Output Parameter:
763eb3f98d2SBarry Smith .   level - number of refinements
764eb3f98d2SBarry Smith 
765eb3f98d2SBarry Smith     Level: developer
766eb3f98d2SBarry Smith 
7676a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
768eb3f98d2SBarry Smith 
769eb3f98d2SBarry Smith @*/
770eb3f98d2SBarry Smith PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
771eb3f98d2SBarry Smith {
772eb3f98d2SBarry Smith   PetscFunctionBegin;
773eb3f98d2SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
774eb3f98d2SBarry Smith   *level = dm->levelup;
775eb3f98d2SBarry Smith   PetscFunctionReturn(0);
776eb3f98d2SBarry Smith }
777eb3f98d2SBarry Smith 
778eb3f98d2SBarry Smith #undef __FUNCT__
77947c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin"
78047c6ae99SBarry Smith /*@
78147c6ae99SBarry Smith     DMGlobalToLocalBegin - Begins updating local vectors from global vector
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith     Neighbor-wise Collective on DM
78447c6ae99SBarry Smith 
78547c6ae99SBarry Smith     Input Parameters:
78647c6ae99SBarry Smith +   dm - the DM object
78747c6ae99SBarry Smith .   g - the global vector
78847c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
78947c6ae99SBarry Smith -   l - the local vector
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith     Level: beginner
79347c6ae99SBarry Smith 
794e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
79547c6ae99SBarry Smith 
79647c6ae99SBarry Smith @*/
7977087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
79847c6ae99SBarry Smith {
79947c6ae99SBarry Smith   PetscErrorCode ierr;
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith   PetscFunctionBegin;
802171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
803843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
80447c6ae99SBarry Smith   PetscFunctionReturn(0);
80547c6ae99SBarry Smith }
80647c6ae99SBarry Smith 
80747c6ae99SBarry Smith #undef __FUNCT__
80847c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd"
80947c6ae99SBarry Smith /*@
81047c6ae99SBarry Smith     DMGlobalToLocalEnd - Ends updating local vectors from global vector
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith     Neighbor-wise Collective on DM
81347c6ae99SBarry Smith 
81447c6ae99SBarry Smith     Input Parameters:
81547c6ae99SBarry Smith +   dm - the DM object
81647c6ae99SBarry Smith .   g - the global vector
81747c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
81847c6ae99SBarry Smith -   l - the local vector
81947c6ae99SBarry Smith 
82047c6ae99SBarry Smith 
82147c6ae99SBarry Smith     Level: beginner
82247c6ae99SBarry Smith 
823e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith @*/
8267087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
82747c6ae99SBarry Smith {
82847c6ae99SBarry Smith   PetscErrorCode ierr;
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith   PetscFunctionBegin;
831171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
832843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
83347c6ae99SBarry Smith   PetscFunctionReturn(0);
83447c6ae99SBarry Smith }
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith #undef __FUNCT__
8379a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalBegin"
83847c6ae99SBarry Smith /*@
8399a42bb27SBarry Smith     DMLocalToGlobalBegin - updates global vectors from local vectors
8409a42bb27SBarry Smith 
8419a42bb27SBarry Smith     Neighbor-wise Collective on DM
8429a42bb27SBarry Smith 
8439a42bb27SBarry Smith     Input Parameters:
8449a42bb27SBarry Smith +   dm - the DM object
845f6813fd5SJed Brown .   l - the local vector
8469a42bb27SBarry 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
8479a42bb27SBarry Smith            base point.
848f6813fd5SJed Brown - - the global vector
8499a42bb27SBarry Smith 
8509a42bb27SBarry 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
8519a42bb27SBarry 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
8529a42bb27SBarry Smith            global array to the final global array with VecAXPY().
8539a42bb27SBarry Smith 
8549a42bb27SBarry Smith     Level: beginner
8559a42bb27SBarry Smith 
856e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
8579a42bb27SBarry Smith 
8589a42bb27SBarry Smith @*/
8597087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
8609a42bb27SBarry Smith {
8619a42bb27SBarry Smith   PetscErrorCode ierr;
8629a42bb27SBarry Smith 
8639a42bb27SBarry Smith   PetscFunctionBegin;
864171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
865843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
8669a42bb27SBarry Smith   PetscFunctionReturn(0);
8679a42bb27SBarry Smith }
8689a42bb27SBarry Smith 
8699a42bb27SBarry Smith #undef __FUNCT__
8709a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalEnd"
8719a42bb27SBarry Smith /*@
8729a42bb27SBarry Smith     DMLocalToGlobalEnd - updates global vectors from local vectors
87347c6ae99SBarry Smith 
87447c6ae99SBarry Smith     Neighbor-wise Collective on DM
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith     Input Parameters:
87747c6ae99SBarry Smith +   dm - the DM object
878f6813fd5SJed Brown .   l - the local vector
87947c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
880f6813fd5SJed Brown -   g - the global vector
88147c6ae99SBarry Smith 
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith     Level: beginner
88447c6ae99SBarry Smith 
885e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
88647c6ae99SBarry Smith 
88747c6ae99SBarry Smith @*/
8887087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
88947c6ae99SBarry Smith {
89047c6ae99SBarry Smith   PetscErrorCode ierr;
89147c6ae99SBarry Smith 
89247c6ae99SBarry Smith   PetscFunctionBegin;
893171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
894843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
89547c6ae99SBarry Smith   PetscFunctionReturn(0);
89647c6ae99SBarry Smith }
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith #undef __FUNCT__
89947c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobianDefault"
90047c6ae99SBarry Smith /*@
90147c6ae99SBarry Smith     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
90247c6ae99SBarry Smith 
90347c6ae99SBarry Smith     Collective on DM
90447c6ae99SBarry Smith 
90547c6ae99SBarry Smith     Input Parameter:
90647c6ae99SBarry Smith +   dm - the DM object
90747c6ae99SBarry Smith .   x - location to compute Jacobian at; may be ignored for linear problems
90847c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
90947c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
91047c6ae99SBarry Smith 
91147c6ae99SBarry Smith     Level: developer
91247c6ae99SBarry Smith 
913e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
91447c6ae99SBarry Smith          DMSetFunction()
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith @*/
9177087cfbeSBarry Smith PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
91847c6ae99SBarry Smith {
91947c6ae99SBarry Smith   PetscErrorCode ierr;
920171400e9SBarry Smith 
92147c6ae99SBarry Smith   PetscFunctionBegin;
922171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
92347c6ae99SBarry Smith   *stflag = SAME_NONZERO_PATTERN;
92447c6ae99SBarry Smith   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
92547c6ae99SBarry Smith   if (A != B) {
92647c6ae99SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92747c6ae99SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92847c6ae99SBarry Smith   }
92947c6ae99SBarry Smith   PetscFunctionReturn(0);
93047c6ae99SBarry Smith }
93147c6ae99SBarry Smith 
93247c6ae99SBarry Smith #undef __FUNCT__
93347c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen"
93447c6ae99SBarry Smith /*@
93547c6ae99SBarry Smith     DMCoarsen - Coarsens a DM object
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith     Collective on DM
93847c6ae99SBarry Smith 
93947c6ae99SBarry Smith     Input Parameter:
94047c6ae99SBarry Smith +   dm - the DM object
94147c6ae99SBarry Smith -   comm - the communicator to contain the new DM object (or PETSC_NULL)
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith     Output Parameter:
94447c6ae99SBarry Smith .   dmc - the coarsened DM
94547c6ae99SBarry Smith 
94647c6ae99SBarry Smith     Level: developer
94747c6ae99SBarry Smith 
948e727c939SJed Brown .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith @*/
9517087cfbeSBarry Smith PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
95247c6ae99SBarry Smith {
95347c6ae99SBarry Smith   PetscErrorCode ierr;
954b17ce1afSJed Brown   DMCoarsenHookLink link;
95547c6ae99SBarry Smith 
95647c6ae99SBarry Smith   PetscFunctionBegin;
957171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
95847c6ae99SBarry Smith   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
95947c6ae99SBarry Smith   (*dmc)->ops->initialguess = dm->ops->initialguess;
96047c6ae99SBarry Smith   (*dmc)->ops->function     = dm->ops->function;
96147c6ae99SBarry Smith   (*dmc)->ops->functionj    = dm->ops->functionj;
96247c6ae99SBarry Smith   if (dm->ops->jacobian != DMComputeJacobianDefault) {
96347c6ae99SBarry Smith     (*dmc)->ops->jacobian     = dm->ops->jacobian;
96447c6ae99SBarry Smith   }
9658cd211a4SJed Brown   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
966644e2e5bSBarry Smith   (*dmc)->ctx       = dm->ctx;
967656b349aSBarry Smith   (*dmc)->leveldown = dm->leveldown + 1;
968b17ce1afSJed Brown   for (link=dm->coarsenhook; link; link=link->next) {
969b17ce1afSJed Brown     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
970b17ce1afSJed Brown   }
971b17ce1afSJed Brown   PetscFunctionReturn(0);
972b17ce1afSJed Brown }
973b17ce1afSJed Brown 
974b17ce1afSJed Brown #undef __FUNCT__
975b17ce1afSJed Brown #define __FUNCT__ "DMCoarsenHookAdd"
976b17ce1afSJed Brown /*@
977b17ce1afSJed Brown    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
978b17ce1afSJed Brown 
979b17ce1afSJed Brown    Logically Collective
980b17ce1afSJed Brown 
981b17ce1afSJed Brown    Input Arguments:
982b17ce1afSJed Brown +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
983b17ce1afSJed Brown .  coarsenhook - function to run when setting up a coarser level
984b17ce1afSJed Brown .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
985b17ce1afSJed Brown -  ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL)
986b17ce1afSJed Brown 
987b17ce1afSJed Brown    Calling sequence of coarsenhook:
988b17ce1afSJed Brown $    coarsenhook(DM fine,DM coarse,void *ctx);
989b17ce1afSJed Brown 
990b17ce1afSJed Brown +  fine - fine level DM
991b17ce1afSJed Brown .  coarse - coarse level DM to restrict problem to
992b17ce1afSJed Brown -  ctx - optional user-defined function context
993b17ce1afSJed Brown 
994b17ce1afSJed Brown    Calling sequence for restricthook:
995b17ce1afSJed Brown $    restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx)
996b17ce1afSJed Brown 
997b17ce1afSJed Brown +  fine - fine level DM
998b17ce1afSJed Brown .  mrestrict - matrix restricting a fine-level solution to the coarse grid
999b17ce1afSJed Brown .  inject - matrix restricting by applying the transpose of injection
1000b17ce1afSJed Brown .  coarse - coarse level DM to update
1001b17ce1afSJed Brown -  ctx - optional user-defined function context
1002b17ce1afSJed Brown 
1003b17ce1afSJed Brown    Level: advanced
1004b17ce1afSJed Brown 
1005b17ce1afSJed Brown    Notes:
1006b17ce1afSJed Brown    This function is only needed if auxiliary data needs to be set up on coarse grids.
1007b17ce1afSJed Brown 
1008b17ce1afSJed Brown    If this function is called multiple times, the hooks will be run in the order they are added.
1009b17ce1afSJed Brown 
1010b17ce1afSJed Brown    The
1011b17ce1afSJed Brown 
1012b17ce1afSJed Brown    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1013b17ce1afSJed Brown    extract the finest level information from its context (instead of from the SNES).
1014b17ce1afSJed Brown 
1015b17ce1afSJed Brown .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1016b17ce1afSJed Brown @*/
1017b17ce1afSJed Brown PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1018b17ce1afSJed Brown {
1019b17ce1afSJed Brown   PetscErrorCode ierr;
1020b17ce1afSJed Brown   DMCoarsenHookLink link,*p;
1021b17ce1afSJed Brown 
1022b17ce1afSJed Brown   PetscFunctionBegin;
1023b17ce1afSJed Brown   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
1024b17ce1afSJed Brown   for (p=&fine->coarsenhook; *p; *p=(*p)->next) {} /* Scan to the end of the current list of hooks */
1025b17ce1afSJed Brown   ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1026b17ce1afSJed Brown   link->coarsenhook = coarsenhook;
1027b17ce1afSJed Brown   link->restricthook = restricthook;
1028b17ce1afSJed Brown   link->ctx = ctx;
1029b17ce1afSJed Brown   *p = link;
1030b17ce1afSJed Brown   PetscFunctionReturn(0);
1031b17ce1afSJed Brown }
1032b17ce1afSJed Brown 
1033b17ce1afSJed Brown #undef __FUNCT__
1034b17ce1afSJed Brown #define __FUNCT__ "DMRestrict"
1035b17ce1afSJed Brown /*@
1036b17ce1afSJed Brown    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1037b17ce1afSJed Brown 
1038b17ce1afSJed Brown    Collective if any hooks are
1039b17ce1afSJed Brown 
1040b17ce1afSJed Brown    Input Arguments:
1041b17ce1afSJed Brown +  fine - finer DM to use as a base
1042b17ce1afSJed Brown .  restrct - restriction matrix, apply using MatRestrict()
1043b17ce1afSJed Brown .  inject - injection matrix, also use MatRestrict()
1044b17ce1afSJed Brown -  coarse - coarer DM to update
1045b17ce1afSJed Brown 
1046b17ce1afSJed Brown    Level: developer
1047b17ce1afSJed Brown 
1048b17ce1afSJed Brown .seealso: DMCoarsenHookAdd(), MatRestrict()
1049b17ce1afSJed Brown @*/
1050b17ce1afSJed Brown PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1051b17ce1afSJed Brown {
1052b17ce1afSJed Brown   PetscErrorCode ierr;
1053b17ce1afSJed Brown   DMCoarsenHookLink link;
1054b17ce1afSJed Brown 
1055b17ce1afSJed Brown   PetscFunctionBegin;
1056b17ce1afSJed Brown   for (link=fine->coarsenhook; link; link=link->next) {
1057b17ce1afSJed Brown     if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);}
1058b17ce1afSJed Brown   }
105947c6ae99SBarry Smith   PetscFunctionReturn(0);
106047c6ae99SBarry Smith }
106147c6ae99SBarry Smith 
106247c6ae99SBarry Smith #undef __FUNCT__
10635fe1f584SPeter Brune #define __FUNCT__ "DMGetCoarsenLevel"
10645fe1f584SPeter Brune /*@
10656a7d9d85SPeter Brune     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
10665fe1f584SPeter Brune 
10675fe1f584SPeter Brune     Not Collective
10685fe1f584SPeter Brune 
10695fe1f584SPeter Brune     Input Parameter:
10705fe1f584SPeter Brune .   dm - the DM object
10715fe1f584SPeter Brune 
10725fe1f584SPeter Brune     Output Parameter:
10736a7d9d85SPeter Brune .   level - number of coarsenings
10745fe1f584SPeter Brune 
10755fe1f584SPeter Brune     Level: developer
10765fe1f584SPeter Brune 
10776a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
10785fe1f584SPeter Brune 
10795fe1f584SPeter Brune @*/
10805fe1f584SPeter Brune PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
10815fe1f584SPeter Brune {
10825fe1f584SPeter Brune   PetscFunctionBegin;
10835fe1f584SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
10845fe1f584SPeter Brune   *level = dm->leveldown;
10855fe1f584SPeter Brune   PetscFunctionReturn(0);
10865fe1f584SPeter Brune }
10875fe1f584SPeter Brune 
10885fe1f584SPeter Brune 
10895fe1f584SPeter Brune 
10905fe1f584SPeter Brune #undef __FUNCT__
109147c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy"
109247c6ae99SBarry Smith /*@C
109347c6ae99SBarry Smith     DMRefineHierarchy - Refines a DM object, all levels at once
109447c6ae99SBarry Smith 
109547c6ae99SBarry Smith     Collective on DM
109647c6ae99SBarry Smith 
109747c6ae99SBarry Smith     Input Parameter:
109847c6ae99SBarry Smith +   dm - the DM object
109947c6ae99SBarry Smith -   nlevels - the number of levels of refinement
110047c6ae99SBarry Smith 
110147c6ae99SBarry Smith     Output Parameter:
110247c6ae99SBarry Smith .   dmf - the refined DM hierarchy
110347c6ae99SBarry Smith 
110447c6ae99SBarry Smith     Level: developer
110547c6ae99SBarry Smith 
1106e727c939SJed Brown .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
110747c6ae99SBarry Smith 
110847c6ae99SBarry Smith @*/
11097087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
111047c6ae99SBarry Smith {
111147c6ae99SBarry Smith   PetscErrorCode ierr;
111247c6ae99SBarry Smith 
111347c6ae99SBarry Smith   PetscFunctionBegin;
1114171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
111547c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
111647c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
111747c6ae99SBarry Smith   if (dm->ops->refinehierarchy) {
111847c6ae99SBarry Smith     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
111947c6ae99SBarry Smith   } else if (dm->ops->refine) {
112047c6ae99SBarry Smith     PetscInt i;
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
112347c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
112447c6ae99SBarry Smith       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
112547c6ae99SBarry Smith     }
112647c6ae99SBarry Smith   } else {
112747c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
112847c6ae99SBarry Smith   }
112947c6ae99SBarry Smith   PetscFunctionReturn(0);
113047c6ae99SBarry Smith }
113147c6ae99SBarry Smith 
113247c6ae99SBarry Smith #undef __FUNCT__
113347c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy"
113447c6ae99SBarry Smith /*@C
113547c6ae99SBarry Smith     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
113647c6ae99SBarry Smith 
113747c6ae99SBarry Smith     Collective on DM
113847c6ae99SBarry Smith 
113947c6ae99SBarry Smith     Input Parameter:
114047c6ae99SBarry Smith +   dm - the DM object
114147c6ae99SBarry Smith -   nlevels - the number of levels of coarsening
114247c6ae99SBarry Smith 
114347c6ae99SBarry Smith     Output Parameter:
114447c6ae99SBarry Smith .   dmc - the coarsened DM hierarchy
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith     Level: developer
114747c6ae99SBarry Smith 
1148e727c939SJed Brown .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
114947c6ae99SBarry Smith 
115047c6ae99SBarry Smith @*/
11517087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
115247c6ae99SBarry Smith {
115347c6ae99SBarry Smith   PetscErrorCode ierr;
115447c6ae99SBarry Smith 
115547c6ae99SBarry Smith   PetscFunctionBegin;
1156171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
115747c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
115847c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
115947c6ae99SBarry Smith   PetscValidPointer(dmc,3);
116047c6ae99SBarry Smith   if (dm->ops->coarsenhierarchy) {
116147c6ae99SBarry Smith     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
116247c6ae99SBarry Smith   } else if (dm->ops->coarsen) {
116347c6ae99SBarry Smith     PetscInt i;
116447c6ae99SBarry Smith 
116547c6ae99SBarry Smith     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
116647c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
116747c6ae99SBarry Smith       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
116847c6ae99SBarry Smith     }
116947c6ae99SBarry Smith   } else {
117047c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
117147c6ae99SBarry Smith   }
117247c6ae99SBarry Smith   PetscFunctionReturn(0);
117347c6ae99SBarry Smith }
117447c6ae99SBarry Smith 
117547c6ae99SBarry Smith #undef __FUNCT__
1176e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates"
117747c6ae99SBarry Smith /*@
1178e727c939SJed Brown    DMCreateAggregates - Gets the aggregates that map between
117947c6ae99SBarry Smith    grids associated with two DMs.
118047c6ae99SBarry Smith 
118147c6ae99SBarry Smith    Collective on DM
118247c6ae99SBarry Smith 
118347c6ae99SBarry Smith    Input Parameters:
118447c6ae99SBarry Smith +  dmc - the coarse grid DM
118547c6ae99SBarry Smith -  dmf - the fine grid DM
118647c6ae99SBarry Smith 
118747c6ae99SBarry Smith    Output Parameters:
118847c6ae99SBarry Smith .  rest - the restriction matrix (transpose of the projection matrix)
118947c6ae99SBarry Smith 
119047c6ae99SBarry Smith    Level: intermediate
119147c6ae99SBarry Smith 
119247c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid
119347c6ae99SBarry Smith 
1194e727c939SJed Brown .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
119547c6ae99SBarry Smith @*/
1196e727c939SJed Brown PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
119747c6ae99SBarry Smith {
119847c6ae99SBarry Smith   PetscErrorCode ierr;
119947c6ae99SBarry Smith 
120047c6ae99SBarry Smith   PetscFunctionBegin;
1201171400e9SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1202171400e9SBarry Smith   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
120347c6ae99SBarry Smith   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
120447c6ae99SBarry Smith   PetscFunctionReturn(0);
120547c6ae99SBarry Smith }
120647c6ae99SBarry Smith 
120747c6ae99SBarry Smith #undef __FUNCT__
12081a266240SBarry Smith #define __FUNCT__ "DMSetApplicationContextDestroy"
12091a266240SBarry Smith /*@C
12101a266240SBarry Smith     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
12111a266240SBarry Smith 
12121a266240SBarry Smith     Not Collective
12131a266240SBarry Smith 
12141a266240SBarry Smith     Input Parameters:
12151a266240SBarry Smith +   dm - the DM object
12161a266240SBarry Smith -   destroy - the destroy function
12171a266240SBarry Smith 
12181a266240SBarry Smith     Level: intermediate
12191a266240SBarry Smith 
1220e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
12211a266240SBarry Smith 
1222f07f9ceaSJed Brown @*/
12231a266240SBarry Smith PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
12241a266240SBarry Smith {
12251a266240SBarry Smith   PetscFunctionBegin;
1226171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12271a266240SBarry Smith   dm->ctxdestroy = destroy;
12281a266240SBarry Smith   PetscFunctionReturn(0);
12291a266240SBarry Smith }
12301a266240SBarry Smith 
12311a266240SBarry Smith #undef __FUNCT__
12321b2093e4SBarry Smith #define __FUNCT__ "DMSetApplicationContext"
1233b07ff414SBarry Smith /*@
12341b2093e4SBarry Smith     DMSetApplicationContext - Set a user context into a DM object
123547c6ae99SBarry Smith 
123647c6ae99SBarry Smith     Not Collective
123747c6ae99SBarry Smith 
123847c6ae99SBarry Smith     Input Parameters:
123947c6ae99SBarry Smith +   dm - the DM object
124047c6ae99SBarry Smith -   ctx - the user context
124147c6ae99SBarry Smith 
124247c6ae99SBarry Smith     Level: intermediate
124347c6ae99SBarry Smith 
1244e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
124547c6ae99SBarry Smith 
124647c6ae99SBarry Smith @*/
12471b2093e4SBarry Smith PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
124847c6ae99SBarry Smith {
124947c6ae99SBarry Smith   PetscFunctionBegin;
1250171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
125147c6ae99SBarry Smith   dm->ctx = ctx;
125247c6ae99SBarry Smith   PetscFunctionReturn(0);
125347c6ae99SBarry Smith }
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith #undef __FUNCT__
12561b2093e4SBarry Smith #define __FUNCT__ "DMGetApplicationContext"
125747c6ae99SBarry Smith /*@
12581b2093e4SBarry Smith     DMGetApplicationContext - Gets a user context from a DM object
125947c6ae99SBarry Smith 
126047c6ae99SBarry Smith     Not Collective
126147c6ae99SBarry Smith 
126247c6ae99SBarry Smith     Input Parameter:
126347c6ae99SBarry Smith .   dm - the DM object
126447c6ae99SBarry Smith 
126547c6ae99SBarry Smith     Output Parameter:
126647c6ae99SBarry Smith .   ctx - the user context
126747c6ae99SBarry Smith 
126847c6ae99SBarry Smith     Level: intermediate
126947c6ae99SBarry Smith 
1270e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
127147c6ae99SBarry Smith 
127247c6ae99SBarry Smith @*/
12731b2093e4SBarry Smith PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
127447c6ae99SBarry Smith {
127547c6ae99SBarry Smith   PetscFunctionBegin;
1276171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12771b2093e4SBarry Smith   *(void**)ctx = dm->ctx;
127847c6ae99SBarry Smith   PetscFunctionReturn(0);
127947c6ae99SBarry Smith }
128047c6ae99SBarry Smith 
128147c6ae99SBarry Smith #undef __FUNCT__
128247c6ae99SBarry Smith #define __FUNCT__ "DMSetInitialGuess"
12837e833e3aSBarry Smith /*@C
128447c6ae99SBarry Smith     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
128547c6ae99SBarry Smith 
128647c6ae99SBarry Smith     Logically Collective on DM
128747c6ae99SBarry Smith 
128847c6ae99SBarry Smith     Input Parameter:
128947c6ae99SBarry Smith +   dm - the DM object to destroy
129047c6ae99SBarry Smith -   f - the function to compute the initial guess
129147c6ae99SBarry Smith 
129247c6ae99SBarry Smith     Level: intermediate
129347c6ae99SBarry Smith 
1294e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
129547c6ae99SBarry Smith 
1296f07f9ceaSJed Brown @*/
12977087cfbeSBarry Smith PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
129847c6ae99SBarry Smith {
129947c6ae99SBarry Smith   PetscFunctionBegin;
1300171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
130147c6ae99SBarry Smith   dm->ops->initialguess = f;
130247c6ae99SBarry Smith   PetscFunctionReturn(0);
130347c6ae99SBarry Smith }
130447c6ae99SBarry Smith 
130547c6ae99SBarry Smith #undef __FUNCT__
130647c6ae99SBarry Smith #define __FUNCT__ "DMSetFunction"
13077e833e3aSBarry Smith /*@C
130847c6ae99SBarry Smith     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
130947c6ae99SBarry Smith 
131047c6ae99SBarry Smith     Logically Collective on DM
131147c6ae99SBarry Smith 
131247c6ae99SBarry Smith     Input Parameter:
131347c6ae99SBarry Smith +   dm - the DM object
131447c6ae99SBarry Smith -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
131547c6ae99SBarry Smith 
131647c6ae99SBarry Smith     Level: intermediate
131747c6ae99SBarry Smith 
131847c6ae99SBarry Smith     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
131947c6ae99SBarry Smith            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
132047c6ae99SBarry Smith 
1321e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
132247c6ae99SBarry Smith          DMSetJacobian()
132347c6ae99SBarry Smith 
1324f07f9ceaSJed Brown @*/
13257087cfbeSBarry Smith PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
132647c6ae99SBarry Smith {
132747c6ae99SBarry Smith   PetscFunctionBegin;
1328171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
132947c6ae99SBarry Smith   dm->ops->function = f;
133047c6ae99SBarry Smith   if (f) {
133147c6ae99SBarry Smith     dm->ops->functionj = f;
133247c6ae99SBarry Smith   }
133347c6ae99SBarry Smith   PetscFunctionReturn(0);
133447c6ae99SBarry Smith }
133547c6ae99SBarry Smith 
133647c6ae99SBarry Smith #undef __FUNCT__
133747c6ae99SBarry Smith #define __FUNCT__ "DMSetJacobian"
13387e833e3aSBarry Smith /*@C
133947c6ae99SBarry Smith     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
134047c6ae99SBarry Smith 
134147c6ae99SBarry Smith     Logically Collective on DM
134247c6ae99SBarry Smith 
134347c6ae99SBarry Smith     Input Parameter:
134447c6ae99SBarry Smith +   dm - the DM object to destroy
134547c6ae99SBarry Smith -   f - the function to compute the matrix entries
134647c6ae99SBarry Smith 
134747c6ae99SBarry Smith     Level: intermediate
134847c6ae99SBarry Smith 
1349e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
135047c6ae99SBarry Smith          DMSetFunction()
135147c6ae99SBarry Smith 
1352f07f9ceaSJed Brown @*/
13537087cfbeSBarry Smith PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
135447c6ae99SBarry Smith {
135547c6ae99SBarry Smith   PetscFunctionBegin;
1356171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
135747c6ae99SBarry Smith   dm->ops->jacobian = f;
135847c6ae99SBarry Smith   PetscFunctionReturn(0);
135947c6ae99SBarry Smith }
136047c6ae99SBarry Smith 
136147c6ae99SBarry Smith #undef __FUNCT__
1362*08da532bSDmitry Karpeev #define __FUNCT__ "DMSetVariableBounds"
1363*08da532bSDmitry Karpeev /*@C
1364*08da532bSDmitry Karpeev     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
1365*08da532bSDmitry Karpeev 
1366*08da532bSDmitry Karpeev     Logically Collective on DM
1367*08da532bSDmitry Karpeev 
1368*08da532bSDmitry Karpeev     Input Parameter:
1369*08da532bSDmitry Karpeev +   dm - the DM object
1370*08da532bSDmitry Karpeev -   f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set)
1371*08da532bSDmitry Karpeev 
1372*08da532bSDmitry Karpeev     Level: intermediate
1373*08da532bSDmitry Karpeev 
1374*08da532bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1375*08da532bSDmitry Karpeev          DMSetJacobian()
1376*08da532bSDmitry Karpeev 
1377*08da532bSDmitry Karpeev @*/
1378*08da532bSDmitry Karpeev PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1379*08da532bSDmitry Karpeev {
1380*08da532bSDmitry Karpeev   PetscFunctionBegin;
1381*08da532bSDmitry Karpeev   dm->ops->computevariablebounds = f;
1382*08da532bSDmitry Karpeev   PetscFunctionReturn(0);
1383*08da532bSDmitry Karpeev }
1384*08da532bSDmitry Karpeev 
1385*08da532bSDmitry Karpeev #undef __FUNCT__
1386*08da532bSDmitry Karpeev #define __FUNCT__ "DMHasVariableBounds"
1387*08da532bSDmitry Karpeev /*@
1388*08da532bSDmitry Karpeev     DMHasVariableBounds - does the DM object have a variable bounds function?
1389*08da532bSDmitry Karpeev 
1390*08da532bSDmitry Karpeev     Not Collective
1391*08da532bSDmitry Karpeev 
1392*08da532bSDmitry Karpeev     Input Parameter:
1393*08da532bSDmitry Karpeev .   dm - the DM object to destroy
1394*08da532bSDmitry Karpeev 
1395*08da532bSDmitry Karpeev     Output Parameter:
1396*08da532bSDmitry Karpeev .   flg - PETSC_TRUE if the variable bounds function exists
1397*08da532bSDmitry Karpeev 
1398*08da532bSDmitry Karpeev     Level: developer
1399*08da532bSDmitry Karpeev 
1400*08da532bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1401*08da532bSDmitry Karpeev 
1402*08da532bSDmitry Karpeev @*/
1403*08da532bSDmitry Karpeev PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
1404*08da532bSDmitry Karpeev {
1405*08da532bSDmitry Karpeev   PetscFunctionBegin;
1406*08da532bSDmitry Karpeev   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
1407*08da532bSDmitry Karpeev   PetscFunctionReturn(0);
1408*08da532bSDmitry Karpeev }
1409*08da532bSDmitry Karpeev 
1410*08da532bSDmitry Karpeev #undef __FUNCT__
1411*08da532bSDmitry Karpeev #define __FUNCT__ "DMComputeVariableBounds"
1412*08da532bSDmitry Karpeev /*@C
1413*08da532bSDmitry Karpeev     DMComputeVariableBounds - compute variable bounds used by SNESVI.
1414*08da532bSDmitry Karpeev 
1415*08da532bSDmitry Karpeev     Logically Collective on DM
1416*08da532bSDmitry Karpeev 
1417*08da532bSDmitry Karpeev     Input Parameters:
1418*08da532bSDmitry Karpeev +   dm - the DM object to destroy
1419*08da532bSDmitry Karpeev -   x  - current solution at which the bounds are computed
1420*08da532bSDmitry Karpeev 
1421*08da532bSDmitry Karpeev     Output parameters:
1422*08da532bSDmitry Karpeev +   xl - lower bound
1423*08da532bSDmitry Karpeev -   xu - upper bound
1424*08da532bSDmitry Karpeev 
1425*08da532bSDmitry Karpeev     Level: intermediate
1426*08da532bSDmitry Karpeev 
1427*08da532bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1428*08da532bSDmitry Karpeev          DMSetFunction(), DMSetVariableBounds()
1429*08da532bSDmitry Karpeev 
1430*08da532bSDmitry Karpeev @*/
1431*08da532bSDmitry Karpeev PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
1432*08da532bSDmitry Karpeev {
1433*08da532bSDmitry Karpeev   PetscErrorCode ierr;
1434*08da532bSDmitry Karpeev   PetscFunctionBegin;
1435*08da532bSDmitry Karpeev   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1436*08da532bSDmitry Karpeev   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
1437*08da532bSDmitry Karpeev   if(dm->ops->computevariablebounds) {
1438*08da532bSDmitry Karpeev     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr);
1439*08da532bSDmitry Karpeev   }
1440*08da532bSDmitry Karpeev   else {
1441*08da532bSDmitry Karpeev     ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr);
1442*08da532bSDmitry Karpeev     ierr = VecSet(xu,SNES_VI_INF);  CHKERRQ(ierr);
1443*08da532bSDmitry Karpeev   }
1444*08da532bSDmitry Karpeev   PetscFunctionReturn(0);
1445*08da532bSDmitry Karpeev }
1446*08da532bSDmitry Karpeev 
1447*08da532bSDmitry Karpeev #undef __FUNCT__
144847c6ae99SBarry Smith #define __FUNCT__ "DMComputeInitialGuess"
144947c6ae99SBarry Smith /*@
145047c6ae99SBarry Smith     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
145147c6ae99SBarry Smith 
145247c6ae99SBarry Smith     Collective on DM
145347c6ae99SBarry Smith 
145447c6ae99SBarry Smith     Input Parameter:
145547c6ae99SBarry Smith +   dm - the DM object to destroy
145647c6ae99SBarry Smith -   x - the vector to hold the initial guess values
145747c6ae99SBarry Smith 
145847c6ae99SBarry Smith     Level: developer
145947c6ae99SBarry Smith 
1460e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
146147c6ae99SBarry Smith 
146247c6ae99SBarry Smith @*/
14637087cfbeSBarry Smith PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
146447c6ae99SBarry Smith {
146547c6ae99SBarry Smith   PetscErrorCode ierr;
146647c6ae99SBarry Smith   PetscFunctionBegin;
1467171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
146847c6ae99SBarry Smith   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
146947c6ae99SBarry Smith   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
147047c6ae99SBarry Smith   PetscFunctionReturn(0);
147147c6ae99SBarry Smith }
147247c6ae99SBarry Smith 
147347c6ae99SBarry Smith #undef __FUNCT__
147447c6ae99SBarry Smith #define __FUNCT__ "DMHasInitialGuess"
147547c6ae99SBarry Smith /*@
147647c6ae99SBarry Smith     DMHasInitialGuess - does the DM object have an initial guess function
147747c6ae99SBarry Smith 
147847c6ae99SBarry Smith     Not Collective
147947c6ae99SBarry Smith 
148047c6ae99SBarry Smith     Input Parameter:
148147c6ae99SBarry Smith .   dm - the DM object to destroy
148247c6ae99SBarry Smith 
148347c6ae99SBarry Smith     Output Parameter:
148447c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
148547c6ae99SBarry Smith 
148647c6ae99SBarry Smith     Level: developer
148747c6ae99SBarry Smith 
1488e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
148947c6ae99SBarry Smith 
149047c6ae99SBarry Smith @*/
14917087cfbeSBarry Smith PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
149247c6ae99SBarry Smith {
149347c6ae99SBarry Smith   PetscFunctionBegin;
1494171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
149547c6ae99SBarry Smith   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
149647c6ae99SBarry Smith   PetscFunctionReturn(0);
149747c6ae99SBarry Smith }
149847c6ae99SBarry Smith 
149947c6ae99SBarry Smith #undef __FUNCT__
150047c6ae99SBarry Smith #define __FUNCT__ "DMHasFunction"
150147c6ae99SBarry Smith /*@
150247c6ae99SBarry Smith     DMHasFunction - does the DM object have a function
150347c6ae99SBarry Smith 
150447c6ae99SBarry Smith     Not Collective
150547c6ae99SBarry Smith 
150647c6ae99SBarry Smith     Input Parameter:
150747c6ae99SBarry Smith .   dm - the DM object to destroy
150847c6ae99SBarry Smith 
150947c6ae99SBarry Smith     Output Parameter:
151047c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
151147c6ae99SBarry Smith 
151247c6ae99SBarry Smith     Level: developer
151347c6ae99SBarry Smith 
1514e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
151547c6ae99SBarry Smith 
151647c6ae99SBarry Smith @*/
15177087cfbeSBarry Smith PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
151847c6ae99SBarry Smith {
151947c6ae99SBarry Smith   PetscFunctionBegin;
1520171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
152147c6ae99SBarry Smith   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
152247c6ae99SBarry Smith   PetscFunctionReturn(0);
152347c6ae99SBarry Smith }
152447c6ae99SBarry Smith 
152547c6ae99SBarry Smith #undef __FUNCT__
152647c6ae99SBarry Smith #define __FUNCT__ "DMHasJacobian"
152747c6ae99SBarry Smith /*@
152847c6ae99SBarry Smith     DMHasJacobian - does the DM object have a matrix function
152947c6ae99SBarry Smith 
153047c6ae99SBarry Smith     Not Collective
153147c6ae99SBarry Smith 
153247c6ae99SBarry Smith     Input Parameter:
153347c6ae99SBarry Smith .   dm - the DM object to destroy
153447c6ae99SBarry Smith 
153547c6ae99SBarry Smith     Output Parameter:
153647c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
153747c6ae99SBarry Smith 
153847c6ae99SBarry Smith     Level: developer
153947c6ae99SBarry Smith 
1540e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
154147c6ae99SBarry Smith 
154247c6ae99SBarry Smith @*/
15437087cfbeSBarry Smith PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
154447c6ae99SBarry Smith {
154547c6ae99SBarry Smith   PetscFunctionBegin;
1546171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
154747c6ae99SBarry Smith   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
154847c6ae99SBarry Smith   PetscFunctionReturn(0);
154947c6ae99SBarry Smith }
155047c6ae99SBarry Smith 
155147c6ae99SBarry Smith #undef  __FUNCT__
1552*08da532bSDmitry Karpeev #define __FUNCT__ "DMSetVec"
1553*08da532bSDmitry Karpeev /*@
1554*08da532bSDmitry Karpeev     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
1555*08da532bSDmitry Karpeev 
1556*08da532bSDmitry Karpeev     Collective on DM
1557*08da532bSDmitry Karpeev 
1558*08da532bSDmitry Karpeev     Input Parameter:
1559*08da532bSDmitry Karpeev +   dm - the DM object
1560*08da532bSDmitry Karpeev -   x - location to compute residual, Jacobian and VI bounds at; will be PETSC_NULL for linear problems.
1561*08da532bSDmitry Karpeev 
1562*08da532bSDmitry Karpeev     Level: developer
1563*08da532bSDmitry Karpeev 
1564*08da532bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1565*08da532bSDmitry Karpeev          DMSetFunction(), DMSetJacobian(), DMSetVariableBounds()
1566*08da532bSDmitry Karpeev 
1567*08da532bSDmitry Karpeev @*/
1568*08da532bSDmitry Karpeev PetscErrorCode  DMSetVec(DM dm,Vec x)
1569*08da532bSDmitry Karpeev {
1570*08da532bSDmitry Karpeev   PetscErrorCode ierr;
1571*08da532bSDmitry Karpeev   PetscFunctionBegin;
1572*08da532bSDmitry Karpeev   if (x) {
1573*08da532bSDmitry Karpeev     if (!dm->x) {
1574*08da532bSDmitry Karpeev       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1575*08da532bSDmitry Karpeev     }
1576*08da532bSDmitry Karpeev     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1577*08da532bSDmitry Karpeev   }
1578*08da532bSDmitry Karpeev   else if(dm->x) {
1579*08da532bSDmitry Karpeev     ierr = VecDestroy(&dm->x);  CHKERRQ(ierr);
1580*08da532bSDmitry Karpeev   }
1581*08da532bSDmitry Karpeev   PetscFunctionReturn(0);
1582*08da532bSDmitry Karpeev }
1583*08da532bSDmitry Karpeev 
1584*08da532bSDmitry Karpeev 
1585*08da532bSDmitry Karpeev #undef __FUNCT__
158647c6ae99SBarry Smith #define __FUNCT__ "DMComputeFunction"
158747c6ae99SBarry Smith /*@
158847c6ae99SBarry Smith     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
158947c6ae99SBarry Smith 
159047c6ae99SBarry Smith     Collective on DM
159147c6ae99SBarry Smith 
159247c6ae99SBarry Smith     Input Parameter:
159347c6ae99SBarry Smith +   dm - the DM object to destroy
159447c6ae99SBarry Smith .   x - the location where the function is evaluationed, may be ignored for linear problems
159547c6ae99SBarry Smith -   b - the vector to hold the right hand side entries
159647c6ae99SBarry Smith 
159747c6ae99SBarry Smith     Level: developer
159847c6ae99SBarry Smith 
1599e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
160047c6ae99SBarry Smith          DMSetJacobian()
160147c6ae99SBarry Smith 
160247c6ae99SBarry Smith @*/
16037087cfbeSBarry Smith PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
160447c6ae99SBarry Smith {
160547c6ae99SBarry Smith   PetscErrorCode ierr;
160647c6ae99SBarry Smith   PetscFunctionBegin;
1607171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
160847c6ae99SBarry Smith   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1609644e2e5bSBarry Smith   PetscStackPush("DM user function");
161047c6ae99SBarry Smith   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1611644e2e5bSBarry Smith   PetscStackPop;
161247c6ae99SBarry Smith   PetscFunctionReturn(0);
161347c6ae99SBarry Smith }
161447c6ae99SBarry Smith 
161547c6ae99SBarry Smith 
1616*08da532bSDmitry Karpeev 
161747c6ae99SBarry Smith #undef __FUNCT__
161847c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobian"
161947c6ae99SBarry Smith /*@
162047c6ae99SBarry Smith     DMComputeJacobian - compute the matrix entries for the solver
162147c6ae99SBarry Smith 
162247c6ae99SBarry Smith     Collective on DM
162347c6ae99SBarry Smith 
162447c6ae99SBarry Smith     Input Parameter:
162547c6ae99SBarry Smith +   dm - the DM object
1626cab2e9ccSBarry Smith .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
162747c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
162847c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
162947c6ae99SBarry Smith 
163047c6ae99SBarry Smith     Level: developer
163147c6ae99SBarry Smith 
1632e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
163347c6ae99SBarry Smith          DMSetFunction()
163447c6ae99SBarry Smith 
163547c6ae99SBarry Smith @*/
16367087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
163747c6ae99SBarry Smith {
163847c6ae99SBarry Smith   PetscErrorCode ierr;
163947c6ae99SBarry Smith 
164047c6ae99SBarry Smith   PetscFunctionBegin;
1641171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
164247c6ae99SBarry Smith   if (!dm->ops->jacobian) {
164347c6ae99SBarry Smith     ISColoring     coloring;
164447c6ae99SBarry Smith     MatFDColoring  fd;
164547c6ae99SBarry Smith 
1646171400e9SBarry Smith     ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
164747c6ae99SBarry Smith     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1648fcfd50ebSBarry Smith     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
164947c6ae99SBarry Smith     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
16500bdded8aSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
16510bdded8aSJed Brown     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
165271cd77b2SBarry Smith 
165347c6ae99SBarry Smith     dm->fd = fd;
165447c6ae99SBarry Smith     dm->ops->jacobian = DMComputeJacobianDefault;
16552533e041SBarry Smith 
165671cd77b2SBarry Smith     /* don't know why this is needed */
165771cd77b2SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
165847c6ae99SBarry Smith   }
165947c6ae99SBarry Smith   if (!x) x = dm->x;
166047c6ae99SBarry Smith   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1661cab2e9ccSBarry Smith 
166271cd77b2SBarry 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 */
1663649052a6SBarry Smith   if (x) {
1664cab2e9ccSBarry Smith     if (!dm->x) {
166571cd77b2SBarry Smith       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1666cab2e9ccSBarry Smith     }
1667cab2e9ccSBarry Smith     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1668649052a6SBarry Smith   }
1669a8248277SBarry Smith   if (A != B) {
1670a8248277SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1671a8248277SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1672a8248277SBarry Smith   }
167347c6ae99SBarry Smith   PetscFunctionReturn(0);
167447c6ae99SBarry Smith }
1675264ace61SBarry Smith 
1676cab2e9ccSBarry Smith 
1677264ace61SBarry Smith PetscFList DMList                       = PETSC_NULL;
1678264ace61SBarry Smith PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1679264ace61SBarry Smith 
1680264ace61SBarry Smith #undef __FUNCT__
1681264ace61SBarry Smith #define __FUNCT__ "DMSetType"
1682264ace61SBarry Smith /*@C
1683264ace61SBarry Smith   DMSetType - Builds a DM, for a particular DM implementation.
1684264ace61SBarry Smith 
1685264ace61SBarry Smith   Collective on DM
1686264ace61SBarry Smith 
1687264ace61SBarry Smith   Input Parameters:
1688264ace61SBarry Smith + dm     - The DM object
1689264ace61SBarry Smith - method - The name of the DM type
1690264ace61SBarry Smith 
1691264ace61SBarry Smith   Options Database Key:
1692264ace61SBarry Smith . -dm_type <type> - Sets the DM type; use -help for a list of available types
1693264ace61SBarry Smith 
1694264ace61SBarry Smith   Notes:
1695e1589f56SBarry Smith   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1696264ace61SBarry Smith 
1697264ace61SBarry Smith   Level: intermediate
1698264ace61SBarry Smith 
1699264ace61SBarry Smith .keywords: DM, set, type
1700264ace61SBarry Smith .seealso: DMGetType(), DMCreate()
1701264ace61SBarry Smith @*/
17027087cfbeSBarry Smith PetscErrorCode  DMSetType(DM dm, const DMType method)
1703264ace61SBarry Smith {
1704264ace61SBarry Smith   PetscErrorCode (*r)(DM);
1705264ace61SBarry Smith   PetscBool      match;
1706264ace61SBarry Smith   PetscErrorCode ierr;
1707264ace61SBarry Smith 
1708264ace61SBarry Smith   PetscFunctionBegin;
1709264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1710264ace61SBarry Smith   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1711264ace61SBarry Smith   if (match) PetscFunctionReturn(0);
1712264ace61SBarry Smith 
1713264ace61SBarry Smith   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
17144b91b6eaSBarry Smith   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1715264ace61SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1716264ace61SBarry Smith 
1717264ace61SBarry Smith   if (dm->ops->destroy) {
1718264ace61SBarry Smith     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1719b5c23020SJed Brown     dm->ops->destroy = PETSC_NULL;
1720264ace61SBarry Smith   }
1721264ace61SBarry Smith   ierr = (*r)(dm);CHKERRQ(ierr);
1722264ace61SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1723264ace61SBarry Smith   PetscFunctionReturn(0);
1724264ace61SBarry Smith }
1725264ace61SBarry Smith 
1726264ace61SBarry Smith #undef __FUNCT__
1727264ace61SBarry Smith #define __FUNCT__ "DMGetType"
1728264ace61SBarry Smith /*@C
1729264ace61SBarry Smith   DMGetType - Gets the DM type name (as a string) from the DM.
1730264ace61SBarry Smith 
1731264ace61SBarry Smith   Not Collective
1732264ace61SBarry Smith 
1733264ace61SBarry Smith   Input Parameter:
1734264ace61SBarry Smith . dm  - The DM
1735264ace61SBarry Smith 
1736264ace61SBarry Smith   Output Parameter:
1737264ace61SBarry Smith . type - The DM type name
1738264ace61SBarry Smith 
1739264ace61SBarry Smith   Level: intermediate
1740264ace61SBarry Smith 
1741264ace61SBarry Smith .keywords: DM, get, type, name
1742264ace61SBarry Smith .seealso: DMSetType(), DMCreate()
1743264ace61SBarry Smith @*/
17447087cfbeSBarry Smith PetscErrorCode  DMGetType(DM dm, const DMType *type)
1745264ace61SBarry Smith {
1746264ace61SBarry Smith   PetscErrorCode ierr;
1747264ace61SBarry Smith 
1748264ace61SBarry Smith   PetscFunctionBegin;
1749264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1750264ace61SBarry Smith   PetscValidCharPointer(type,2);
1751264ace61SBarry Smith   if (!DMRegisterAllCalled) {
1752264ace61SBarry Smith     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1753264ace61SBarry Smith   }
1754264ace61SBarry Smith   *type = ((PetscObject)dm)->type_name;
1755264ace61SBarry Smith   PetscFunctionReturn(0);
1756264ace61SBarry Smith }
1757264ace61SBarry Smith 
175867a56275SMatthew G Knepley #undef __FUNCT__
175967a56275SMatthew G Knepley #define __FUNCT__ "DMConvert"
176067a56275SMatthew G Knepley /*@C
176167a56275SMatthew G Knepley   DMConvert - Converts a DM to another DM, either of the same or different type.
176267a56275SMatthew G Knepley 
176367a56275SMatthew G Knepley   Collective on DM
176467a56275SMatthew G Knepley 
176567a56275SMatthew G Knepley   Input Parameters:
176667a56275SMatthew G Knepley + dm - the DM
176767a56275SMatthew G Knepley - newtype - new DM type (use "same" for the same type)
176867a56275SMatthew G Knepley 
176967a56275SMatthew G Knepley   Output Parameter:
177067a56275SMatthew G Knepley . M - pointer to new DM
177167a56275SMatthew G Knepley 
177267a56275SMatthew G Knepley   Notes:
177367a56275SMatthew G Knepley   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
177467a56275SMatthew G Knepley   the MPI communicator of the generated DM is always the same as the communicator
177567a56275SMatthew G Knepley   of the input DM.
177667a56275SMatthew G Knepley 
177767a56275SMatthew G Knepley   Level: intermediate
177867a56275SMatthew G Knepley 
177967a56275SMatthew G Knepley .seealso: DMCreate()
178067a56275SMatthew G Knepley @*/
178167a56275SMatthew G Knepley PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
178267a56275SMatthew G Knepley {
178367a56275SMatthew G Knepley   DM             B;
178467a56275SMatthew G Knepley   char           convname[256];
178567a56275SMatthew G Knepley   PetscBool      sametype, issame;
178667a56275SMatthew G Knepley   PetscErrorCode ierr;
178767a56275SMatthew G Knepley 
178867a56275SMatthew G Knepley   PetscFunctionBegin;
178967a56275SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
179067a56275SMatthew G Knepley   PetscValidType(dm,1);
179167a56275SMatthew G Knepley   PetscValidPointer(M,3);
179267a56275SMatthew G Knepley   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
179367a56275SMatthew G Knepley   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
179467a56275SMatthew G Knepley   {
179567a56275SMatthew G Knepley     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
179667a56275SMatthew G Knepley 
179767a56275SMatthew G Knepley     /*
179867a56275SMatthew G Knepley        Order of precedence:
179967a56275SMatthew G Knepley        1) See if a specialized converter is known to the current DM.
180067a56275SMatthew G Knepley        2) See if a specialized converter is known to the desired DM class.
180167a56275SMatthew G Knepley        3) See if a good general converter is registered for the desired class
180267a56275SMatthew G Knepley        4) See if a good general converter is known for the current matrix.
180367a56275SMatthew G Knepley        5) Use a really basic converter.
180467a56275SMatthew G Knepley     */
180567a56275SMatthew G Knepley 
180667a56275SMatthew G Knepley     /* 1) See if a specialized converter is known to the current DM and the desired class */
180767a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
180867a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
180967a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
181067a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
181167a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
181267a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
181367a56275SMatthew G Knepley     if (conv) goto foundconv;
181467a56275SMatthew G Knepley 
181567a56275SMatthew G Knepley     /* 2)  See if a specialized converter is known to the desired DM class. */
181667a56275SMatthew G Knepley     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
181767a56275SMatthew G Knepley     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
181867a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
181967a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
182067a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
182167a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
182267a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
182367a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
182467a56275SMatthew G Knepley     if (conv) {
1825fcfd50ebSBarry Smith       ierr = DMDestroy(&B);CHKERRQ(ierr);
182667a56275SMatthew G Knepley       goto foundconv;
182767a56275SMatthew G Knepley     }
182867a56275SMatthew G Knepley 
182967a56275SMatthew G Knepley #if 0
183067a56275SMatthew G Knepley     /* 3) See if a good general converter is registered for the desired class */
183167a56275SMatthew G Knepley     conv = B->ops->convertfrom;
1832fcfd50ebSBarry Smith     ierr = DMDestroy(&B);CHKERRQ(ierr);
183367a56275SMatthew G Knepley     if (conv) goto foundconv;
183467a56275SMatthew G Knepley 
183567a56275SMatthew G Knepley     /* 4) See if a good general converter is known for the current matrix */
183667a56275SMatthew G Knepley     if (dm->ops->convert) {
183767a56275SMatthew G Knepley       conv = dm->ops->convert;
183867a56275SMatthew G Knepley     }
183967a56275SMatthew G Knepley     if (conv) goto foundconv;
184067a56275SMatthew G Knepley #endif
184167a56275SMatthew G Knepley 
184267a56275SMatthew G Knepley     /* 5) Use a really basic converter. */
184367a56275SMatthew G Knepley     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
184467a56275SMatthew G Knepley 
184567a56275SMatthew G Knepley     foundconv:
184667a56275SMatthew G Knepley     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
184767a56275SMatthew G Knepley     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
184867a56275SMatthew G Knepley     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
184967a56275SMatthew G Knepley   }
185067a56275SMatthew G Knepley   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
185167a56275SMatthew G Knepley   PetscFunctionReturn(0);
185267a56275SMatthew G Knepley }
1853264ace61SBarry Smith 
1854264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
1855264ace61SBarry Smith 
1856264ace61SBarry Smith #undef __FUNCT__
1857264ace61SBarry Smith #define __FUNCT__ "DMRegister"
1858264ace61SBarry Smith /*@C
1859264ace61SBarry Smith   DMRegister - See DMRegisterDynamic()
1860264ace61SBarry Smith 
1861264ace61SBarry Smith   Level: advanced
1862264ace61SBarry Smith @*/
18637087cfbeSBarry Smith PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1864264ace61SBarry Smith {
1865264ace61SBarry Smith   char fullname[PETSC_MAX_PATH_LEN];
1866264ace61SBarry Smith   PetscErrorCode ierr;
1867264ace61SBarry Smith 
1868264ace61SBarry Smith   PetscFunctionBegin;
1869264ace61SBarry Smith   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1870264ace61SBarry Smith   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1871264ace61SBarry Smith   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1872264ace61SBarry Smith   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1873264ace61SBarry Smith   PetscFunctionReturn(0);
1874264ace61SBarry Smith }
1875264ace61SBarry Smith 
1876264ace61SBarry Smith 
1877264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
1878264ace61SBarry Smith #undef __FUNCT__
1879264ace61SBarry Smith #define __FUNCT__ "DMRegisterDestroy"
1880264ace61SBarry Smith /*@C
1881264ace61SBarry Smith    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1882264ace61SBarry Smith 
1883264ace61SBarry Smith    Not Collective
1884264ace61SBarry Smith 
1885264ace61SBarry Smith    Level: advanced
1886264ace61SBarry Smith 
1887264ace61SBarry Smith .keywords: DM, register, destroy
1888264ace61SBarry Smith .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1889264ace61SBarry Smith @*/
18907087cfbeSBarry Smith PetscErrorCode  DMRegisterDestroy(void)
1891264ace61SBarry Smith {
1892264ace61SBarry Smith   PetscErrorCode ierr;
1893264ace61SBarry Smith 
1894264ace61SBarry Smith   PetscFunctionBegin;
1895264ace61SBarry Smith   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1896264ace61SBarry Smith   DMRegisterAllCalled = PETSC_FALSE;
1897264ace61SBarry Smith   PetscFunctionReturn(0);
1898264ace61SBarry Smith }
189923f975d1SBarry Smith 
190023f975d1SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1901c6db04a5SJed Brown #include <mex.h>
190223f975d1SBarry Smith 
19033014e516SBarry Smith typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
190423f975d1SBarry Smith 
190523f975d1SBarry Smith #undef __FUNCT__
190623f975d1SBarry Smith #define __FUNCT__ "DMComputeFunction_Matlab"
190723f975d1SBarry Smith /*
190823f975d1SBarry Smith    DMComputeFunction_Matlab - Calls the function that has been set with
190923f975d1SBarry Smith                          DMSetFunctionMatlab().
191023f975d1SBarry Smith 
191123f975d1SBarry Smith    For linear problems x is null
191223f975d1SBarry Smith 
191323f975d1SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
191423f975d1SBarry Smith */
19157087cfbeSBarry Smith PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
191623f975d1SBarry Smith {
191723f975d1SBarry Smith   PetscErrorCode    ierr;
191823f975d1SBarry Smith   DMMatlabContext   *sctx;
191923f975d1SBarry Smith   int               nlhs = 1,nrhs = 4;
192023f975d1SBarry Smith   mxArray	    *plhs[1],*prhs[4];
192123f975d1SBarry Smith   long long int     lx = 0,ly = 0,ls = 0;
192223f975d1SBarry Smith 
192323f975d1SBarry Smith   PetscFunctionBegin;
192423f975d1SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
192523f975d1SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
192623f975d1SBarry Smith   PetscCheckSameComm(dm,1,y,3);
192723f975d1SBarry Smith 
192823f975d1SBarry Smith   /* call Matlab function in ctx with arguments x and y */
19291b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
193023f975d1SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
193123f975d1SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
19323014e516SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
193323f975d1SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
193423f975d1SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
193523f975d1SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)ly);
193623f975d1SBarry Smith   prhs[3] =  mxCreateString(sctx->funcname);
1937b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
193823f975d1SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
193923f975d1SBarry Smith   mxDestroyArray(prhs[0]);
194023f975d1SBarry Smith   mxDestroyArray(prhs[1]);
194123f975d1SBarry Smith   mxDestroyArray(prhs[2]);
194223f975d1SBarry Smith   mxDestroyArray(prhs[3]);
194323f975d1SBarry Smith   mxDestroyArray(plhs[0]);
194423f975d1SBarry Smith   PetscFunctionReturn(0);
194523f975d1SBarry Smith }
194623f975d1SBarry Smith 
194723f975d1SBarry Smith 
194823f975d1SBarry Smith #undef __FUNCT__
194923f975d1SBarry Smith #define __FUNCT__ "DMSetFunctionMatlab"
195023f975d1SBarry Smith /*
195123f975d1SBarry Smith    DMSetFunctionMatlab - Sets the function evaluation routine
195223f975d1SBarry Smith 
195323f975d1SBarry Smith */
19547087cfbeSBarry Smith PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
195523f975d1SBarry Smith {
195623f975d1SBarry Smith   PetscErrorCode    ierr;
195723f975d1SBarry Smith   DMMatlabContext   *sctx;
195823f975d1SBarry Smith 
195923f975d1SBarry Smith   PetscFunctionBegin;
1960171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
196123f975d1SBarry Smith   /* currently sctx is memory bleed */
19621b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
19633014e516SBarry Smith   if (!sctx) {
196423f975d1SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
19653014e516SBarry Smith   }
196623f975d1SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
19671b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
196823f975d1SBarry Smith   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
196923f975d1SBarry Smith   PetscFunctionReturn(0);
197023f975d1SBarry Smith }
19713014e516SBarry Smith 
19723014e516SBarry Smith #undef __FUNCT__
19733014e516SBarry Smith #define __FUNCT__ "DMComputeJacobian_Matlab"
19743014e516SBarry Smith /*
19753014e516SBarry Smith    DMComputeJacobian_Matlab - Calls the function that has been set with
19763014e516SBarry Smith                          DMSetJacobianMatlab().
19773014e516SBarry Smith 
19783014e516SBarry Smith    For linear problems x is null
19793014e516SBarry Smith 
19803014e516SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
19813014e516SBarry Smith */
19827087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
19833014e516SBarry Smith {
19843014e516SBarry Smith   PetscErrorCode    ierr;
19853014e516SBarry Smith   DMMatlabContext   *sctx;
19863014e516SBarry Smith   int               nlhs = 2,nrhs = 5;
19873014e516SBarry Smith   mxArray	    *plhs[2],*prhs[5];
19883014e516SBarry Smith   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
19893014e516SBarry Smith 
19903014e516SBarry Smith   PetscFunctionBegin;
19913014e516SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
19923014e516SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
19933014e516SBarry Smith 
1994e3c5b3baSBarry Smith   /* call MATLAB function in ctx with arguments x, A, and B */
19951b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
19963014e516SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
19973014e516SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
19983014e516SBarry Smith   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
19993014e516SBarry Smith   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
20003014e516SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
20013014e516SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
20023014e516SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lA);
20033014e516SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lB);
20043014e516SBarry Smith   prhs[4] =  mxCreateString(sctx->jacname);
2005b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
2006c980e822SBarry Smith   *str    =  (MatStructure) mxGetScalar(plhs[0]);
2007c088a8dcSBarry Smith   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
20083014e516SBarry Smith   mxDestroyArray(prhs[0]);
20093014e516SBarry Smith   mxDestroyArray(prhs[1]);
20103014e516SBarry Smith   mxDestroyArray(prhs[2]);
20113014e516SBarry Smith   mxDestroyArray(prhs[3]);
20123014e516SBarry Smith   mxDestroyArray(prhs[4]);
20133014e516SBarry Smith   mxDestroyArray(plhs[0]);
20143014e516SBarry Smith   mxDestroyArray(plhs[1]);
20153014e516SBarry Smith   PetscFunctionReturn(0);
20163014e516SBarry Smith }
20173014e516SBarry Smith 
20183014e516SBarry Smith 
20193014e516SBarry Smith #undef __FUNCT__
20203014e516SBarry Smith #define __FUNCT__ "DMSetJacobianMatlab"
20213014e516SBarry Smith /*
20223014e516SBarry Smith    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
20233014e516SBarry Smith 
20243014e516SBarry Smith */
20257087cfbeSBarry Smith PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
20263014e516SBarry Smith {
20273014e516SBarry Smith   PetscErrorCode    ierr;
20283014e516SBarry Smith   DMMatlabContext   *sctx;
20293014e516SBarry Smith 
20303014e516SBarry Smith   PetscFunctionBegin;
2031171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
20323014e516SBarry Smith   /* currently sctx is memory bleed */
20331b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
20343014e516SBarry Smith   if (!sctx) {
20353014e516SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
20363014e516SBarry Smith   }
20373014e516SBarry Smith   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
20381b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
20393014e516SBarry Smith   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
20403014e516SBarry Smith   PetscFunctionReturn(0);
20413014e516SBarry Smith }
204223f975d1SBarry Smith #endif
2043b859378eSBarry Smith 
2044b859378eSBarry Smith #undef __FUNCT__
2045b859378eSBarry Smith #define __FUNCT__ "DMLoad"
2046b859378eSBarry Smith /*@C
2047b859378eSBarry Smith   DMLoad - Loads a DM that has been stored in binary or HDF5 format
2048b859378eSBarry Smith   with DMView().
2049b859378eSBarry Smith 
2050b859378eSBarry Smith   Collective on PetscViewer
2051b859378eSBarry Smith 
2052b859378eSBarry Smith   Input Parameters:
2053b859378eSBarry Smith + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2054b859378eSBarry Smith            some related function before a call to DMLoad().
2055b859378eSBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2056b859378eSBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2057b859378eSBarry Smith 
2058b859378eSBarry Smith    Level: intermediate
2059b859378eSBarry Smith 
2060b859378eSBarry Smith   Notes:
2061b859378eSBarry Smith   Defaults to the DM DA.
2062b859378eSBarry Smith 
2063b859378eSBarry Smith   Notes for advanced users:
2064b859378eSBarry Smith   Most users should not need to know the details of the binary storage
2065b859378eSBarry Smith   format, since DMLoad() and DMView() completely hide these details.
2066b859378eSBarry Smith   But for anyone who's interested, the standard binary matrix storage
2067b859378eSBarry Smith   format is
2068b859378eSBarry Smith .vb
2069b859378eSBarry Smith      has not yet been determined
2070b859378eSBarry Smith .ve
2071b859378eSBarry Smith 
2072b859378eSBarry Smith    In addition, PETSc automatically does the byte swapping for
2073b859378eSBarry Smith machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
2074b859378eSBarry Smith linux, Windows and the paragon; thus if you write your own binary
2075b859378eSBarry Smith read/write routines you have to swap the bytes; see PetscBinaryRead()
2076b859378eSBarry Smith and PetscBinaryWrite() to see how this may be done.
2077b859378eSBarry Smith 
2078b859378eSBarry Smith   Concepts: vector^loading from file
2079b859378eSBarry Smith 
2080b859378eSBarry Smith .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2081b859378eSBarry Smith @*/
2082b859378eSBarry Smith PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2083b859378eSBarry Smith {
2084b859378eSBarry Smith   PetscErrorCode ierr;
2085b859378eSBarry Smith 
2086b859378eSBarry Smith   PetscFunctionBegin;
2087b859378eSBarry Smith   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2088b859378eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2089b859378eSBarry Smith 
2090b859378eSBarry Smith   if (!((PetscObject)newdm)->type_name) {
2091b859378eSBarry Smith     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
2092b859378eSBarry Smith   }
2093b859378eSBarry Smith   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2094b859378eSBarry Smith   PetscFunctionReturn(0);
2095b859378eSBarry Smith }
2096b859378eSBarry Smith 
20977da65231SMatthew G Knepley /******************************** FEM Support **********************************/
20987da65231SMatthew G Knepley 
20997da65231SMatthew G Knepley #undef __FUNCT__
21007da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellVector"
21017da65231SMatthew G Knepley PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) {
21021d47ebbbSSatish Balay   PetscInt       f;
21031b30c384SMatthew G Knepley   PetscErrorCode ierr;
21041b30c384SMatthew G Knepley 
21057da65231SMatthew G Knepley   PetscFunctionBegin;
210674778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
21071d47ebbbSSatish Balay   for(f = 0; f < len; ++f) {
210874778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
21097da65231SMatthew G Knepley   }
21107da65231SMatthew G Knepley   PetscFunctionReturn(0);
21117da65231SMatthew G Knepley }
21127da65231SMatthew G Knepley 
21137da65231SMatthew G Knepley #undef __FUNCT__
21147da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellMatrix"
21157da65231SMatthew G Knepley PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) {
21161b30c384SMatthew G Knepley   PetscInt       f, g;
21177da65231SMatthew G Knepley   PetscErrorCode ierr;
21187da65231SMatthew G Knepley 
21197da65231SMatthew G Knepley   PetscFunctionBegin;
212074778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
21211d47ebbbSSatish Balay   for(f = 0; f < rows; ++f) {
212274778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
21231d47ebbbSSatish Balay     for(g = 0; g < cols; ++g) {
212474778d6cSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
21257da65231SMatthew G Knepley     }
212674778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
21277da65231SMatthew G Knepley   }
21287da65231SMatthew G Knepley   PetscFunctionReturn(0);
21297da65231SMatthew G Knepley }
2130