xref: /petsc/src/dm/interface/dm.c (revision 521d9a4c1fa86e77f87c59897acde77f37c691c3)
108da532bSDmitry Karpeev #include <petscsnes.h>
2b45d2f2cSJed Brown #include <petsc-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 
44e7c4fc90SDmitry Karpeev 
45a89ea682SMatthew G Knepley   v->workSize     = 0;
46a89ea682SMatthew G Knepley   v->workArray    = PETSC_NULL;
471411c6eeSJed Brown   v->ltogmap      = PETSC_NULL;
481411c6eeSJed Brown   v->ltogmapb     = PETSC_NULL;
491411c6eeSJed Brown   v->bs           = 1;
50171400e9SBarry Smith   v->coloringtype = IS_COLORING_GLOBAL;
511411c6eeSJed Brown 
521411c6eeSJed Brown   *dm = v;
53a4121054SBarry Smith   PetscFunctionReturn(0);
54a4121054SBarry Smith }
55a4121054SBarry Smith 
56a4121054SBarry Smith 
57a4121054SBarry Smith #undef __FUNCT__
589a42bb27SBarry Smith #define __FUNCT__ "DMSetVecType"
599a42bb27SBarry Smith /*@C
60564755cdSBarry Smith        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
619a42bb27SBarry Smith 
62aa219208SBarry Smith    Logically Collective on DMDA
639a42bb27SBarry Smith 
649a42bb27SBarry Smith    Input Parameter:
659a42bb27SBarry Smith +  da - initial distributed array
668154be41SBarry Smith .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
679a42bb27SBarry Smith 
689a42bb27SBarry Smith    Options Database:
69dd85299cSBarry Smith .   -dm_vec_type ctype
709a42bb27SBarry Smith 
719a42bb27SBarry Smith    Level: intermediate
729a42bb27SBarry Smith 
73aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
749a42bb27SBarry Smith @*/
757087cfbeSBarry Smith PetscErrorCode  DMSetVecType(DM da,const VecType ctype)
769a42bb27SBarry Smith {
779a42bb27SBarry Smith   PetscErrorCode ierr;
789a42bb27SBarry Smith 
799a42bb27SBarry Smith   PetscFunctionBegin;
809a42bb27SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
819a42bb27SBarry Smith   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
829a42bb27SBarry Smith   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
839a42bb27SBarry Smith   PetscFunctionReturn(0);
849a42bb27SBarry Smith }
859a42bb27SBarry Smith 
869a42bb27SBarry Smith #undef __FUNCT__
87*521d9a4cSLisandro Dalcin #define __FUNCT__ "DMSetMatType"
88*521d9a4cSLisandro Dalcin /*@C
89*521d9a4cSLisandro Dalcin        DMSetMatType - Sets the type of matrix created with DMCreateMatrix()
90*521d9a4cSLisandro Dalcin 
91*521d9a4cSLisandro Dalcin    Logically Collective on DM
92*521d9a4cSLisandro Dalcin 
93*521d9a4cSLisandro Dalcin    Input Parameter:
94*521d9a4cSLisandro Dalcin +  dm - the DM context
95*521d9a4cSLisandro Dalcin .  ctype - the matrix type
96*521d9a4cSLisandro Dalcin 
97*521d9a4cSLisandro Dalcin    Options Database:
98*521d9a4cSLisandro Dalcin .   -dm_mat_type ctype
99*521d9a4cSLisandro Dalcin 
100*521d9a4cSLisandro Dalcin    Level: intermediate
101*521d9a4cSLisandro Dalcin 
102*521d9a4cSLisandro Dalcin .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType
103*521d9a4cSLisandro Dalcin @*/
104*521d9a4cSLisandro Dalcin PetscErrorCode  DMSetMatType(DM dm,const MatType ctype)
105*521d9a4cSLisandro Dalcin {
106*521d9a4cSLisandro Dalcin   PetscErrorCode ierr;
107*521d9a4cSLisandro Dalcin   PetscFunctionBegin;
108*521d9a4cSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
109*521d9a4cSLisandro Dalcin   ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
110*521d9a4cSLisandro Dalcin   ierr = PetscStrallocpy(ctype,&dm->mattype);CHKERRQ(ierr);
111*521d9a4cSLisandro Dalcin   PetscFunctionReturn(0);
112*521d9a4cSLisandro Dalcin }
113*521d9a4cSLisandro Dalcin 
114*521d9a4cSLisandro Dalcin #undef __FUNCT__
1159a42bb27SBarry Smith #define __FUNCT__ "DMSetOptionsPrefix"
1169a42bb27SBarry Smith /*@C
1179a42bb27SBarry Smith    DMSetOptionsPrefix - Sets the prefix used for searching for all
118aa219208SBarry Smith    DMDA options in the database.
1199a42bb27SBarry Smith 
120aa219208SBarry Smith    Logically Collective on DMDA
1219a42bb27SBarry Smith 
1229a42bb27SBarry Smith    Input Parameter:
123aa219208SBarry Smith +  da - the DMDA context
1249a42bb27SBarry Smith -  prefix - the prefix to prepend to all option names
1259a42bb27SBarry Smith 
1269a42bb27SBarry Smith    Notes:
1279a42bb27SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1289a42bb27SBarry Smith    The first character of all runtime options is AUTOMATICALLY the hyphen.
1299a42bb27SBarry Smith 
1309a42bb27SBarry Smith    Level: advanced
1319a42bb27SBarry Smith 
132aa219208SBarry Smith .keywords: DMDA, set, options, prefix, database
1339a42bb27SBarry Smith 
1349a42bb27SBarry Smith .seealso: DMSetFromOptions()
1359a42bb27SBarry Smith @*/
1367087cfbeSBarry Smith PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
1379a42bb27SBarry Smith {
1389a42bb27SBarry Smith   PetscErrorCode ierr;
1399a42bb27SBarry Smith 
1409a42bb27SBarry Smith   PetscFunctionBegin;
1419a42bb27SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1429a42bb27SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
1439a42bb27SBarry Smith   PetscFunctionReturn(0);
1449a42bb27SBarry Smith }
1459a42bb27SBarry Smith 
1469a42bb27SBarry Smith #undef __FUNCT__
14747c6ae99SBarry Smith #define __FUNCT__ "DMDestroy"
14847c6ae99SBarry Smith /*@
149aa219208SBarry Smith     DMDestroy - Destroys a vector packer or DMDA.
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith     Collective on DM
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith     Input Parameter:
15447c6ae99SBarry Smith .   dm - the DM object to destroy
15547c6ae99SBarry Smith 
15647c6ae99SBarry Smith     Level: developer
15747c6ae99SBarry Smith 
158e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith @*/
161fcfd50ebSBarry Smith PetscErrorCode  DMDestroy(DM *dm)
16247c6ae99SBarry Smith {
163732e2eb9SMatthew G Knepley   PetscInt       i, cnt = 0;
164b17ce1afSJed Brown   DMCoarsenHookLink link,next;
16547c6ae99SBarry Smith   PetscErrorCode ierr;
16647c6ae99SBarry Smith 
16747c6ae99SBarry Smith   PetscFunctionBegin;
1686bf464f9SBarry Smith   if (!*dm) PetscFunctionReturn(0);
1696bf464f9SBarry Smith   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
17087e657c6SBarry Smith 
17187e657c6SBarry Smith   /* count all the circular references of DM and its contained Vecs */
172732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1736bf464f9SBarry Smith     if ((*dm)->localin[i])  {cnt++;}
1746bf464f9SBarry Smith     if ((*dm)->globalin[i]) {cnt++;}
175732e2eb9SMatthew G Knepley   }
17671cd77b2SBarry Smith   if ((*dm)->x) {
17771cd77b2SBarry Smith     PetscObject obj;
17871cd77b2SBarry Smith     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
17971cd77b2SBarry Smith     if (obj == (PetscObject)*dm) cnt++;
18071cd77b2SBarry Smith   }
181732e2eb9SMatthew G Knepley 
1826bf464f9SBarry Smith   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
183732e2eb9SMatthew G Knepley   /*
184732e2eb9SMatthew G Knepley      Need this test because the dm references the vectors that
185732e2eb9SMatthew G Knepley      reference the dm, so destroying the dm calls destroy on the
186732e2eb9SMatthew G Knepley      vectors that cause another destroy on the dm
187732e2eb9SMatthew G Knepley   */
1886bf464f9SBarry Smith   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
1896bf464f9SBarry Smith   ((PetscObject) (*dm))->refct = 0;
190732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1916bf464f9SBarry Smith     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
1926bf464f9SBarry Smith     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
193732e2eb9SMatthew G Knepley   }
1941a266240SBarry Smith 
195b17ce1afSJed Brown   /* Destroy the list of hooks */
196b17ce1afSJed Brown   for (link=(*dm)->coarsenhook; link; link=next) {
197b17ce1afSJed Brown     next = link->next;
198b17ce1afSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
199b17ce1afSJed Brown   }
200b17ce1afSJed Brown   (*dm)->coarsenhook = PETSC_NULL;
201b17ce1afSJed Brown 
2021a266240SBarry Smith   if ((*dm)->ctx && (*dm)->ctxdestroy) {
2031a266240SBarry Smith     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
2041a266240SBarry Smith   }
20587e657c6SBarry Smith   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
20671cd77b2SBarry Smith   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
2074dcab191SBarry Smith   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
2086bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
2096bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
2106bf464f9SBarry Smith   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
211073dac72SJed Brown   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
212a89ea682SMatthew G Knepley   ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr);
213732e2eb9SMatthew G Knepley   /* if memory was published with AMS then destroy it */
2146bf464f9SBarry Smith   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
215732e2eb9SMatthew G Knepley 
2166bf464f9SBarry Smith   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
2176bf464f9SBarry Smith   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
218732e2eb9SMatthew G Knepley   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
21947c6ae99SBarry Smith   PetscFunctionReturn(0);
22047c6ae99SBarry Smith }
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith #undef __FUNCT__
223d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp"
224d7bf68aeSBarry Smith /*@
225d7bf68aeSBarry Smith     DMSetUp - sets up the data structures inside a DM object
226d7bf68aeSBarry Smith 
227d7bf68aeSBarry Smith     Collective on DM
228d7bf68aeSBarry Smith 
229d7bf68aeSBarry Smith     Input Parameter:
230d7bf68aeSBarry Smith .   dm - the DM object to setup
231d7bf68aeSBarry Smith 
232d7bf68aeSBarry Smith     Level: developer
233d7bf68aeSBarry Smith 
234e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
235d7bf68aeSBarry Smith 
236d7bf68aeSBarry Smith @*/
2377087cfbeSBarry Smith PetscErrorCode  DMSetUp(DM dm)
238d7bf68aeSBarry Smith {
239d7bf68aeSBarry Smith   PetscErrorCode ierr;
240d7bf68aeSBarry Smith 
241d7bf68aeSBarry Smith   PetscFunctionBegin;
242171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2438387afaaSJed Brown   if (dm->setupcalled) PetscFunctionReturn(0);
244d7bf68aeSBarry Smith   if (dm->ops->setup) {
245d7bf68aeSBarry Smith     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
246d7bf68aeSBarry Smith   }
2478387afaaSJed Brown   dm->setupcalled = PETSC_TRUE;
248d7bf68aeSBarry Smith   PetscFunctionReturn(0);
249d7bf68aeSBarry Smith }
250d7bf68aeSBarry Smith 
251d7bf68aeSBarry Smith #undef __FUNCT__
252d7bf68aeSBarry Smith #define __FUNCT__ "DMSetFromOptions"
253d7bf68aeSBarry Smith /*@
254d7bf68aeSBarry Smith     DMSetFromOptions - sets parameters in a DM from the options database
255d7bf68aeSBarry Smith 
256d7bf68aeSBarry Smith     Collective on DM
257d7bf68aeSBarry Smith 
258d7bf68aeSBarry Smith     Input Parameter:
259d7bf68aeSBarry Smith .   dm - the DM object to set options for
260d7bf68aeSBarry Smith 
261732e2eb9SMatthew G Knepley     Options Database:
262dd85299cSBarry Smith +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
263dd85299cSBarry Smith .   -dm_vec_type <type>  type of vector to create inside DM
264171400e9SBarry Smith .   -dm_mat_type <type>  type of matrix to create inside DM
265171400e9SBarry Smith -   -dm_coloring_type <global or ghosted>
266732e2eb9SMatthew G Knepley 
267d7bf68aeSBarry Smith     Level: developer
268d7bf68aeSBarry Smith 
269e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
270d7bf68aeSBarry Smith 
271d7bf68aeSBarry Smith @*/
2727087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions(DM dm)
273d7bf68aeSBarry Smith {
27467ad5babSMatthew G Knepley   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg;
275d7bf68aeSBarry Smith   PetscErrorCode ierr;
276f9ba7244SBarry Smith   char           typeName[256] = MATAIJ;
277d7bf68aeSBarry Smith 
278d7bf68aeSBarry Smith   PetscFunctionBegin;
279171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2803194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
28182fcb398SMatthew G Knepley     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
282c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr);
283c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr);
28467ad5babSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr);
285073dac72SJed 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);
286f9ba7244SBarry Smith     ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
287f9ba7244SBarry Smith     if (flg) {
288f9ba7244SBarry Smith       ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
289f9ba7244SBarry Smith     }
290*521d9a4cSLisandro Dalcin     ierr = PetscOptionsList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype?dm->mattype:typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr);
291073dac72SJed Brown     if (flg) {
292*521d9a4cSLisandro Dalcin       ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr);
293073dac72SJed Brown     }
2941b89239cSHong 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);
295f9ba7244SBarry Smith     if (dm->ops->setfromoptions) {
296f9ba7244SBarry Smith       ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
297f9ba7244SBarry Smith     }
298f9ba7244SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
299f9ba7244SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
30082fcb398SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
30182fcb398SMatthew G Knepley   if (flg1) {
30282fcb398SMatthew G Knepley     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
30382fcb398SMatthew G Knepley   }
304c4721b0eSMatthew G Knepley   if (flg2) {
305c4721b0eSMatthew G Knepley     PetscViewer viewer;
306c4721b0eSMatthew G Knepley 
307c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
308c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
309c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
310c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
311c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
312c4721b0eSMatthew G Knepley   }
313c4721b0eSMatthew G Knepley   if (flg3) {
314c4721b0eSMatthew G Knepley     PetscViewer viewer;
315c4721b0eSMatthew G Knepley 
316c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
317c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
318c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
319c4721b0eSMatthew G Knepley     ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr);
320c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
321c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
322c4721b0eSMatthew G Knepley   }
32367ad5babSMatthew G Knepley   if (flg4) {
32467ad5babSMatthew G Knepley     PetscViewer viewer;
32567ad5babSMatthew G Knepley 
32667ad5babSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
32767ad5babSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
32867ad5babSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr);
32967ad5babSMatthew G Knepley     ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr);
33067ad5babSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
33167ad5babSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
33267ad5babSMatthew G Knepley   }
333d7bf68aeSBarry Smith   PetscFunctionReturn(0);
334d7bf68aeSBarry Smith }
335d7bf68aeSBarry Smith 
336d7bf68aeSBarry Smith #undef __FUNCT__
33747c6ae99SBarry Smith #define __FUNCT__ "DMView"
338fc9bc008SSatish Balay /*@C
339aa219208SBarry Smith     DMView - Views a vector packer or DMDA.
34047c6ae99SBarry Smith 
34147c6ae99SBarry Smith     Collective on DM
34247c6ae99SBarry Smith 
34347c6ae99SBarry Smith     Input Parameter:
34447c6ae99SBarry Smith +   dm - the DM object to view
34547c6ae99SBarry Smith -   v - the viewer
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith     Level: developer
34847c6ae99SBarry Smith 
349e727c939SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith @*/
3527087cfbeSBarry Smith PetscErrorCode  DMView(DM dm,PetscViewer v)
35347c6ae99SBarry Smith {
35447c6ae99SBarry Smith   PetscErrorCode ierr;
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith   PetscFunctionBegin;
357171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3583014e516SBarry Smith  if (!v) {
3593014e516SBarry Smith     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
3603014e516SBarry Smith   }
3610c010503SBarry Smith   if (dm->ops->view) {
3620c010503SBarry Smith     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
36347c6ae99SBarry Smith   }
36447c6ae99SBarry Smith   PetscFunctionReturn(0);
36547c6ae99SBarry Smith }
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith #undef __FUNCT__
36847c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector"
36947c6ae99SBarry Smith /*@
370aa219208SBarry Smith     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
37147c6ae99SBarry Smith 
37247c6ae99SBarry Smith     Collective on DM
37347c6ae99SBarry Smith 
37447c6ae99SBarry Smith     Input Parameter:
37547c6ae99SBarry Smith .   dm - the DM object
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith     Output Parameter:
37847c6ae99SBarry Smith .   vec - the global vector
37947c6ae99SBarry Smith 
380073dac72SJed Brown     Level: beginner
38147c6ae99SBarry Smith 
382e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
38347c6ae99SBarry Smith 
38447c6ae99SBarry Smith @*/
3857087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector(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->createglobalvector)(dm,vec);CHKERRQ(ierr);
39247c6ae99SBarry Smith   PetscFunctionReturn(0);
39347c6ae99SBarry Smith }
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith #undef __FUNCT__
39647c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector"
39747c6ae99SBarry Smith /*@
398aa219208SBarry Smith     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
39947c6ae99SBarry Smith 
40047c6ae99SBarry Smith     Not Collective
40147c6ae99SBarry Smith 
40247c6ae99SBarry Smith     Input Parameter:
40347c6ae99SBarry Smith .   dm - the DM object
40447c6ae99SBarry Smith 
40547c6ae99SBarry Smith     Output Parameter:
40647c6ae99SBarry Smith .   vec - the local vector
40747c6ae99SBarry Smith 
408073dac72SJed Brown     Level: beginner
40947c6ae99SBarry Smith 
410e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
41147c6ae99SBarry Smith 
41247c6ae99SBarry Smith @*/
4137087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
41447c6ae99SBarry Smith {
41547c6ae99SBarry Smith   PetscErrorCode ierr;
41647c6ae99SBarry Smith 
41747c6ae99SBarry Smith   PetscFunctionBegin;
418171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
41947c6ae99SBarry Smith   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
42047c6ae99SBarry Smith   PetscFunctionReturn(0);
42147c6ae99SBarry Smith }
42247c6ae99SBarry Smith 
42347c6ae99SBarry Smith #undef __FUNCT__
4241411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping"
4251411c6eeSJed Brown /*@
4261411c6eeSJed Brown    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
4271411c6eeSJed Brown 
4281411c6eeSJed Brown    Collective on DM
4291411c6eeSJed Brown 
4301411c6eeSJed Brown    Input Parameter:
4311411c6eeSJed Brown .  dm - the DM that provides the mapping
4321411c6eeSJed Brown 
4331411c6eeSJed Brown    Output Parameter:
4341411c6eeSJed Brown .  ltog - the mapping
4351411c6eeSJed Brown 
4361411c6eeSJed Brown    Level: intermediate
4371411c6eeSJed Brown 
4381411c6eeSJed Brown    Notes:
4391411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMapping() or
4401411c6eeSJed Brown    MatSetLocalToGlobalMapping().
4411411c6eeSJed Brown 
4421411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
4431411c6eeSJed Brown @*/
4447087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
4451411c6eeSJed Brown {
4461411c6eeSJed Brown   PetscErrorCode ierr;
4471411c6eeSJed Brown 
4481411c6eeSJed Brown   PetscFunctionBegin;
4491411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4501411c6eeSJed Brown   PetscValidPointer(ltog,2);
4511411c6eeSJed Brown   if (!dm->ltogmap) {
4521411c6eeSJed Brown     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
4531411c6eeSJed Brown     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
4541411c6eeSJed Brown   }
4551411c6eeSJed Brown   *ltog = dm->ltogmap;
4561411c6eeSJed Brown   PetscFunctionReturn(0);
4571411c6eeSJed Brown }
4581411c6eeSJed Brown 
4591411c6eeSJed Brown #undef __FUNCT__
4601411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
4611411c6eeSJed Brown /*@
4621411c6eeSJed Brown    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
4631411c6eeSJed Brown 
4641411c6eeSJed Brown    Collective on DM
4651411c6eeSJed Brown 
4661411c6eeSJed Brown    Input Parameter:
4671411c6eeSJed Brown .  da - the distributed array that provides the mapping
4681411c6eeSJed Brown 
4691411c6eeSJed Brown    Output Parameter:
4701411c6eeSJed Brown .  ltog - the block mapping
4711411c6eeSJed Brown 
4721411c6eeSJed Brown    Level: intermediate
4731411c6eeSJed Brown 
4741411c6eeSJed Brown    Notes:
4751411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
4761411c6eeSJed Brown    MatSetLocalToGlobalMappingBlock().
4771411c6eeSJed Brown 
4781411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
4791411c6eeSJed Brown @*/
4807087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
4811411c6eeSJed Brown {
4821411c6eeSJed Brown   PetscErrorCode ierr;
4831411c6eeSJed Brown 
4841411c6eeSJed Brown   PetscFunctionBegin;
4851411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4861411c6eeSJed Brown   PetscValidPointer(ltog,2);
4871411c6eeSJed Brown   if (!dm->ltogmapb) {
4881411c6eeSJed Brown     PetscInt bs;
4891411c6eeSJed Brown     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
4901411c6eeSJed Brown     if (bs > 1) {
4911411c6eeSJed Brown       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
4921411c6eeSJed Brown       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
4931411c6eeSJed Brown     } else {
4941411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
4951411c6eeSJed Brown       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
4961411c6eeSJed Brown     }
4971411c6eeSJed Brown   }
4981411c6eeSJed Brown   *ltog = dm->ltogmapb;
4991411c6eeSJed Brown   PetscFunctionReturn(0);
5001411c6eeSJed Brown }
5011411c6eeSJed Brown 
5021411c6eeSJed Brown #undef __FUNCT__
5031411c6eeSJed Brown #define __FUNCT__ "DMGetBlockSize"
5041411c6eeSJed Brown /*@
5051411c6eeSJed Brown    DMGetBlockSize - Gets the inherent block size associated with a DM
5061411c6eeSJed Brown 
5071411c6eeSJed Brown    Not Collective
5081411c6eeSJed Brown 
5091411c6eeSJed Brown    Input Parameter:
5101411c6eeSJed Brown .  dm - the DM with block structure
5111411c6eeSJed Brown 
5121411c6eeSJed Brown    Output Parameter:
5131411c6eeSJed Brown .  bs - the block size, 1 implies no exploitable block structure
5141411c6eeSJed Brown 
5151411c6eeSJed Brown    Level: intermediate
5161411c6eeSJed Brown 
5171411c6eeSJed Brown .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
5181411c6eeSJed Brown @*/
5197087cfbeSBarry Smith PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
5201411c6eeSJed Brown {
5211411c6eeSJed Brown   PetscFunctionBegin;
5221411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5231411c6eeSJed Brown   PetscValidPointer(bs,2);
5241411c6eeSJed 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");
5251411c6eeSJed Brown   *bs = dm->bs;
5261411c6eeSJed Brown   PetscFunctionReturn(0);
5271411c6eeSJed Brown }
5281411c6eeSJed Brown 
5291411c6eeSJed Brown #undef __FUNCT__
530e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation"
53147c6ae99SBarry Smith /*@
532e727c939SJed Brown     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith     Collective on DM
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith     Input Parameter:
53747c6ae99SBarry Smith +   dm1 - the DM object
53847c6ae99SBarry Smith -   dm2 - the second, finer DM object
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith     Output Parameter:
54147c6ae99SBarry Smith +  mat - the interpolation
54247c6ae99SBarry Smith -  vec - the scaling (optional)
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith     Level: developer
54547c6ae99SBarry Smith 
54685afcc9aSBarry 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
54785afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
548d52bd9f3SBarry Smith 
549d52bd9f3SBarry Smith         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
550d52bd9f3SBarry Smith         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
55185afcc9aSBarry Smith 
55285afcc9aSBarry Smith 
553e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
55447c6ae99SBarry Smith 
55547c6ae99SBarry Smith @*/
556e727c939SJed Brown PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
55747c6ae99SBarry Smith {
55847c6ae99SBarry Smith   PetscErrorCode ierr;
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith   PetscFunctionBegin;
561171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
562171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
56325296bd5SBarry Smith   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
56447c6ae99SBarry Smith   PetscFunctionReturn(0);
56547c6ae99SBarry Smith }
56647c6ae99SBarry Smith 
56747c6ae99SBarry Smith #undef __FUNCT__
568e727c939SJed Brown #define __FUNCT__ "DMCreateInjection"
56947c6ae99SBarry Smith /*@
570e727c939SJed Brown     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith     Collective on DM
57347c6ae99SBarry Smith 
57447c6ae99SBarry Smith     Input Parameter:
57547c6ae99SBarry Smith +   dm1 - the DM object
57647c6ae99SBarry Smith -   dm2 - the second, finer DM object
57747c6ae99SBarry Smith 
57847c6ae99SBarry Smith     Output Parameter:
57947c6ae99SBarry Smith .   ctx - the injection
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith     Level: developer
58247c6ae99SBarry Smith 
58385afcc9aSBarry 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
58485afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
58585afcc9aSBarry Smith 
586e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith @*/
589e727c939SJed Brown PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
59047c6ae99SBarry Smith {
59147c6ae99SBarry Smith   PetscErrorCode ierr;
59247c6ae99SBarry Smith 
59347c6ae99SBarry Smith   PetscFunctionBegin;
594171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
595171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
59647c6ae99SBarry Smith   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
59747c6ae99SBarry Smith   PetscFunctionReturn(0);
59847c6ae99SBarry Smith }
59947c6ae99SBarry Smith 
60047c6ae99SBarry Smith #undef __FUNCT__
601e727c939SJed Brown #define __FUNCT__ "DMCreateColoring"
602d1e2c406SBarry Smith /*@C
603e727c939SJed Brown     DMCreateColoring - Gets coloring for a DMDA or DMComposite
60447c6ae99SBarry Smith 
60547c6ae99SBarry Smith     Collective on DM
60647c6ae99SBarry Smith 
60747c6ae99SBarry Smith     Input Parameter:
60847c6ae99SBarry Smith +   dm - the DM object
60947c6ae99SBarry Smith .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
61047c6ae99SBarry Smith -   matype - either MATAIJ or MATBAIJ
61147c6ae99SBarry Smith 
61247c6ae99SBarry Smith     Output Parameter:
61347c6ae99SBarry Smith .   coloring - the coloring
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith     Level: developer
61647c6ae99SBarry Smith 
617e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
61847c6ae99SBarry Smith 
619aab9d709SJed Brown @*/
620e727c939SJed Brown PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
62147c6ae99SBarry Smith {
62247c6ae99SBarry Smith   PetscErrorCode ierr;
62347c6ae99SBarry Smith 
62447c6ae99SBarry Smith   PetscFunctionBegin;
625171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
62647c6ae99SBarry Smith   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
62747c6ae99SBarry Smith   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
62847c6ae99SBarry Smith   PetscFunctionReturn(0);
62947c6ae99SBarry Smith }
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith #undef __FUNCT__
632950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix"
63347c6ae99SBarry Smith /*@C
634950540a4SJed Brown     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
63547c6ae99SBarry Smith 
63647c6ae99SBarry Smith     Collective on DM
63747c6ae99SBarry Smith 
63847c6ae99SBarry Smith     Input Parameter:
63947c6ae99SBarry Smith +   dm - the DM object
64047c6ae99SBarry Smith -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
64194013140SBarry Smith             any type which inherits from one of these (such as MATAIJ)
64247c6ae99SBarry Smith 
64347c6ae99SBarry Smith     Output Parameter:
64447c6ae99SBarry Smith .   mat - the empty Jacobian
64547c6ae99SBarry Smith 
646073dac72SJed Brown     Level: beginner
64747c6ae99SBarry Smith 
64894013140SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
64994013140SBarry Smith        do not need to do it yourself.
65094013140SBarry Smith 
65194013140SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
652aa219208SBarry Smith        the nonzero pattern call DMDASetMatPreallocateOnly()
65394013140SBarry Smith 
65494013140SBarry 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
65594013140SBarry Smith        internally by PETSc.
65694013140SBarry Smith 
65794013140SBarry Smith        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
658aa219208SBarry Smith        the indices for the global numbering for DMDAs which is complicated.
65994013140SBarry Smith 
660e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
66147c6ae99SBarry Smith 
662aab9d709SJed Brown @*/
663950540a4SJed Brown PetscErrorCode  DMCreateMatrix(DM dm,const MatType mtype,Mat *mat)
66447c6ae99SBarry Smith {
66547c6ae99SBarry Smith   PetscErrorCode ierr;
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith   PetscFunctionBegin;
668171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
669235683edSBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
670235683edSBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
671235683edSBarry Smith #endif
672c7b7c8a4SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
673c7b7c8a4SJed Brown   PetscValidPointer(mat,3);
674073dac72SJed Brown   if (dm->mattype) {
67525296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
676073dac72SJed Brown   } else {
67725296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
678c7b7c8a4SJed Brown   }
67947c6ae99SBarry Smith   PetscFunctionReturn(0);
68047c6ae99SBarry Smith }
68147c6ae99SBarry Smith 
68247c6ae99SBarry Smith #undef __FUNCT__
683732e2eb9SMatthew G Knepley #define __FUNCT__ "DMSetMatrixPreallocateOnly"
684732e2eb9SMatthew G Knepley /*@
685950540a4SJed Brown   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
686732e2eb9SMatthew G Knepley     preallocated but the nonzero structure and zero values will not be set.
687732e2eb9SMatthew G Knepley 
688732e2eb9SMatthew G Knepley   Logically Collective on DMDA
689732e2eb9SMatthew G Knepley 
690732e2eb9SMatthew G Knepley   Input Parameter:
691732e2eb9SMatthew G Knepley + dm - the DM
692732e2eb9SMatthew G Knepley - only - PETSC_TRUE if only want preallocation
693732e2eb9SMatthew G Knepley 
694732e2eb9SMatthew G Knepley   Level: developer
695950540a4SJed Brown .seealso DMCreateMatrix()
696732e2eb9SMatthew G Knepley @*/
697732e2eb9SMatthew G Knepley PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
698732e2eb9SMatthew G Knepley {
699732e2eb9SMatthew G Knepley   PetscFunctionBegin;
700732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
701732e2eb9SMatthew G Knepley   dm->prealloc_only = only;
702732e2eb9SMatthew G Knepley   PetscFunctionReturn(0);
703732e2eb9SMatthew G Knepley }
704732e2eb9SMatthew G Knepley 
705732e2eb9SMatthew G Knepley #undef __FUNCT__
706a89ea682SMatthew G Knepley #define __FUNCT__ "DMGetWorkArray"
707a89ea682SMatthew G Knepley /*@C
708a89ea682SMatthew G Knepley   DMGetWorkArray - Gets a work array guaranteed to be at least the input size
709a89ea682SMatthew G Knepley 
710a89ea682SMatthew G Knepley   Not Collective
711a89ea682SMatthew G Knepley 
712a89ea682SMatthew G Knepley   Input Parameters:
713a89ea682SMatthew G Knepley + dm - the DM object
714a89ea682SMatthew G Knepley - size - The minium size
715a89ea682SMatthew G Knepley 
716a89ea682SMatthew G Knepley   Output Parameter:
717a89ea682SMatthew G Knepley . array - the work array
718a89ea682SMatthew G Knepley 
719a89ea682SMatthew G Knepley   Level: developer
720a89ea682SMatthew G Knepley 
721a89ea682SMatthew G Knepley .seealso DMDestroy(), DMCreate()
722a89ea682SMatthew G Knepley @*/
723a89ea682SMatthew G Knepley PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array)
724a89ea682SMatthew G Knepley {
725a89ea682SMatthew G Knepley   PetscErrorCode ierr;
726a89ea682SMatthew G Knepley 
727a89ea682SMatthew G Knepley   PetscFunctionBegin;
728a89ea682SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
729a89ea682SMatthew G Knepley   PetscValidPointer(array,3);
730a89ea682SMatthew G Knepley   if (size > dm->workSize) {
731a89ea682SMatthew G Knepley     dm->workSize = size;
732a89ea682SMatthew G Knepley     ierr = PetscFree(dm->workArray);CHKERRQ(ierr);
733a89ea682SMatthew G Knepley     ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr);
734a89ea682SMatthew G Knepley   }
735a89ea682SMatthew G Knepley   *array = dm->workArray;
736a89ea682SMatthew G Knepley   PetscFunctionReturn(0);
737a89ea682SMatthew G Knepley }
738a89ea682SMatthew G Knepley 
7394d343eeaSMatthew G Knepley #undef __FUNCT__
740e7c4fc90SDmitry Karpeev #define __FUNCT__ "DMCreateDecompositionDM"
741e7c4fc90SDmitry Karpeev /*@C
742e7c4fc90SDmitry Karpeev   DMCreateDecompositionDM - creates a DM that defines a decomposition of the original DM.
743e7c4fc90SDmitry Karpeev 
744e7c4fc90SDmitry Karpeev   Not Collective
745e7c4fc90SDmitry Karpeev 
746e7c4fc90SDmitry Karpeev   Input Parameters:
747e7c4fc90SDmitry Karpeev + dm   - the DM object
748e7c4fc90SDmitry Karpeev - name - the name of the decomposition
749e7c4fc90SDmitry Karpeev 
750e7c4fc90SDmitry Karpeev   Output Parameter:
751e7c4fc90SDmitry Karpeev . ddm  - the decomposition DM (PETSC_NULL, if no such decomposition is known)
752e7c4fc90SDmitry Karpeev 
753e7c4fc90SDmitry Karpeev   Level: advanced
754e7c4fc90SDmitry Karpeev 
755e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMCreate(), DMCreateDecomposition()
756e7c4fc90SDmitry Karpeev @*/
757e7c4fc90SDmitry Karpeev PetscErrorCode DMCreateDecompositionDM(DM dm, const char* name, DM *ddm)
758e7c4fc90SDmitry Karpeev {
759e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
760e7c4fc90SDmitry Karpeev 
761e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
762e7c4fc90SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
763e7c4fc90SDmitry Karpeev   PetscValidCharPointer(name,2);
764e7c4fc90SDmitry Karpeev   PetscValidPointer(ddm,3);
765e7c4fc90SDmitry Karpeev   if(!dm->ops->createdecompositiondm) {
766e7c4fc90SDmitry Karpeev     *ddm = PETSC_NULL;
767e7c4fc90SDmitry Karpeev   }
768e7c4fc90SDmitry Karpeev   else {
769e7c4fc90SDmitry Karpeev     ierr = (*dm->ops->createdecompositiondm)(dm,name,ddm); CHKERRQ(ierr);
770e7c4fc90SDmitry Karpeev   }
771e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
772e7c4fc90SDmitry Karpeev }
773e7c4fc90SDmitry Karpeev 
774e7c4fc90SDmitry Karpeev #undef __FUNCT__
7754d343eeaSMatthew G Knepley #define __FUNCT__ "DMCreateFieldIS"
7764f3b5142SJed Brown /*@C
7774d343eeaSMatthew G Knepley   DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field
7784d343eeaSMatthew G Knepley 
7794d343eeaSMatthew G Knepley   Not collective
7804d343eeaSMatthew G Knepley 
7814d343eeaSMatthew G Knepley   Input Parameter:
7824d343eeaSMatthew G Knepley . dm - the DM object
7834d343eeaSMatthew G Knepley 
7844d343eeaSMatthew G Knepley   Output Parameters:
78521c9b008SJed Brown + numFields - The number of fields (or PETSC_NULL if not requested)
78621c9b008SJed Brown . names     - The name for each field (or PETSC_NULL if not requested)
78721c9b008SJed Brown - fields    - The global indices for each field (or PETSC_NULL if not requested)
7884d343eeaSMatthew G Knepley 
7894d343eeaSMatthew G Knepley   Level: intermediate
7904d343eeaSMatthew G Knepley 
79121c9b008SJed Brown   Notes:
79221c9b008SJed Brown   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
79321c9b008SJed Brown   PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with
79421c9b008SJed Brown   PetscFree().
79521c9b008SJed Brown 
7964d343eeaSMatthew G Knepley .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
7974d343eeaSMatthew G Knepley @*/
79821c9b008SJed Brown PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***names, IS **fields)
7994d343eeaSMatthew G Knepley {
8004d343eeaSMatthew G Knepley   PetscErrorCode ierr;
8014d343eeaSMatthew G Knepley 
8024d343eeaSMatthew G Knepley   PetscFunctionBegin;
8034d343eeaSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8044d343eeaSMatthew G Knepley   if (numFields) {PetscValidPointer(numFields,2);}
8054d343eeaSMatthew G Knepley   if (names) {PetscValidPointer(names,3);}
8064d343eeaSMatthew G Knepley   if (fields) {PetscValidPointer(fields,4);}
8074d343eeaSMatthew G Knepley   ierr = (*dm->ops->createfieldis)(dm, numFields, names, fields);CHKERRQ(ierr);
8084d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
8094d343eeaSMatthew G Knepley }
8104d343eeaSMatthew G Knepley 
8115fe1f584SPeter Brune 
812a89ea682SMatthew G Knepley #undef __FUNCT__
813e7c4fc90SDmitry Karpeev #define __FUNCT__ "DMCreateDecomposition"
814e7c4fc90SDmitry Karpeev /*@C
815e7c4fc90SDmitry Karpeev   DMCreateDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems:
816e7c4fc90SDmitry Karpeev                           each IS contains the global indices of the dofs of the corresponding subproblem.
817e7c4fc90SDmitry Karpeev                           The optional list of DMs define the DM for each subproblem.
818e7c4fc90SDmitry Karpeev                           Generalizes DMCreateFieldIS().
819e7c4fc90SDmitry Karpeev 
820e7c4fc90SDmitry Karpeev   Not collective
821e7c4fc90SDmitry Karpeev 
822e7c4fc90SDmitry Karpeev   Input Parameter:
823e7c4fc90SDmitry Karpeev . dm - the DM object
824e7c4fc90SDmitry Karpeev 
825e7c4fc90SDmitry Karpeev   Output Parameters:
826e7c4fc90SDmitry Karpeev + len       - The number of subproblems in the decomposition (or PETSC_NULL if not requested)
827e7c4fc90SDmitry Karpeev . namelist  - The name for each subproblem (or PETSC_NULL if not requested)
828e7c4fc90SDmitry Karpeev . islist    - The global indices for each subproblem (or PETSC_NULL if not requested)
829e7c4fc90SDmitry Karpeev - dmlist    - The DMs for each subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined)
830e7c4fc90SDmitry Karpeev 
831e7c4fc90SDmitry Karpeev   Level: intermediate
832e7c4fc90SDmitry Karpeev 
833e7c4fc90SDmitry Karpeev   Notes:
834e7c4fc90SDmitry Karpeev   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
835e7c4fc90SDmitry Karpeev   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
836e7c4fc90SDmitry Karpeev   and all of the arrays should be freed with PetscFree().
837e7c4fc90SDmitry Karpeev 
838e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
839e7c4fc90SDmitry Karpeev @*/
840e7c4fc90SDmitry Karpeev PetscErrorCode DMCreateDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
841e7c4fc90SDmitry Karpeev {
842e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
843e7c4fc90SDmitry Karpeev 
844e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
845e7c4fc90SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
846e7c4fc90SDmitry Karpeev   if (len) {PetscValidPointer(len,2);}
847e7c4fc90SDmitry Karpeev   if (namelist) {PetscValidPointer(namelist,3);}
848e7c4fc90SDmitry Karpeev   if (islist) {PetscValidPointer(islist,4);}
849e7c4fc90SDmitry Karpeev   if (dmlist) {PetscValidPointer(dmlist,5);}
850e7c4fc90SDmitry Karpeev   if(!dm->ops->createdecomposition) {
851e7c4fc90SDmitry Karpeev     ierr = (*dm->ops->createfieldis)(dm, len, namelist, islist);CHKERRQ(ierr);
852e7c4fc90SDmitry Karpeev     /* By default there are no DMs associated with subproblems. */
853e7c4fc90SDmitry Karpeev     if(dmlist) *dmlist = PETSC_NULL;
854e7c4fc90SDmitry Karpeev   }
855e7c4fc90SDmitry Karpeev   else {
856e7c4fc90SDmitry Karpeev     ierr = (*dm->ops->createdecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr);
857e7c4fc90SDmitry Karpeev   }
858e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
859e7c4fc90SDmitry Karpeev }
860e7c4fc90SDmitry Karpeev 
861e7c4fc90SDmitry Karpeev 
862e7c4fc90SDmitry Karpeev #undef __FUNCT__
86347c6ae99SBarry Smith #define __FUNCT__ "DMRefine"
86447c6ae99SBarry Smith /*@
86547c6ae99SBarry Smith   DMRefine - Refines a DM object
86647c6ae99SBarry Smith 
86747c6ae99SBarry Smith   Collective on DM
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith   Input Parameter:
87047c6ae99SBarry Smith + dm   - the DM object
87191d95f02SJed Brown - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
87247c6ae99SBarry Smith 
87347c6ae99SBarry Smith   Output Parameter:
874ae0a1c52SMatthew G Knepley . dmf - the refined DM, or PETSC_NULL
875ae0a1c52SMatthew G Knepley 
876ae0a1c52SMatthew G Knepley   Note: If no refinement was done, the return value is PETSC_NULL
87747c6ae99SBarry Smith 
87847c6ae99SBarry Smith   Level: developer
87947c6ae99SBarry Smith 
880e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
88147c6ae99SBarry Smith @*/
8827087cfbeSBarry Smith PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
88347c6ae99SBarry Smith {
88447c6ae99SBarry Smith   PetscErrorCode ierr;
88547c6ae99SBarry Smith 
88647c6ae99SBarry Smith   PetscFunctionBegin;
887732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
88847c6ae99SBarry Smith   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
8894057135bSMatthew G Knepley   if (*dmf) {
89043842a1eSJed Brown     (*dmf)->ops->creatematrix = dm->ops->creatematrix;
891644e2e5bSBarry Smith     (*dmf)->ops->initialguess = dm->ops->initialguess;
892644e2e5bSBarry Smith     (*dmf)->ops->function     = dm->ops->function;
893644e2e5bSBarry Smith     (*dmf)->ops->functionj    = dm->ops->functionj;
894644e2e5bSBarry Smith     if (dm->ops->jacobian != DMComputeJacobianDefault) {
895644e2e5bSBarry Smith       (*dmf)->ops->jacobian     = dm->ops->jacobian;
896644e2e5bSBarry Smith     }
8978cd211a4SJed Brown     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
898644e2e5bSBarry Smith     (*dmf)->ctx     = dm->ctx;
899656b349aSBarry Smith     (*dmf)->levelup = dm->levelup + 1;
9004057135bSMatthew G Knepley   }
90147c6ae99SBarry Smith   PetscFunctionReturn(0);
90247c6ae99SBarry Smith }
90347c6ae99SBarry Smith 
90447c6ae99SBarry Smith #undef __FUNCT__
905eb3f98d2SBarry Smith #define __FUNCT__ "DMGetRefineLevel"
906eb3f98d2SBarry Smith /*@
907eb3f98d2SBarry Smith     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
908eb3f98d2SBarry Smith 
909eb3f98d2SBarry Smith     Not Collective
910eb3f98d2SBarry Smith 
911eb3f98d2SBarry Smith     Input Parameter:
912eb3f98d2SBarry Smith .   dm - the DM object
913eb3f98d2SBarry Smith 
914eb3f98d2SBarry Smith     Output Parameter:
915eb3f98d2SBarry Smith .   level - number of refinements
916eb3f98d2SBarry Smith 
917eb3f98d2SBarry Smith     Level: developer
918eb3f98d2SBarry Smith 
9196a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
920eb3f98d2SBarry Smith 
921eb3f98d2SBarry Smith @*/
922eb3f98d2SBarry Smith PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
923eb3f98d2SBarry Smith {
924eb3f98d2SBarry Smith   PetscFunctionBegin;
925eb3f98d2SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
926eb3f98d2SBarry Smith   *level = dm->levelup;
927eb3f98d2SBarry Smith   PetscFunctionReturn(0);
928eb3f98d2SBarry Smith }
929eb3f98d2SBarry Smith 
930eb3f98d2SBarry Smith #undef __FUNCT__
93147c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin"
93247c6ae99SBarry Smith /*@
93347c6ae99SBarry Smith     DMGlobalToLocalBegin - Begins updating local vectors from global vector
93447c6ae99SBarry Smith 
93547c6ae99SBarry Smith     Neighbor-wise Collective on DM
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith     Input Parameters:
93847c6ae99SBarry Smith +   dm - the DM object
93947c6ae99SBarry Smith .   g - the global vector
94047c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
94147c6ae99SBarry Smith -   l - the local vector
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith 
94447c6ae99SBarry Smith     Level: beginner
94547c6ae99SBarry Smith 
946e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith @*/
9497087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
95047c6ae99SBarry Smith {
95147c6ae99SBarry Smith   PetscErrorCode ierr;
95247c6ae99SBarry Smith 
95347c6ae99SBarry Smith   PetscFunctionBegin;
954171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
955843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
95647c6ae99SBarry Smith   PetscFunctionReturn(0);
95747c6ae99SBarry Smith }
95847c6ae99SBarry Smith 
95947c6ae99SBarry Smith #undef __FUNCT__
96047c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd"
96147c6ae99SBarry Smith /*@
96247c6ae99SBarry Smith     DMGlobalToLocalEnd - Ends updating local vectors from global vector
96347c6ae99SBarry Smith 
96447c6ae99SBarry Smith     Neighbor-wise Collective on DM
96547c6ae99SBarry Smith 
96647c6ae99SBarry Smith     Input Parameters:
96747c6ae99SBarry Smith +   dm - the DM object
96847c6ae99SBarry Smith .   g - the global vector
96947c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
97047c6ae99SBarry Smith -   l - the local vector
97147c6ae99SBarry Smith 
97247c6ae99SBarry Smith 
97347c6ae99SBarry Smith     Level: beginner
97447c6ae99SBarry Smith 
975e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith @*/
9787087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
97947c6ae99SBarry Smith {
98047c6ae99SBarry Smith   PetscErrorCode ierr;
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith   PetscFunctionBegin;
983171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
984843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
98547c6ae99SBarry Smith   PetscFunctionReturn(0);
98647c6ae99SBarry Smith }
98747c6ae99SBarry Smith 
98847c6ae99SBarry Smith #undef __FUNCT__
9899a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalBegin"
99047c6ae99SBarry Smith /*@
9919a42bb27SBarry Smith     DMLocalToGlobalBegin - updates global vectors from local vectors
9929a42bb27SBarry Smith 
9939a42bb27SBarry Smith     Neighbor-wise Collective on DM
9949a42bb27SBarry Smith 
9959a42bb27SBarry Smith     Input Parameters:
9969a42bb27SBarry Smith +   dm - the DM object
997f6813fd5SJed Brown .   l - the local vector
9989a42bb27SBarry 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
9999a42bb27SBarry Smith            base point.
1000f6813fd5SJed Brown - - the global vector
10019a42bb27SBarry Smith 
10029a42bb27SBarry 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
10039a42bb27SBarry 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
10049a42bb27SBarry Smith            global array to the final global array with VecAXPY().
10059a42bb27SBarry Smith 
10069a42bb27SBarry Smith     Level: beginner
10079a42bb27SBarry Smith 
1008e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
10099a42bb27SBarry Smith 
10109a42bb27SBarry Smith @*/
10117087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
10129a42bb27SBarry Smith {
10139a42bb27SBarry Smith   PetscErrorCode ierr;
10149a42bb27SBarry Smith 
10159a42bb27SBarry Smith   PetscFunctionBegin;
1016171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1017843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
10189a42bb27SBarry Smith   PetscFunctionReturn(0);
10199a42bb27SBarry Smith }
10209a42bb27SBarry Smith 
10219a42bb27SBarry Smith #undef __FUNCT__
10229a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalEnd"
10239a42bb27SBarry Smith /*@
10249a42bb27SBarry Smith     DMLocalToGlobalEnd - updates global vectors from local vectors
102547c6ae99SBarry Smith 
102647c6ae99SBarry Smith     Neighbor-wise Collective on DM
102747c6ae99SBarry Smith 
102847c6ae99SBarry Smith     Input Parameters:
102947c6ae99SBarry Smith +   dm - the DM object
1030f6813fd5SJed Brown .   l - the local vector
103147c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
1032f6813fd5SJed Brown -   g - the global vector
103347c6ae99SBarry Smith 
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith     Level: beginner
103647c6ae99SBarry Smith 
1037e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
103847c6ae99SBarry Smith 
103947c6ae99SBarry Smith @*/
10407087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
104147c6ae99SBarry Smith {
104247c6ae99SBarry Smith   PetscErrorCode ierr;
104347c6ae99SBarry Smith 
104447c6ae99SBarry Smith   PetscFunctionBegin;
1045171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1046843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
104747c6ae99SBarry Smith   PetscFunctionReturn(0);
104847c6ae99SBarry Smith }
104947c6ae99SBarry Smith 
105047c6ae99SBarry Smith #undef __FUNCT__
105147c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobianDefault"
105247c6ae99SBarry Smith /*@
105347c6ae99SBarry Smith     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
105447c6ae99SBarry Smith 
105547c6ae99SBarry Smith     Collective on DM
105647c6ae99SBarry Smith 
105747c6ae99SBarry Smith     Input Parameter:
105847c6ae99SBarry Smith +   dm - the DM object
105947c6ae99SBarry Smith .   x - location to compute Jacobian at; may be ignored for linear problems
106047c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
106147c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
106247c6ae99SBarry Smith 
106347c6ae99SBarry Smith     Level: developer
106447c6ae99SBarry Smith 
1065e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
106647c6ae99SBarry Smith          DMSetFunction()
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith @*/
10697087cfbeSBarry Smith PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
107047c6ae99SBarry Smith {
107147c6ae99SBarry Smith   PetscErrorCode ierr;
1072171400e9SBarry Smith 
107347c6ae99SBarry Smith   PetscFunctionBegin;
1074171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
107547c6ae99SBarry Smith   *stflag = SAME_NONZERO_PATTERN;
107647c6ae99SBarry Smith   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
107747c6ae99SBarry Smith   if (A != B) {
107847c6ae99SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107947c6ae99SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108047c6ae99SBarry Smith   }
108147c6ae99SBarry Smith   PetscFunctionReturn(0);
108247c6ae99SBarry Smith }
108347c6ae99SBarry Smith 
108447c6ae99SBarry Smith #undef __FUNCT__
108547c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen"
108647c6ae99SBarry Smith /*@
108747c6ae99SBarry Smith     DMCoarsen - Coarsens a DM object
108847c6ae99SBarry Smith 
108947c6ae99SBarry Smith     Collective on DM
109047c6ae99SBarry Smith 
109147c6ae99SBarry Smith     Input Parameter:
109247c6ae99SBarry Smith +   dm - the DM object
109391d95f02SJed Brown -   comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
109447c6ae99SBarry Smith 
109547c6ae99SBarry Smith     Output Parameter:
109647c6ae99SBarry Smith .   dmc - the coarsened DM
109747c6ae99SBarry Smith 
109847c6ae99SBarry Smith     Level: developer
109947c6ae99SBarry Smith 
1100e727c939SJed Brown .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
110147c6ae99SBarry Smith 
110247c6ae99SBarry Smith @*/
11037087cfbeSBarry Smith PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
110447c6ae99SBarry Smith {
110547c6ae99SBarry Smith   PetscErrorCode ierr;
1106b17ce1afSJed Brown   DMCoarsenHookLink link;
110747c6ae99SBarry Smith 
110847c6ae99SBarry Smith   PetscFunctionBegin;
1109171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
111047c6ae99SBarry Smith   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
111143842a1eSJed Brown   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
111247c6ae99SBarry Smith   (*dmc)->ops->initialguess = dm->ops->initialguess;
111347c6ae99SBarry Smith   (*dmc)->ops->function     = dm->ops->function;
111447c6ae99SBarry Smith   (*dmc)->ops->functionj    = dm->ops->functionj;
111547c6ae99SBarry Smith   if (dm->ops->jacobian != DMComputeJacobianDefault) {
111647c6ae99SBarry Smith     (*dmc)->ops->jacobian     = dm->ops->jacobian;
111747c6ae99SBarry Smith   }
11188cd211a4SJed Brown   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1119644e2e5bSBarry Smith   (*dmc)->ctx       = dm->ctx;
1120656b349aSBarry Smith   (*dmc)->leveldown = dm->leveldown + 1;
1121b17ce1afSJed Brown   for (link=dm->coarsenhook; link; link=link->next) {
1122b17ce1afSJed Brown     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1123b17ce1afSJed Brown   }
1124b17ce1afSJed Brown   PetscFunctionReturn(0);
1125b17ce1afSJed Brown }
1126b17ce1afSJed Brown 
1127b17ce1afSJed Brown #undef __FUNCT__
1128b17ce1afSJed Brown #define __FUNCT__ "DMCoarsenHookAdd"
1129b17ce1afSJed Brown /*@
1130b17ce1afSJed Brown    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1131b17ce1afSJed Brown 
1132b17ce1afSJed Brown    Logically Collective
1133b17ce1afSJed Brown 
1134b17ce1afSJed Brown    Input Arguments:
1135b17ce1afSJed Brown +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1136b17ce1afSJed Brown .  coarsenhook - function to run when setting up a coarser level
1137b17ce1afSJed Brown .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
1138b17ce1afSJed Brown -  ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL)
1139b17ce1afSJed Brown 
1140b17ce1afSJed Brown    Calling sequence of coarsenhook:
1141b17ce1afSJed Brown $    coarsenhook(DM fine,DM coarse,void *ctx);
1142b17ce1afSJed Brown 
1143b17ce1afSJed Brown +  fine - fine level DM
1144b17ce1afSJed Brown .  coarse - coarse level DM to restrict problem to
1145b17ce1afSJed Brown -  ctx - optional user-defined function context
1146b17ce1afSJed Brown 
1147b17ce1afSJed Brown    Calling sequence for restricthook:
1148b17ce1afSJed Brown $    restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx)
1149b17ce1afSJed Brown 
1150b17ce1afSJed Brown +  fine - fine level DM
1151b17ce1afSJed Brown .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1152b17ce1afSJed Brown .  inject - matrix restricting by applying the transpose of injection
1153b17ce1afSJed Brown .  coarse - coarse level DM to update
1154b17ce1afSJed Brown -  ctx - optional user-defined function context
1155b17ce1afSJed Brown 
1156b17ce1afSJed Brown    Level: advanced
1157b17ce1afSJed Brown 
1158b17ce1afSJed Brown    Notes:
1159b17ce1afSJed Brown    This function is only needed if auxiliary data needs to be set up on coarse grids.
1160b17ce1afSJed Brown 
1161b17ce1afSJed Brown    If this function is called multiple times, the hooks will be run in the order they are added.
1162b17ce1afSJed Brown 
1163b17ce1afSJed Brown    The
1164b17ce1afSJed Brown 
1165b17ce1afSJed Brown    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1166b17ce1afSJed Brown    extract the finest level information from its context (instead of from the SNES).
1167b17ce1afSJed Brown 
1168b17ce1afSJed Brown .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1169b17ce1afSJed Brown @*/
1170b17ce1afSJed Brown PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1171b17ce1afSJed Brown {
1172b17ce1afSJed Brown   PetscErrorCode ierr;
1173b17ce1afSJed Brown   DMCoarsenHookLink link,*p;
1174b17ce1afSJed Brown 
1175b17ce1afSJed Brown   PetscFunctionBegin;
1176b17ce1afSJed Brown   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
11776bfea28cSJed Brown   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1178b17ce1afSJed Brown   ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1179b17ce1afSJed Brown   link->coarsenhook = coarsenhook;
1180b17ce1afSJed Brown   link->restricthook = restricthook;
1181b17ce1afSJed Brown   link->ctx = ctx;
11826cab3a1bSJed Brown   link->next = PETSC_NULL;
1183b17ce1afSJed Brown   *p = link;
1184b17ce1afSJed Brown   PetscFunctionReturn(0);
1185b17ce1afSJed Brown }
1186b17ce1afSJed Brown 
1187b17ce1afSJed Brown #undef __FUNCT__
1188b17ce1afSJed Brown #define __FUNCT__ "DMRestrict"
1189b17ce1afSJed Brown /*@
1190b17ce1afSJed Brown    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1191b17ce1afSJed Brown 
1192b17ce1afSJed Brown    Collective if any hooks are
1193b17ce1afSJed Brown 
1194b17ce1afSJed Brown    Input Arguments:
1195b17ce1afSJed Brown +  fine - finer DM to use as a base
1196b17ce1afSJed Brown .  restrct - restriction matrix, apply using MatRestrict()
1197b17ce1afSJed Brown .  inject - injection matrix, also use MatRestrict()
1198b17ce1afSJed Brown -  coarse - coarer DM to update
1199b17ce1afSJed Brown 
1200b17ce1afSJed Brown    Level: developer
1201b17ce1afSJed Brown 
1202b17ce1afSJed Brown .seealso: DMCoarsenHookAdd(), MatRestrict()
1203b17ce1afSJed Brown @*/
1204b17ce1afSJed Brown PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1205b17ce1afSJed Brown {
1206b17ce1afSJed Brown   PetscErrorCode ierr;
1207b17ce1afSJed Brown   DMCoarsenHookLink link;
1208b17ce1afSJed Brown 
1209b17ce1afSJed Brown   PetscFunctionBegin;
1210b17ce1afSJed Brown   for (link=fine->coarsenhook; link; link=link->next) {
1211b17ce1afSJed Brown     if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);}
1212b17ce1afSJed Brown   }
121347c6ae99SBarry Smith   PetscFunctionReturn(0);
121447c6ae99SBarry Smith }
121547c6ae99SBarry Smith 
121647c6ae99SBarry Smith #undef __FUNCT__
12175fe1f584SPeter Brune #define __FUNCT__ "DMGetCoarsenLevel"
12185fe1f584SPeter Brune /*@
12196a7d9d85SPeter Brune     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
12205fe1f584SPeter Brune 
12215fe1f584SPeter Brune     Not Collective
12225fe1f584SPeter Brune 
12235fe1f584SPeter Brune     Input Parameter:
12245fe1f584SPeter Brune .   dm - the DM object
12255fe1f584SPeter Brune 
12265fe1f584SPeter Brune     Output Parameter:
12276a7d9d85SPeter Brune .   level - number of coarsenings
12285fe1f584SPeter Brune 
12295fe1f584SPeter Brune     Level: developer
12305fe1f584SPeter Brune 
12316a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
12325fe1f584SPeter Brune 
12335fe1f584SPeter Brune @*/
12345fe1f584SPeter Brune PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
12355fe1f584SPeter Brune {
12365fe1f584SPeter Brune   PetscFunctionBegin;
12375fe1f584SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12385fe1f584SPeter Brune   *level = dm->leveldown;
12395fe1f584SPeter Brune   PetscFunctionReturn(0);
12405fe1f584SPeter Brune }
12415fe1f584SPeter Brune 
12425fe1f584SPeter Brune 
12435fe1f584SPeter Brune 
12445fe1f584SPeter Brune #undef __FUNCT__
124547c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy"
124647c6ae99SBarry Smith /*@C
124747c6ae99SBarry Smith     DMRefineHierarchy - Refines a DM object, all levels at once
124847c6ae99SBarry Smith 
124947c6ae99SBarry Smith     Collective on DM
125047c6ae99SBarry Smith 
125147c6ae99SBarry Smith     Input Parameter:
125247c6ae99SBarry Smith +   dm - the DM object
125347c6ae99SBarry Smith -   nlevels - the number of levels of refinement
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith     Output Parameter:
125647c6ae99SBarry Smith .   dmf - the refined DM hierarchy
125747c6ae99SBarry Smith 
125847c6ae99SBarry Smith     Level: developer
125947c6ae99SBarry Smith 
1260e727c939SJed Brown .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
126147c6ae99SBarry Smith 
126247c6ae99SBarry Smith @*/
12637087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
126447c6ae99SBarry Smith {
126547c6ae99SBarry Smith   PetscErrorCode ierr;
126647c6ae99SBarry Smith 
126747c6ae99SBarry Smith   PetscFunctionBegin;
1268171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
126947c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
127047c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
127147c6ae99SBarry Smith   if (dm->ops->refinehierarchy) {
127247c6ae99SBarry Smith     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
127347c6ae99SBarry Smith   } else if (dm->ops->refine) {
127447c6ae99SBarry Smith     PetscInt i;
127547c6ae99SBarry Smith 
127647c6ae99SBarry Smith     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
127747c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
127847c6ae99SBarry Smith       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
127947c6ae99SBarry Smith     }
128047c6ae99SBarry Smith   } else {
128147c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
128247c6ae99SBarry Smith   }
128347c6ae99SBarry Smith   PetscFunctionReturn(0);
128447c6ae99SBarry Smith }
128547c6ae99SBarry Smith 
128647c6ae99SBarry Smith #undef __FUNCT__
128747c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy"
128847c6ae99SBarry Smith /*@C
128947c6ae99SBarry Smith     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
129047c6ae99SBarry Smith 
129147c6ae99SBarry Smith     Collective on DM
129247c6ae99SBarry Smith 
129347c6ae99SBarry Smith     Input Parameter:
129447c6ae99SBarry Smith +   dm - the DM object
129547c6ae99SBarry Smith -   nlevels - the number of levels of coarsening
129647c6ae99SBarry Smith 
129747c6ae99SBarry Smith     Output Parameter:
129847c6ae99SBarry Smith .   dmc - the coarsened DM hierarchy
129947c6ae99SBarry Smith 
130047c6ae99SBarry Smith     Level: developer
130147c6ae99SBarry Smith 
1302e727c939SJed Brown .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
130347c6ae99SBarry Smith 
130447c6ae99SBarry Smith @*/
13057087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
130647c6ae99SBarry Smith {
130747c6ae99SBarry Smith   PetscErrorCode ierr;
130847c6ae99SBarry Smith 
130947c6ae99SBarry Smith   PetscFunctionBegin;
1310171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
131147c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
131247c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
131347c6ae99SBarry Smith   PetscValidPointer(dmc,3);
131447c6ae99SBarry Smith   if (dm->ops->coarsenhierarchy) {
131547c6ae99SBarry Smith     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
131647c6ae99SBarry Smith   } else if (dm->ops->coarsen) {
131747c6ae99SBarry Smith     PetscInt i;
131847c6ae99SBarry Smith 
131947c6ae99SBarry Smith     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
132047c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
132147c6ae99SBarry Smith       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
132247c6ae99SBarry Smith     }
132347c6ae99SBarry Smith   } else {
132447c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
132547c6ae99SBarry Smith   }
132647c6ae99SBarry Smith   PetscFunctionReturn(0);
132747c6ae99SBarry Smith }
132847c6ae99SBarry Smith 
132947c6ae99SBarry Smith #undef __FUNCT__
1330e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates"
133147c6ae99SBarry Smith /*@
1332e727c939SJed Brown    DMCreateAggregates - Gets the aggregates that map between
133347c6ae99SBarry Smith    grids associated with two DMs.
133447c6ae99SBarry Smith 
133547c6ae99SBarry Smith    Collective on DM
133647c6ae99SBarry Smith 
133747c6ae99SBarry Smith    Input Parameters:
133847c6ae99SBarry Smith +  dmc - the coarse grid DM
133947c6ae99SBarry Smith -  dmf - the fine grid DM
134047c6ae99SBarry Smith 
134147c6ae99SBarry Smith    Output Parameters:
134247c6ae99SBarry Smith .  rest - the restriction matrix (transpose of the projection matrix)
134347c6ae99SBarry Smith 
134447c6ae99SBarry Smith    Level: intermediate
134547c6ae99SBarry Smith 
134647c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid
134747c6ae99SBarry Smith 
1348e727c939SJed Brown .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
134947c6ae99SBarry Smith @*/
1350e727c939SJed Brown PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
135147c6ae99SBarry Smith {
135247c6ae99SBarry Smith   PetscErrorCode ierr;
135347c6ae99SBarry Smith 
135447c6ae99SBarry Smith   PetscFunctionBegin;
1355171400e9SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1356171400e9SBarry Smith   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
135747c6ae99SBarry Smith   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
135847c6ae99SBarry Smith   PetscFunctionReturn(0);
135947c6ae99SBarry Smith }
136047c6ae99SBarry Smith 
136147c6ae99SBarry Smith #undef __FUNCT__
13621a266240SBarry Smith #define __FUNCT__ "DMSetApplicationContextDestroy"
13631a266240SBarry Smith /*@C
13641a266240SBarry Smith     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
13651a266240SBarry Smith 
13661a266240SBarry Smith     Not Collective
13671a266240SBarry Smith 
13681a266240SBarry Smith     Input Parameters:
13691a266240SBarry Smith +   dm - the DM object
13701a266240SBarry Smith -   destroy - the destroy function
13711a266240SBarry Smith 
13721a266240SBarry Smith     Level: intermediate
13731a266240SBarry Smith 
1374e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
13751a266240SBarry Smith 
1376f07f9ceaSJed Brown @*/
13771a266240SBarry Smith PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
13781a266240SBarry Smith {
13791a266240SBarry Smith   PetscFunctionBegin;
1380171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13811a266240SBarry Smith   dm->ctxdestroy = destroy;
13821a266240SBarry Smith   PetscFunctionReturn(0);
13831a266240SBarry Smith }
13841a266240SBarry Smith 
13851a266240SBarry Smith #undef __FUNCT__
13861b2093e4SBarry Smith #define __FUNCT__ "DMSetApplicationContext"
1387b07ff414SBarry Smith /*@
13881b2093e4SBarry Smith     DMSetApplicationContext - Set a user context into a DM object
138947c6ae99SBarry Smith 
139047c6ae99SBarry Smith     Not Collective
139147c6ae99SBarry Smith 
139247c6ae99SBarry Smith     Input Parameters:
139347c6ae99SBarry Smith +   dm - the DM object
139447c6ae99SBarry Smith -   ctx - the user context
139547c6ae99SBarry Smith 
139647c6ae99SBarry Smith     Level: intermediate
139747c6ae99SBarry Smith 
1398e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
139947c6ae99SBarry Smith 
140047c6ae99SBarry Smith @*/
14011b2093e4SBarry Smith PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
140247c6ae99SBarry Smith {
140347c6ae99SBarry Smith   PetscFunctionBegin;
1404171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
140547c6ae99SBarry Smith   dm->ctx = ctx;
140647c6ae99SBarry Smith   PetscFunctionReturn(0);
140747c6ae99SBarry Smith }
140847c6ae99SBarry Smith 
140947c6ae99SBarry Smith #undef __FUNCT__
14101b2093e4SBarry Smith #define __FUNCT__ "DMGetApplicationContext"
141147c6ae99SBarry Smith /*@
14121b2093e4SBarry Smith     DMGetApplicationContext - Gets a user context from a DM object
141347c6ae99SBarry Smith 
141447c6ae99SBarry Smith     Not Collective
141547c6ae99SBarry Smith 
141647c6ae99SBarry Smith     Input Parameter:
141747c6ae99SBarry Smith .   dm - the DM object
141847c6ae99SBarry Smith 
141947c6ae99SBarry Smith     Output Parameter:
142047c6ae99SBarry Smith .   ctx - the user context
142147c6ae99SBarry Smith 
142247c6ae99SBarry Smith     Level: intermediate
142347c6ae99SBarry Smith 
1424e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
142547c6ae99SBarry Smith 
142647c6ae99SBarry Smith @*/
14271b2093e4SBarry Smith PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
142847c6ae99SBarry Smith {
142947c6ae99SBarry Smith   PetscFunctionBegin;
1430171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
14311b2093e4SBarry Smith   *(void**)ctx = dm->ctx;
143247c6ae99SBarry Smith   PetscFunctionReturn(0);
143347c6ae99SBarry Smith }
143447c6ae99SBarry Smith 
143547c6ae99SBarry Smith #undef __FUNCT__
143647c6ae99SBarry Smith #define __FUNCT__ "DMSetInitialGuess"
14377e833e3aSBarry Smith /*@C
143847c6ae99SBarry Smith     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
143947c6ae99SBarry Smith 
144047c6ae99SBarry Smith     Logically Collective on DM
144147c6ae99SBarry Smith 
144247c6ae99SBarry Smith     Input Parameter:
144347c6ae99SBarry Smith +   dm - the DM object to destroy
144447c6ae99SBarry Smith -   f - the function to compute the initial guess
144547c6ae99SBarry Smith 
144647c6ae99SBarry Smith     Level: intermediate
144747c6ae99SBarry Smith 
1448e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
144947c6ae99SBarry Smith 
1450f07f9ceaSJed Brown @*/
14517087cfbeSBarry Smith PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
145247c6ae99SBarry Smith {
145347c6ae99SBarry Smith   PetscFunctionBegin;
1454171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
145547c6ae99SBarry Smith   dm->ops->initialguess = f;
145647c6ae99SBarry Smith   PetscFunctionReturn(0);
145747c6ae99SBarry Smith }
145847c6ae99SBarry Smith 
145947c6ae99SBarry Smith #undef __FUNCT__
146047c6ae99SBarry Smith #define __FUNCT__ "DMSetFunction"
14617e833e3aSBarry Smith /*@C
146247c6ae99SBarry Smith     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
146347c6ae99SBarry Smith 
146447c6ae99SBarry Smith     Logically Collective on DM
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith     Input Parameter:
146747c6ae99SBarry Smith +   dm - the DM object
146847c6ae99SBarry Smith -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
146947c6ae99SBarry Smith 
147047c6ae99SBarry Smith     Level: intermediate
147147c6ae99SBarry Smith 
147247c6ae99SBarry Smith     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
147347c6ae99SBarry Smith            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
147447c6ae99SBarry Smith 
1475e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
147647c6ae99SBarry Smith          DMSetJacobian()
147747c6ae99SBarry Smith 
1478f07f9ceaSJed Brown @*/
14797087cfbeSBarry Smith PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
148047c6ae99SBarry Smith {
148147c6ae99SBarry Smith   PetscFunctionBegin;
1482171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
148347c6ae99SBarry Smith   dm->ops->function = f;
148447c6ae99SBarry Smith   if (f) {
148547c6ae99SBarry Smith     dm->ops->functionj = f;
148647c6ae99SBarry Smith   }
148747c6ae99SBarry Smith   PetscFunctionReturn(0);
148847c6ae99SBarry Smith }
148947c6ae99SBarry Smith 
149047c6ae99SBarry Smith #undef __FUNCT__
149147c6ae99SBarry Smith #define __FUNCT__ "DMSetJacobian"
14927e833e3aSBarry Smith /*@C
149347c6ae99SBarry Smith     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
149447c6ae99SBarry Smith 
149547c6ae99SBarry Smith     Logically Collective on DM
149647c6ae99SBarry Smith 
149747c6ae99SBarry Smith     Input Parameter:
149847c6ae99SBarry Smith +   dm - the DM object to destroy
149947c6ae99SBarry Smith -   f - the function to compute the matrix entries
150047c6ae99SBarry Smith 
150147c6ae99SBarry Smith     Level: intermediate
150247c6ae99SBarry Smith 
1503e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
150447c6ae99SBarry Smith          DMSetFunction()
150547c6ae99SBarry Smith 
1506f07f9ceaSJed Brown @*/
15077087cfbeSBarry Smith PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
150847c6ae99SBarry Smith {
150947c6ae99SBarry Smith   PetscFunctionBegin;
1510171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
151147c6ae99SBarry Smith   dm->ops->jacobian = f;
151247c6ae99SBarry Smith   PetscFunctionReturn(0);
151347c6ae99SBarry Smith }
151447c6ae99SBarry Smith 
151547c6ae99SBarry Smith #undef __FUNCT__
151608da532bSDmitry Karpeev #define __FUNCT__ "DMSetVariableBounds"
151708da532bSDmitry Karpeev /*@C
151808da532bSDmitry Karpeev     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
151908da532bSDmitry Karpeev 
152008da532bSDmitry Karpeev     Logically Collective on DM
152108da532bSDmitry Karpeev 
152208da532bSDmitry Karpeev     Input Parameter:
152308da532bSDmitry Karpeev +   dm - the DM object
152408da532bSDmitry Karpeev -   f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set)
152508da532bSDmitry Karpeev 
152608da532bSDmitry Karpeev     Level: intermediate
152708da532bSDmitry Karpeev 
1528e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
152908da532bSDmitry Karpeev          DMSetJacobian()
153008da532bSDmitry Karpeev 
153108da532bSDmitry Karpeev @*/
153208da532bSDmitry Karpeev PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
153308da532bSDmitry Karpeev {
153408da532bSDmitry Karpeev   PetscFunctionBegin;
153508da532bSDmitry Karpeev   dm->ops->computevariablebounds = f;
153608da532bSDmitry Karpeev   PetscFunctionReturn(0);
153708da532bSDmitry Karpeev }
153808da532bSDmitry Karpeev 
153908da532bSDmitry Karpeev #undef __FUNCT__
154008da532bSDmitry Karpeev #define __FUNCT__ "DMHasVariableBounds"
154108da532bSDmitry Karpeev /*@
154208da532bSDmitry Karpeev     DMHasVariableBounds - does the DM object have a variable bounds function?
154308da532bSDmitry Karpeev 
154408da532bSDmitry Karpeev     Not Collective
154508da532bSDmitry Karpeev 
154608da532bSDmitry Karpeev     Input Parameter:
154708da532bSDmitry Karpeev .   dm - the DM object to destroy
154808da532bSDmitry Karpeev 
154908da532bSDmitry Karpeev     Output Parameter:
155008da532bSDmitry Karpeev .   flg - PETSC_TRUE if the variable bounds function exists
155108da532bSDmitry Karpeev 
155208da532bSDmitry Karpeev     Level: developer
155308da532bSDmitry Karpeev 
1554e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
155508da532bSDmitry Karpeev 
155608da532bSDmitry Karpeev @*/
155708da532bSDmitry Karpeev PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
155808da532bSDmitry Karpeev {
155908da532bSDmitry Karpeev   PetscFunctionBegin;
156008da532bSDmitry Karpeev   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
156108da532bSDmitry Karpeev   PetscFunctionReturn(0);
156208da532bSDmitry Karpeev }
156308da532bSDmitry Karpeev 
156408da532bSDmitry Karpeev #undef __FUNCT__
156508da532bSDmitry Karpeev #define __FUNCT__ "DMComputeVariableBounds"
156608da532bSDmitry Karpeev /*@C
156708da532bSDmitry Karpeev     DMComputeVariableBounds - compute variable bounds used by SNESVI.
156808da532bSDmitry Karpeev 
156908da532bSDmitry Karpeev     Logically Collective on DM
157008da532bSDmitry Karpeev 
157108da532bSDmitry Karpeev     Input Parameters:
157208da532bSDmitry Karpeev +   dm - the DM object to destroy
157308da532bSDmitry Karpeev -   x  - current solution at which the bounds are computed
157408da532bSDmitry Karpeev 
157508da532bSDmitry Karpeev     Output parameters:
157608da532bSDmitry Karpeev +   xl - lower bound
157708da532bSDmitry Karpeev -   xu - upper bound
157808da532bSDmitry Karpeev 
157908da532bSDmitry Karpeev     Level: intermediate
158008da532bSDmitry Karpeev 
1581e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
158208da532bSDmitry Karpeev          DMSetFunction(), DMSetVariableBounds()
158308da532bSDmitry Karpeev 
158408da532bSDmitry Karpeev @*/
158508da532bSDmitry Karpeev PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
158608da532bSDmitry Karpeev {
158708da532bSDmitry Karpeev   PetscErrorCode ierr;
158808da532bSDmitry Karpeev   PetscFunctionBegin;
158908da532bSDmitry Karpeev   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
159008da532bSDmitry Karpeev   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
159108da532bSDmitry Karpeev   if(dm->ops->computevariablebounds) {
159208da532bSDmitry Karpeev     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr);
159308da532bSDmitry Karpeev   }
159408da532bSDmitry Karpeev   else {
159508da532bSDmitry Karpeev     ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr);
159608da532bSDmitry Karpeev     ierr = VecSet(xu,SNES_VI_INF);  CHKERRQ(ierr);
159708da532bSDmitry Karpeev   }
159808da532bSDmitry Karpeev   PetscFunctionReturn(0);
159908da532bSDmitry Karpeev }
160008da532bSDmitry Karpeev 
160108da532bSDmitry Karpeev #undef __FUNCT__
160247c6ae99SBarry Smith #define __FUNCT__ "DMComputeInitialGuess"
160347c6ae99SBarry Smith /*@
160447c6ae99SBarry Smith     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
160547c6ae99SBarry Smith 
160647c6ae99SBarry Smith     Collective on DM
160747c6ae99SBarry Smith 
160847c6ae99SBarry Smith     Input Parameter:
160947c6ae99SBarry Smith +   dm - the DM object to destroy
161047c6ae99SBarry Smith -   x - the vector to hold the initial guess values
161147c6ae99SBarry Smith 
161247c6ae99SBarry Smith     Level: developer
161347c6ae99SBarry Smith 
1614e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
161547c6ae99SBarry Smith 
161647c6ae99SBarry Smith @*/
16177087cfbeSBarry Smith PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
161847c6ae99SBarry Smith {
161947c6ae99SBarry Smith   PetscErrorCode ierr;
162047c6ae99SBarry Smith   PetscFunctionBegin;
1621171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
162247c6ae99SBarry Smith   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
162347c6ae99SBarry Smith   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
162447c6ae99SBarry Smith   PetscFunctionReturn(0);
162547c6ae99SBarry Smith }
162647c6ae99SBarry Smith 
162747c6ae99SBarry Smith #undef __FUNCT__
162847c6ae99SBarry Smith #define __FUNCT__ "DMHasInitialGuess"
162947c6ae99SBarry Smith /*@
163047c6ae99SBarry Smith     DMHasInitialGuess - does the DM object have an initial guess function
163147c6ae99SBarry Smith 
163247c6ae99SBarry Smith     Not Collective
163347c6ae99SBarry Smith 
163447c6ae99SBarry Smith     Input Parameter:
163547c6ae99SBarry Smith .   dm - the DM object to destroy
163647c6ae99SBarry Smith 
163747c6ae99SBarry Smith     Output Parameter:
163847c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
163947c6ae99SBarry Smith 
164047c6ae99SBarry Smith     Level: developer
164147c6ae99SBarry Smith 
1642e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
164347c6ae99SBarry Smith 
164447c6ae99SBarry Smith @*/
16457087cfbeSBarry Smith PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
164647c6ae99SBarry Smith {
164747c6ae99SBarry Smith   PetscFunctionBegin;
1648171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
164947c6ae99SBarry Smith   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
165047c6ae99SBarry Smith   PetscFunctionReturn(0);
165147c6ae99SBarry Smith }
165247c6ae99SBarry Smith 
165347c6ae99SBarry Smith #undef __FUNCT__
165447c6ae99SBarry Smith #define __FUNCT__ "DMHasFunction"
165547c6ae99SBarry Smith /*@
165647c6ae99SBarry Smith     DMHasFunction - does the DM object have a function
165747c6ae99SBarry Smith 
165847c6ae99SBarry Smith     Not Collective
165947c6ae99SBarry Smith 
166047c6ae99SBarry Smith     Input Parameter:
166147c6ae99SBarry Smith .   dm - the DM object to destroy
166247c6ae99SBarry Smith 
166347c6ae99SBarry Smith     Output Parameter:
166447c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
166547c6ae99SBarry Smith 
166647c6ae99SBarry Smith     Level: developer
166747c6ae99SBarry Smith 
1668e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
166947c6ae99SBarry Smith 
167047c6ae99SBarry Smith @*/
16717087cfbeSBarry Smith PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
167247c6ae99SBarry Smith {
167347c6ae99SBarry Smith   PetscFunctionBegin;
1674171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
167547c6ae99SBarry Smith   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
167647c6ae99SBarry Smith   PetscFunctionReturn(0);
167747c6ae99SBarry Smith }
167847c6ae99SBarry Smith 
167947c6ae99SBarry Smith #undef __FUNCT__
168047c6ae99SBarry Smith #define __FUNCT__ "DMHasJacobian"
168147c6ae99SBarry Smith /*@
168247c6ae99SBarry Smith     DMHasJacobian - does the DM object have a matrix function
168347c6ae99SBarry Smith 
168447c6ae99SBarry Smith     Not Collective
168547c6ae99SBarry Smith 
168647c6ae99SBarry Smith     Input Parameter:
168747c6ae99SBarry Smith .   dm - the DM object to destroy
168847c6ae99SBarry Smith 
168947c6ae99SBarry Smith     Output Parameter:
169047c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
169147c6ae99SBarry Smith 
169247c6ae99SBarry Smith     Level: developer
169347c6ae99SBarry Smith 
1694e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
169547c6ae99SBarry Smith 
169647c6ae99SBarry Smith @*/
16977087cfbeSBarry Smith PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
169847c6ae99SBarry Smith {
169947c6ae99SBarry Smith   PetscFunctionBegin;
1700171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
170147c6ae99SBarry Smith   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
170247c6ae99SBarry Smith   PetscFunctionReturn(0);
170347c6ae99SBarry Smith }
170447c6ae99SBarry Smith 
170547c6ae99SBarry Smith #undef  __FUNCT__
170608da532bSDmitry Karpeev #define __FUNCT__ "DMSetVec"
1707748fac09SDmitry Karpeev /*@C
170808da532bSDmitry Karpeev     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
170908da532bSDmitry Karpeev 
171008da532bSDmitry Karpeev     Collective on DM
171108da532bSDmitry Karpeev 
171208da532bSDmitry Karpeev     Input Parameter:
171308da532bSDmitry Karpeev +   dm - the DM object
1714e88d7f4bSDmitry Karpeev -   x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems.
171508da532bSDmitry Karpeev 
171608da532bSDmitry Karpeev     Level: developer
171708da532bSDmitry Karpeev 
1718e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
171908da532bSDmitry Karpeev          DMSetFunction(), DMSetJacobian(), DMSetVariableBounds()
172008da532bSDmitry Karpeev 
172108da532bSDmitry Karpeev @*/
172208da532bSDmitry Karpeev PetscErrorCode  DMSetVec(DM dm,Vec x)
172308da532bSDmitry Karpeev {
172408da532bSDmitry Karpeev   PetscErrorCode ierr;
172508da532bSDmitry Karpeev   PetscFunctionBegin;
172608da532bSDmitry Karpeev   if (x) {
172708da532bSDmitry Karpeev     if (!dm->x) {
172808da532bSDmitry Karpeev       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
172908da532bSDmitry Karpeev     }
173008da532bSDmitry Karpeev     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
173108da532bSDmitry Karpeev   }
173208da532bSDmitry Karpeev   else if(dm->x) {
173308da532bSDmitry Karpeev     ierr = VecDestroy(&dm->x);  CHKERRQ(ierr);
173408da532bSDmitry Karpeev   }
173508da532bSDmitry Karpeev   PetscFunctionReturn(0);
173608da532bSDmitry Karpeev }
173708da532bSDmitry Karpeev 
173808da532bSDmitry Karpeev 
173908da532bSDmitry Karpeev #undef __FUNCT__
174047c6ae99SBarry Smith #define __FUNCT__ "DMComputeFunction"
174147c6ae99SBarry Smith /*@
174247c6ae99SBarry Smith     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
174347c6ae99SBarry Smith 
174447c6ae99SBarry Smith     Collective on DM
174547c6ae99SBarry Smith 
174647c6ae99SBarry Smith     Input Parameter:
174747c6ae99SBarry Smith +   dm - the DM object to destroy
174847c6ae99SBarry Smith .   x - the location where the function is evaluationed, may be ignored for linear problems
174947c6ae99SBarry Smith -   b - the vector to hold the right hand side entries
175047c6ae99SBarry Smith 
175147c6ae99SBarry Smith     Level: developer
175247c6ae99SBarry Smith 
1753e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
175447c6ae99SBarry Smith          DMSetJacobian()
175547c6ae99SBarry Smith 
175647c6ae99SBarry Smith @*/
17577087cfbeSBarry Smith PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
175847c6ae99SBarry Smith {
175947c6ae99SBarry Smith   PetscErrorCode ierr;
176047c6ae99SBarry Smith   PetscFunctionBegin;
1761171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
176247c6ae99SBarry Smith   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1763644e2e5bSBarry Smith   PetscStackPush("DM user function");
176447c6ae99SBarry Smith   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1765644e2e5bSBarry Smith   PetscStackPop;
176647c6ae99SBarry Smith   PetscFunctionReturn(0);
176747c6ae99SBarry Smith }
176847c6ae99SBarry Smith 
176947c6ae99SBarry Smith 
177008da532bSDmitry Karpeev 
177147c6ae99SBarry Smith #undef __FUNCT__
177247c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobian"
177347c6ae99SBarry Smith /*@
177447c6ae99SBarry Smith     DMComputeJacobian - compute the matrix entries for the solver
177547c6ae99SBarry Smith 
177647c6ae99SBarry Smith     Collective on DM
177747c6ae99SBarry Smith 
177847c6ae99SBarry Smith     Input Parameter:
177947c6ae99SBarry Smith +   dm - the DM object
1780cab2e9ccSBarry Smith .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
178147c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
178247c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
178347c6ae99SBarry Smith 
178447c6ae99SBarry Smith     Level: developer
178547c6ae99SBarry Smith 
1786e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
178747c6ae99SBarry Smith          DMSetFunction()
178847c6ae99SBarry Smith 
178947c6ae99SBarry Smith @*/
17907087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
179147c6ae99SBarry Smith {
179247c6ae99SBarry Smith   PetscErrorCode ierr;
179347c6ae99SBarry Smith 
179447c6ae99SBarry Smith   PetscFunctionBegin;
1795171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
179647c6ae99SBarry Smith   if (!dm->ops->jacobian) {
179747c6ae99SBarry Smith     ISColoring     coloring;
179847c6ae99SBarry Smith     MatFDColoring  fd;
179947c6ae99SBarry Smith 
1800171400e9SBarry Smith     ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
180147c6ae99SBarry Smith     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1802fcfd50ebSBarry Smith     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
180347c6ae99SBarry Smith     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
18040bdded8aSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
18050bdded8aSJed Brown     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
180671cd77b2SBarry Smith 
180747c6ae99SBarry Smith     dm->fd = fd;
180847c6ae99SBarry Smith     dm->ops->jacobian = DMComputeJacobianDefault;
18092533e041SBarry Smith 
181071cd77b2SBarry Smith     /* don't know why this is needed */
181171cd77b2SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
181247c6ae99SBarry Smith   }
181347c6ae99SBarry Smith   if (!x) x = dm->x;
181447c6ae99SBarry Smith   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1815cab2e9ccSBarry Smith 
181671cd77b2SBarry 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 */
1817649052a6SBarry Smith   if (x) {
1818cab2e9ccSBarry Smith     if (!dm->x) {
181971cd77b2SBarry Smith       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1820cab2e9ccSBarry Smith     }
1821cab2e9ccSBarry Smith     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1822649052a6SBarry Smith   }
1823a8248277SBarry Smith   if (A != B) {
1824a8248277SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1825a8248277SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1826a8248277SBarry Smith   }
182747c6ae99SBarry Smith   PetscFunctionReturn(0);
182847c6ae99SBarry Smith }
1829264ace61SBarry Smith 
1830cab2e9ccSBarry Smith 
1831264ace61SBarry Smith PetscFList DMList                       = PETSC_NULL;
1832264ace61SBarry Smith PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1833264ace61SBarry Smith 
1834264ace61SBarry Smith #undef __FUNCT__
1835264ace61SBarry Smith #define __FUNCT__ "DMSetType"
1836264ace61SBarry Smith /*@C
1837264ace61SBarry Smith   DMSetType - Builds a DM, for a particular DM implementation.
1838264ace61SBarry Smith 
1839264ace61SBarry Smith   Collective on DM
1840264ace61SBarry Smith 
1841264ace61SBarry Smith   Input Parameters:
1842264ace61SBarry Smith + dm     - The DM object
1843264ace61SBarry Smith - method - The name of the DM type
1844264ace61SBarry Smith 
1845264ace61SBarry Smith   Options Database Key:
1846264ace61SBarry Smith . -dm_type <type> - Sets the DM type; use -help for a list of available types
1847264ace61SBarry Smith 
1848264ace61SBarry Smith   Notes:
1849e1589f56SBarry Smith   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1850264ace61SBarry Smith 
1851264ace61SBarry Smith   Level: intermediate
1852264ace61SBarry Smith 
1853264ace61SBarry Smith .keywords: DM, set, type
1854264ace61SBarry Smith .seealso: DMGetType(), DMCreate()
1855264ace61SBarry Smith @*/
18567087cfbeSBarry Smith PetscErrorCode  DMSetType(DM dm, const DMType method)
1857264ace61SBarry Smith {
1858264ace61SBarry Smith   PetscErrorCode (*r)(DM);
1859264ace61SBarry Smith   PetscBool      match;
1860264ace61SBarry Smith   PetscErrorCode ierr;
1861264ace61SBarry Smith 
1862264ace61SBarry Smith   PetscFunctionBegin;
1863264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1864264ace61SBarry Smith   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1865264ace61SBarry Smith   if (match) PetscFunctionReturn(0);
1866264ace61SBarry Smith 
1867264ace61SBarry Smith   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
18684b91b6eaSBarry Smith   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1869264ace61SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1870264ace61SBarry Smith 
1871264ace61SBarry Smith   if (dm->ops->destroy) {
1872264ace61SBarry Smith     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1873b5c23020SJed Brown     dm->ops->destroy = PETSC_NULL;
1874264ace61SBarry Smith   }
1875264ace61SBarry Smith   ierr = (*r)(dm);CHKERRQ(ierr);
1876264ace61SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1877264ace61SBarry Smith   PetscFunctionReturn(0);
1878264ace61SBarry Smith }
1879264ace61SBarry Smith 
1880264ace61SBarry Smith #undef __FUNCT__
1881264ace61SBarry Smith #define __FUNCT__ "DMGetType"
1882264ace61SBarry Smith /*@C
1883264ace61SBarry Smith   DMGetType - Gets the DM type name (as a string) from the DM.
1884264ace61SBarry Smith 
1885264ace61SBarry Smith   Not Collective
1886264ace61SBarry Smith 
1887264ace61SBarry Smith   Input Parameter:
1888264ace61SBarry Smith . dm  - The DM
1889264ace61SBarry Smith 
1890264ace61SBarry Smith   Output Parameter:
1891264ace61SBarry Smith . type - The DM type name
1892264ace61SBarry Smith 
1893264ace61SBarry Smith   Level: intermediate
1894264ace61SBarry Smith 
1895264ace61SBarry Smith .keywords: DM, get, type, name
1896264ace61SBarry Smith .seealso: DMSetType(), DMCreate()
1897264ace61SBarry Smith @*/
18987087cfbeSBarry Smith PetscErrorCode  DMGetType(DM dm, const DMType *type)
1899264ace61SBarry Smith {
1900264ace61SBarry Smith   PetscErrorCode ierr;
1901264ace61SBarry Smith 
1902264ace61SBarry Smith   PetscFunctionBegin;
1903264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1904264ace61SBarry Smith   PetscValidCharPointer(type,2);
1905264ace61SBarry Smith   if (!DMRegisterAllCalled) {
1906264ace61SBarry Smith     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1907264ace61SBarry Smith   }
1908264ace61SBarry Smith   *type = ((PetscObject)dm)->type_name;
1909264ace61SBarry Smith   PetscFunctionReturn(0);
1910264ace61SBarry Smith }
1911264ace61SBarry Smith 
191267a56275SMatthew G Knepley #undef __FUNCT__
191367a56275SMatthew G Knepley #define __FUNCT__ "DMConvert"
191467a56275SMatthew G Knepley /*@C
191567a56275SMatthew G Knepley   DMConvert - Converts a DM to another DM, either of the same or different type.
191667a56275SMatthew G Knepley 
191767a56275SMatthew G Knepley   Collective on DM
191867a56275SMatthew G Knepley 
191967a56275SMatthew G Knepley   Input Parameters:
192067a56275SMatthew G Knepley + dm - the DM
192167a56275SMatthew G Knepley - newtype - new DM type (use "same" for the same type)
192267a56275SMatthew G Knepley 
192367a56275SMatthew G Knepley   Output Parameter:
192467a56275SMatthew G Knepley . M - pointer to new DM
192567a56275SMatthew G Knepley 
192667a56275SMatthew G Knepley   Notes:
192767a56275SMatthew G Knepley   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
192867a56275SMatthew G Knepley   the MPI communicator of the generated DM is always the same as the communicator
192967a56275SMatthew G Knepley   of the input DM.
193067a56275SMatthew G Knepley 
193167a56275SMatthew G Knepley   Level: intermediate
193267a56275SMatthew G Knepley 
193367a56275SMatthew G Knepley .seealso: DMCreate()
193467a56275SMatthew G Knepley @*/
193567a56275SMatthew G Knepley PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
193667a56275SMatthew G Knepley {
193767a56275SMatthew G Knepley   DM             B;
193867a56275SMatthew G Knepley   char           convname[256];
193967a56275SMatthew G Knepley   PetscBool      sametype, issame;
194067a56275SMatthew G Knepley   PetscErrorCode ierr;
194167a56275SMatthew G Knepley 
194267a56275SMatthew G Knepley   PetscFunctionBegin;
194367a56275SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
194467a56275SMatthew G Knepley   PetscValidType(dm,1);
194567a56275SMatthew G Knepley   PetscValidPointer(M,3);
194667a56275SMatthew G Knepley   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
194767a56275SMatthew G Knepley   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
194867a56275SMatthew G Knepley   {
194967a56275SMatthew G Knepley     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
195067a56275SMatthew G Knepley 
195167a56275SMatthew G Knepley     /*
195267a56275SMatthew G Knepley        Order of precedence:
195367a56275SMatthew G Knepley        1) See if a specialized converter is known to the current DM.
195467a56275SMatthew G Knepley        2) See if a specialized converter is known to the desired DM class.
195567a56275SMatthew G Knepley        3) See if a good general converter is registered for the desired class
195667a56275SMatthew G Knepley        4) See if a good general converter is known for the current matrix.
195767a56275SMatthew G Knepley        5) Use a really basic converter.
195867a56275SMatthew G Knepley     */
195967a56275SMatthew G Knepley 
196067a56275SMatthew G Knepley     /* 1) See if a specialized converter is known to the current DM and the desired class */
196167a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
196267a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
196367a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
196467a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
196567a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
196667a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
196767a56275SMatthew G Knepley     if (conv) goto foundconv;
196867a56275SMatthew G Knepley 
196967a56275SMatthew G Knepley     /* 2)  See if a specialized converter is known to the desired DM class. */
197067a56275SMatthew G Knepley     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
197167a56275SMatthew G Knepley     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
197267a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
197367a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
197467a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
197567a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
197667a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
197767a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
197867a56275SMatthew G Knepley     if (conv) {
1979fcfd50ebSBarry Smith       ierr = DMDestroy(&B);CHKERRQ(ierr);
198067a56275SMatthew G Knepley       goto foundconv;
198167a56275SMatthew G Knepley     }
198267a56275SMatthew G Knepley 
198367a56275SMatthew G Knepley #if 0
198467a56275SMatthew G Knepley     /* 3) See if a good general converter is registered for the desired class */
198567a56275SMatthew G Knepley     conv = B->ops->convertfrom;
1986fcfd50ebSBarry Smith     ierr = DMDestroy(&B);CHKERRQ(ierr);
198767a56275SMatthew G Knepley     if (conv) goto foundconv;
198867a56275SMatthew G Knepley 
198967a56275SMatthew G Knepley     /* 4) See if a good general converter is known for the current matrix */
199067a56275SMatthew G Knepley     if (dm->ops->convert) {
199167a56275SMatthew G Knepley       conv = dm->ops->convert;
199267a56275SMatthew G Knepley     }
199367a56275SMatthew G Knepley     if (conv) goto foundconv;
199467a56275SMatthew G Knepley #endif
199567a56275SMatthew G Knepley 
199667a56275SMatthew G Knepley     /* 5) Use a really basic converter. */
199767a56275SMatthew G Knepley     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
199867a56275SMatthew G Knepley 
199967a56275SMatthew G Knepley     foundconv:
200067a56275SMatthew G Knepley     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
200167a56275SMatthew G Knepley     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
200267a56275SMatthew G Knepley     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
200367a56275SMatthew G Knepley   }
200467a56275SMatthew G Knepley   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
200567a56275SMatthew G Knepley   PetscFunctionReturn(0);
200667a56275SMatthew G Knepley }
2007264ace61SBarry Smith 
2008264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
2009264ace61SBarry Smith 
2010264ace61SBarry Smith #undef __FUNCT__
2011264ace61SBarry Smith #define __FUNCT__ "DMRegister"
2012264ace61SBarry Smith /*@C
2013264ace61SBarry Smith   DMRegister - See DMRegisterDynamic()
2014264ace61SBarry Smith 
2015264ace61SBarry Smith   Level: advanced
2016264ace61SBarry Smith @*/
20177087cfbeSBarry Smith PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
2018264ace61SBarry Smith {
2019264ace61SBarry Smith   char fullname[PETSC_MAX_PATH_LEN];
2020264ace61SBarry Smith   PetscErrorCode ierr;
2021264ace61SBarry Smith 
2022264ace61SBarry Smith   PetscFunctionBegin;
2023264ace61SBarry Smith   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
2024264ace61SBarry Smith   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
2025264ace61SBarry Smith   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
2026264ace61SBarry Smith   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
2027264ace61SBarry Smith   PetscFunctionReturn(0);
2028264ace61SBarry Smith }
2029264ace61SBarry Smith 
2030264ace61SBarry Smith 
2031264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
2032264ace61SBarry Smith #undef __FUNCT__
2033264ace61SBarry Smith #define __FUNCT__ "DMRegisterDestroy"
2034264ace61SBarry Smith /*@C
2035264ace61SBarry Smith    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
2036264ace61SBarry Smith 
2037264ace61SBarry Smith    Not Collective
2038264ace61SBarry Smith 
2039264ace61SBarry Smith    Level: advanced
2040264ace61SBarry Smith 
2041264ace61SBarry Smith .keywords: DM, register, destroy
2042264ace61SBarry Smith .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
2043264ace61SBarry Smith @*/
20447087cfbeSBarry Smith PetscErrorCode  DMRegisterDestroy(void)
2045264ace61SBarry Smith {
2046264ace61SBarry Smith   PetscErrorCode ierr;
2047264ace61SBarry Smith 
2048264ace61SBarry Smith   PetscFunctionBegin;
2049264ace61SBarry Smith   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
2050264ace61SBarry Smith   DMRegisterAllCalled = PETSC_FALSE;
2051264ace61SBarry Smith   PetscFunctionReturn(0);
2052264ace61SBarry Smith }
205323f975d1SBarry Smith 
205423f975d1SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
2055c6db04a5SJed Brown #include <mex.h>
205623f975d1SBarry Smith 
20573014e516SBarry Smith typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
205823f975d1SBarry Smith 
205923f975d1SBarry Smith #undef __FUNCT__
206023f975d1SBarry Smith #define __FUNCT__ "DMComputeFunction_Matlab"
206123f975d1SBarry Smith /*
206223f975d1SBarry Smith    DMComputeFunction_Matlab - Calls the function that has been set with
206323f975d1SBarry Smith                          DMSetFunctionMatlab().
206423f975d1SBarry Smith 
206523f975d1SBarry Smith    For linear problems x is null
206623f975d1SBarry Smith 
206723f975d1SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
206823f975d1SBarry Smith */
20697087cfbeSBarry Smith PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
207023f975d1SBarry Smith {
207123f975d1SBarry Smith   PetscErrorCode    ierr;
207223f975d1SBarry Smith   DMMatlabContext   *sctx;
207323f975d1SBarry Smith   int               nlhs = 1,nrhs = 4;
207423f975d1SBarry Smith   mxArray	    *plhs[1],*prhs[4];
207523f975d1SBarry Smith   long long int     lx = 0,ly = 0,ls = 0;
207623f975d1SBarry Smith 
207723f975d1SBarry Smith   PetscFunctionBegin;
207823f975d1SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
207923f975d1SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
208023f975d1SBarry Smith   PetscCheckSameComm(dm,1,y,3);
208123f975d1SBarry Smith 
208223f975d1SBarry Smith   /* call Matlab function in ctx with arguments x and y */
20831b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
208423f975d1SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
208523f975d1SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
20863014e516SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
208723f975d1SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
208823f975d1SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
208923f975d1SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)ly);
209023f975d1SBarry Smith   prhs[3] =  mxCreateString(sctx->funcname);
2091b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
209223f975d1SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
209323f975d1SBarry Smith   mxDestroyArray(prhs[0]);
209423f975d1SBarry Smith   mxDestroyArray(prhs[1]);
209523f975d1SBarry Smith   mxDestroyArray(prhs[2]);
209623f975d1SBarry Smith   mxDestroyArray(prhs[3]);
209723f975d1SBarry Smith   mxDestroyArray(plhs[0]);
209823f975d1SBarry Smith   PetscFunctionReturn(0);
209923f975d1SBarry Smith }
210023f975d1SBarry Smith 
210123f975d1SBarry Smith 
210223f975d1SBarry Smith #undef __FUNCT__
210323f975d1SBarry Smith #define __FUNCT__ "DMSetFunctionMatlab"
210423f975d1SBarry Smith /*
210523f975d1SBarry Smith    DMSetFunctionMatlab - Sets the function evaluation routine
210623f975d1SBarry Smith 
210723f975d1SBarry Smith */
21087087cfbeSBarry Smith PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
210923f975d1SBarry Smith {
211023f975d1SBarry Smith   PetscErrorCode    ierr;
211123f975d1SBarry Smith   DMMatlabContext   *sctx;
211223f975d1SBarry Smith 
211323f975d1SBarry Smith   PetscFunctionBegin;
2114171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
211523f975d1SBarry Smith   /* currently sctx is memory bleed */
21161b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
21173014e516SBarry Smith   if (!sctx) {
211823f975d1SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
21193014e516SBarry Smith   }
212023f975d1SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
21211b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
212223f975d1SBarry Smith   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
212323f975d1SBarry Smith   PetscFunctionReturn(0);
212423f975d1SBarry Smith }
21253014e516SBarry Smith 
21263014e516SBarry Smith #undef __FUNCT__
21273014e516SBarry Smith #define __FUNCT__ "DMComputeJacobian_Matlab"
21283014e516SBarry Smith /*
21293014e516SBarry Smith    DMComputeJacobian_Matlab - Calls the function that has been set with
21303014e516SBarry Smith                          DMSetJacobianMatlab().
21313014e516SBarry Smith 
21323014e516SBarry Smith    For linear problems x is null
21333014e516SBarry Smith 
21343014e516SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
21353014e516SBarry Smith */
21367087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
21373014e516SBarry Smith {
21383014e516SBarry Smith   PetscErrorCode    ierr;
21393014e516SBarry Smith   DMMatlabContext   *sctx;
21403014e516SBarry Smith   int               nlhs = 2,nrhs = 5;
21413014e516SBarry Smith   mxArray	    *plhs[2],*prhs[5];
21423014e516SBarry Smith   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
21433014e516SBarry Smith 
21443014e516SBarry Smith   PetscFunctionBegin;
21453014e516SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
21463014e516SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
21473014e516SBarry Smith 
2148e3c5b3baSBarry Smith   /* call MATLAB function in ctx with arguments x, A, and B */
21491b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
21503014e516SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
21513014e516SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
21523014e516SBarry Smith   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
21533014e516SBarry Smith   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
21543014e516SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
21553014e516SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
21563014e516SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lA);
21573014e516SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lB);
21583014e516SBarry Smith   prhs[4] =  mxCreateString(sctx->jacname);
2159b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
2160c980e822SBarry Smith   *str    =  (MatStructure) mxGetScalar(plhs[0]);
2161c088a8dcSBarry Smith   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
21623014e516SBarry Smith   mxDestroyArray(prhs[0]);
21633014e516SBarry Smith   mxDestroyArray(prhs[1]);
21643014e516SBarry Smith   mxDestroyArray(prhs[2]);
21653014e516SBarry Smith   mxDestroyArray(prhs[3]);
21663014e516SBarry Smith   mxDestroyArray(prhs[4]);
21673014e516SBarry Smith   mxDestroyArray(plhs[0]);
21683014e516SBarry Smith   mxDestroyArray(plhs[1]);
21693014e516SBarry Smith   PetscFunctionReturn(0);
21703014e516SBarry Smith }
21713014e516SBarry Smith 
21723014e516SBarry Smith 
21733014e516SBarry Smith #undef __FUNCT__
21743014e516SBarry Smith #define __FUNCT__ "DMSetJacobianMatlab"
21753014e516SBarry Smith /*
21763014e516SBarry Smith    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
21773014e516SBarry Smith 
21783014e516SBarry Smith */
21797087cfbeSBarry Smith PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
21803014e516SBarry Smith {
21813014e516SBarry Smith   PetscErrorCode    ierr;
21823014e516SBarry Smith   DMMatlabContext   *sctx;
21833014e516SBarry Smith 
21843014e516SBarry Smith   PetscFunctionBegin;
2185171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
21863014e516SBarry Smith   /* currently sctx is memory bleed */
21871b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
21883014e516SBarry Smith   if (!sctx) {
21893014e516SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
21903014e516SBarry Smith   }
21913014e516SBarry Smith   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
21921b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
21933014e516SBarry Smith   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
21943014e516SBarry Smith   PetscFunctionReturn(0);
21953014e516SBarry Smith }
219623f975d1SBarry Smith #endif
2197b859378eSBarry Smith 
2198b859378eSBarry Smith #undef __FUNCT__
2199b859378eSBarry Smith #define __FUNCT__ "DMLoad"
2200b859378eSBarry Smith /*@C
2201b859378eSBarry Smith   DMLoad - Loads a DM that has been stored in binary or HDF5 format
2202b859378eSBarry Smith   with DMView().
2203b859378eSBarry Smith 
2204b859378eSBarry Smith   Collective on PetscViewer
2205b859378eSBarry Smith 
2206b859378eSBarry Smith   Input Parameters:
2207b859378eSBarry Smith + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2208b859378eSBarry Smith            some related function before a call to DMLoad().
2209b859378eSBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2210b859378eSBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2211b859378eSBarry Smith 
2212b859378eSBarry Smith    Level: intermediate
2213b859378eSBarry Smith 
2214b859378eSBarry Smith   Notes:
2215b859378eSBarry Smith   Defaults to the DM DA.
2216b859378eSBarry Smith 
2217b859378eSBarry Smith   Notes for advanced users:
2218b859378eSBarry Smith   Most users should not need to know the details of the binary storage
2219b859378eSBarry Smith   format, since DMLoad() and DMView() completely hide these details.
2220b859378eSBarry Smith   But for anyone who's interested, the standard binary matrix storage
2221b859378eSBarry Smith   format is
2222b859378eSBarry Smith .vb
2223b859378eSBarry Smith      has not yet been determined
2224b859378eSBarry Smith .ve
2225b859378eSBarry Smith 
2226b859378eSBarry Smith    In addition, PETSc automatically does the byte swapping for
2227b859378eSBarry Smith machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
2228b859378eSBarry Smith linux, Windows and the paragon; thus if you write your own binary
2229b859378eSBarry Smith read/write routines you have to swap the bytes; see PetscBinaryRead()
2230b859378eSBarry Smith and PetscBinaryWrite() to see how this may be done.
2231b859378eSBarry Smith 
2232b859378eSBarry Smith   Concepts: vector^loading from file
2233b859378eSBarry Smith 
2234b859378eSBarry Smith .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2235b859378eSBarry Smith @*/
2236b859378eSBarry Smith PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2237b859378eSBarry Smith {
2238b859378eSBarry Smith   PetscErrorCode ierr;
2239b859378eSBarry Smith 
2240b859378eSBarry Smith   PetscFunctionBegin;
2241b859378eSBarry Smith   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2242b859378eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2243b859378eSBarry Smith 
2244b859378eSBarry Smith   if (!((PetscObject)newdm)->type_name) {
2245b859378eSBarry Smith     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
2246b859378eSBarry Smith   }
2247b859378eSBarry Smith   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2248b859378eSBarry Smith   PetscFunctionReturn(0);
2249b859378eSBarry Smith }
2250b859378eSBarry Smith 
22517da65231SMatthew G Knepley /******************************** FEM Support **********************************/
22527da65231SMatthew G Knepley 
22537da65231SMatthew G Knepley #undef __FUNCT__
22547da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellVector"
22557da65231SMatthew G Knepley PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) {
22561d47ebbbSSatish Balay   PetscInt       f;
22571b30c384SMatthew G Knepley   PetscErrorCode ierr;
22581b30c384SMatthew G Knepley 
22597da65231SMatthew G Knepley   PetscFunctionBegin;
226074778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
22611d47ebbbSSatish Balay   for(f = 0; f < len; ++f) {
226274778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
22637da65231SMatthew G Knepley   }
22647da65231SMatthew G Knepley   PetscFunctionReturn(0);
22657da65231SMatthew G Knepley }
22667da65231SMatthew G Knepley 
22677da65231SMatthew G Knepley #undef __FUNCT__
22687da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellMatrix"
22697da65231SMatthew G Knepley PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) {
22701b30c384SMatthew G Knepley   PetscInt       f, g;
22717da65231SMatthew G Knepley   PetscErrorCode ierr;
22727da65231SMatthew G Knepley 
22737da65231SMatthew G Knepley   PetscFunctionBegin;
227474778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
22751d47ebbbSSatish Balay   for(f = 0; f < rows; ++f) {
227674778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
22771d47ebbbSSatish Balay     for(g = 0; g < cols; ++g) {
227874778d6cSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
22797da65231SMatthew G Knepley     }
228074778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
22817da65231SMatthew G Knepley   }
22827da65231SMatthew G Knepley   PetscFunctionReturn(0);
22837da65231SMatthew G Knepley }
2284e7c4fc90SDmitry Karpeev 
2285