xref: /petsc/src/dm/interface/dm.c (revision ae0a1c52d261b093567d7f5318c3432dadf7dc1d)
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__
879a42bb27SBarry Smith #define __FUNCT__ "DMSetOptionsPrefix"
889a42bb27SBarry Smith /*@C
899a42bb27SBarry Smith    DMSetOptionsPrefix - Sets the prefix used for searching for all
90aa219208SBarry Smith    DMDA options in the database.
919a42bb27SBarry Smith 
92aa219208SBarry Smith    Logically Collective on DMDA
939a42bb27SBarry Smith 
949a42bb27SBarry Smith    Input Parameter:
95aa219208SBarry Smith +  da - the DMDA context
969a42bb27SBarry Smith -  prefix - the prefix to prepend to all option names
979a42bb27SBarry Smith 
989a42bb27SBarry Smith    Notes:
999a42bb27SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1009a42bb27SBarry Smith    The first character of all runtime options is AUTOMATICALLY the hyphen.
1019a42bb27SBarry Smith 
1029a42bb27SBarry Smith    Level: advanced
1039a42bb27SBarry Smith 
104aa219208SBarry Smith .keywords: DMDA, set, options, prefix, database
1059a42bb27SBarry Smith 
1069a42bb27SBarry Smith .seealso: DMSetFromOptions()
1079a42bb27SBarry Smith @*/
1087087cfbeSBarry Smith PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
1099a42bb27SBarry Smith {
1109a42bb27SBarry Smith   PetscErrorCode ierr;
1119a42bb27SBarry Smith 
1129a42bb27SBarry Smith   PetscFunctionBegin;
1139a42bb27SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1149a42bb27SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
1159a42bb27SBarry Smith   PetscFunctionReturn(0);
1169a42bb27SBarry Smith }
1179a42bb27SBarry Smith 
1189a42bb27SBarry Smith #undef __FUNCT__
11947c6ae99SBarry Smith #define __FUNCT__ "DMDestroy"
12047c6ae99SBarry Smith /*@
121aa219208SBarry Smith     DMDestroy - Destroys a vector packer or DMDA.
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith     Collective on DM
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith     Input Parameter:
12647c6ae99SBarry Smith .   dm - the DM object to destroy
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith     Level: developer
12947c6ae99SBarry Smith 
130e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith @*/
133fcfd50ebSBarry Smith PetscErrorCode  DMDestroy(DM *dm)
13447c6ae99SBarry Smith {
135732e2eb9SMatthew G Knepley   PetscInt       i, cnt = 0;
136b17ce1afSJed Brown   DMCoarsenHookLink link,next;
13747c6ae99SBarry Smith   PetscErrorCode ierr;
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith   PetscFunctionBegin;
1406bf464f9SBarry Smith   if (!*dm) PetscFunctionReturn(0);
1416bf464f9SBarry Smith   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
14287e657c6SBarry Smith 
14387e657c6SBarry Smith   /* count all the circular references of DM and its contained Vecs */
144732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1456bf464f9SBarry Smith     if ((*dm)->localin[i])  {cnt++;}
1466bf464f9SBarry Smith     if ((*dm)->globalin[i]) {cnt++;}
147732e2eb9SMatthew G Knepley   }
14871cd77b2SBarry Smith   if ((*dm)->x) {
14971cd77b2SBarry Smith     PetscObject obj;
15071cd77b2SBarry Smith     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
15171cd77b2SBarry Smith     if (obj == (PetscObject)*dm) cnt++;
15271cd77b2SBarry Smith   }
153732e2eb9SMatthew G Knepley 
1546bf464f9SBarry Smith   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
155732e2eb9SMatthew G Knepley   /*
156732e2eb9SMatthew G Knepley      Need this test because the dm references the vectors that
157732e2eb9SMatthew G Knepley      reference the dm, so destroying the dm calls destroy on the
158732e2eb9SMatthew G Knepley      vectors that cause another destroy on the dm
159732e2eb9SMatthew G Knepley   */
1606bf464f9SBarry Smith   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
1616bf464f9SBarry Smith   ((PetscObject) (*dm))->refct = 0;
162732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
1636bf464f9SBarry Smith     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
1646bf464f9SBarry Smith     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
165732e2eb9SMatthew G Knepley   }
1661a266240SBarry Smith 
167b17ce1afSJed Brown   /* Destroy the list of hooks */
168b17ce1afSJed Brown   for (link=(*dm)->coarsenhook; link; link=next) {
169b17ce1afSJed Brown     next = link->next;
170b17ce1afSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
171b17ce1afSJed Brown   }
172b17ce1afSJed Brown   (*dm)->coarsenhook = PETSC_NULL;
173b17ce1afSJed Brown 
1741a266240SBarry Smith   if ((*dm)->ctx && (*dm)->ctxdestroy) {
1751a266240SBarry Smith     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
1761a266240SBarry Smith   }
17787e657c6SBarry Smith   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
17871cd77b2SBarry Smith   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
1794dcab191SBarry Smith   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
1806bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
1816bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
1826bf464f9SBarry Smith   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
183073dac72SJed Brown   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
184a89ea682SMatthew G Knepley   ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr);
185732e2eb9SMatthew G Knepley   /* if memory was published with AMS then destroy it */
1866bf464f9SBarry Smith   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
187732e2eb9SMatthew G Knepley 
1886bf464f9SBarry Smith   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
1896bf464f9SBarry Smith   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
190732e2eb9SMatthew G Knepley   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
19147c6ae99SBarry Smith   PetscFunctionReturn(0);
19247c6ae99SBarry Smith }
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith #undef __FUNCT__
195d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp"
196d7bf68aeSBarry Smith /*@
197d7bf68aeSBarry Smith     DMSetUp - sets up the data structures inside a DM object
198d7bf68aeSBarry Smith 
199d7bf68aeSBarry Smith     Collective on DM
200d7bf68aeSBarry Smith 
201d7bf68aeSBarry Smith     Input Parameter:
202d7bf68aeSBarry Smith .   dm - the DM object to setup
203d7bf68aeSBarry Smith 
204d7bf68aeSBarry Smith     Level: developer
205d7bf68aeSBarry Smith 
206e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
207d7bf68aeSBarry Smith 
208d7bf68aeSBarry Smith @*/
2097087cfbeSBarry Smith PetscErrorCode  DMSetUp(DM dm)
210d7bf68aeSBarry Smith {
211d7bf68aeSBarry Smith   PetscErrorCode ierr;
212d7bf68aeSBarry Smith 
213d7bf68aeSBarry Smith   PetscFunctionBegin;
214171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2158387afaaSJed Brown   if (dm->setupcalled) PetscFunctionReturn(0);
216d7bf68aeSBarry Smith   if (dm->ops->setup) {
217d7bf68aeSBarry Smith     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
218d7bf68aeSBarry Smith   }
2198387afaaSJed Brown   dm->setupcalled = PETSC_TRUE;
220d7bf68aeSBarry Smith   PetscFunctionReturn(0);
221d7bf68aeSBarry Smith }
222d7bf68aeSBarry Smith 
223d7bf68aeSBarry Smith #undef __FUNCT__
224d7bf68aeSBarry Smith #define __FUNCT__ "DMSetFromOptions"
225d7bf68aeSBarry Smith /*@
226d7bf68aeSBarry Smith     DMSetFromOptions - sets parameters in a DM from the options database
227d7bf68aeSBarry Smith 
228d7bf68aeSBarry Smith     Collective on DM
229d7bf68aeSBarry Smith 
230d7bf68aeSBarry Smith     Input Parameter:
231d7bf68aeSBarry Smith .   dm - the DM object to set options for
232d7bf68aeSBarry Smith 
233732e2eb9SMatthew G Knepley     Options Database:
234dd85299cSBarry Smith +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
235dd85299cSBarry Smith .   -dm_vec_type <type>  type of vector to create inside DM
236171400e9SBarry Smith .   -dm_mat_type <type>  type of matrix to create inside DM
237171400e9SBarry Smith -   -dm_coloring_type <global or ghosted>
238732e2eb9SMatthew G Knepley 
239d7bf68aeSBarry Smith     Level: developer
240d7bf68aeSBarry Smith 
241e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
242d7bf68aeSBarry Smith 
243d7bf68aeSBarry Smith @*/
2447087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions(DM dm)
245d7bf68aeSBarry Smith {
24667ad5babSMatthew G Knepley   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg;
247d7bf68aeSBarry Smith   PetscErrorCode ierr;
248f9ba7244SBarry Smith   char           typeName[256] = MATAIJ;
249d7bf68aeSBarry Smith 
250d7bf68aeSBarry Smith   PetscFunctionBegin;
251171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2523194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
25382fcb398SMatthew G Knepley     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
254c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr);
255c4721b0eSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr);
25667ad5babSMatthew G Knepley     ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr);
257073dac72SJed 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);
258f9ba7244SBarry Smith     ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
259f9ba7244SBarry Smith     if (flg) {
260f9ba7244SBarry Smith       ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
261f9ba7244SBarry Smith     }
262f9ba7244SBarry Smith     ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr);
263073dac72SJed Brown     if (flg) {
264073dac72SJed Brown       ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
265f9ba7244SBarry Smith       ierr = PetscStrallocpy(typeName,&dm->mattype);CHKERRQ(ierr);
266073dac72SJed Brown     }
2671b89239cSHong 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);
268f9ba7244SBarry Smith     if (dm->ops->setfromoptions) {
269f9ba7244SBarry Smith       ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
270f9ba7244SBarry Smith     }
271f9ba7244SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
272f9ba7244SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
27382fcb398SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
27482fcb398SMatthew G Knepley   if (flg1) {
27582fcb398SMatthew G Knepley     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
27682fcb398SMatthew G Knepley   }
277c4721b0eSMatthew G Knepley   if (flg2) {
278c4721b0eSMatthew G Knepley     PetscViewer viewer;
279c4721b0eSMatthew G Knepley 
280c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
281c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
282c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
283c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
284c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
285c4721b0eSMatthew G Knepley   }
286c4721b0eSMatthew G Knepley   if (flg3) {
287c4721b0eSMatthew G Knepley     PetscViewer viewer;
288c4721b0eSMatthew G Knepley 
289c4721b0eSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
290c4721b0eSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
291c4721b0eSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
292c4721b0eSMatthew G Knepley     ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr);
293c4721b0eSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
294c4721b0eSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
295c4721b0eSMatthew G Knepley   }
29667ad5babSMatthew G Knepley   if (flg4) {
29767ad5babSMatthew G Knepley     PetscViewer viewer;
29867ad5babSMatthew G Knepley 
29967ad5babSMatthew G Knepley     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
30067ad5babSMatthew G Knepley     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
30167ad5babSMatthew G Knepley     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr);
30267ad5babSMatthew G Knepley     ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr);
30367ad5babSMatthew G Knepley     ierr = DMView(dm, viewer);CHKERRQ(ierr);
30467ad5babSMatthew G Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
30567ad5babSMatthew G Knepley   }
306d7bf68aeSBarry Smith   PetscFunctionReturn(0);
307d7bf68aeSBarry Smith }
308d7bf68aeSBarry Smith 
309d7bf68aeSBarry Smith #undef __FUNCT__
31047c6ae99SBarry Smith #define __FUNCT__ "DMView"
311fc9bc008SSatish Balay /*@C
312aa219208SBarry Smith     DMView - Views a vector packer or DMDA.
31347c6ae99SBarry Smith 
31447c6ae99SBarry Smith     Collective on DM
31547c6ae99SBarry Smith 
31647c6ae99SBarry Smith     Input Parameter:
31747c6ae99SBarry Smith +   dm - the DM object to view
31847c6ae99SBarry Smith -   v - the viewer
31947c6ae99SBarry Smith 
32047c6ae99SBarry Smith     Level: developer
32147c6ae99SBarry Smith 
322e727c939SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith @*/
3257087cfbeSBarry Smith PetscErrorCode  DMView(DM dm,PetscViewer v)
32647c6ae99SBarry Smith {
32747c6ae99SBarry Smith   PetscErrorCode ierr;
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith   PetscFunctionBegin;
330171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3313014e516SBarry Smith  if (!v) {
3323014e516SBarry Smith     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
3333014e516SBarry Smith   }
3340c010503SBarry Smith   if (dm->ops->view) {
3350c010503SBarry Smith     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
33647c6ae99SBarry Smith   }
33747c6ae99SBarry Smith   PetscFunctionReturn(0);
33847c6ae99SBarry Smith }
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith #undef __FUNCT__
34147c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector"
34247c6ae99SBarry Smith /*@
343aa219208SBarry Smith     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
34447c6ae99SBarry Smith 
34547c6ae99SBarry Smith     Collective on DM
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith     Input Parameter:
34847c6ae99SBarry Smith .   dm - the DM object
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith     Output Parameter:
35147c6ae99SBarry Smith .   vec - the global vector
35247c6ae99SBarry Smith 
353073dac72SJed Brown     Level: beginner
35447c6ae99SBarry Smith 
355e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith @*/
3587087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
35947c6ae99SBarry Smith {
36047c6ae99SBarry Smith   PetscErrorCode ierr;
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith   PetscFunctionBegin;
363171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36447c6ae99SBarry Smith   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
36547c6ae99SBarry Smith   PetscFunctionReturn(0);
36647c6ae99SBarry Smith }
36747c6ae99SBarry Smith 
36847c6ae99SBarry Smith #undef __FUNCT__
36947c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector"
37047c6ae99SBarry Smith /*@
371aa219208SBarry Smith     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
37247c6ae99SBarry Smith 
37347c6ae99SBarry Smith     Not Collective
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith     Input Parameter:
37647c6ae99SBarry Smith .   dm - the DM object
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith     Output Parameter:
37947c6ae99SBarry Smith .   vec - the local vector
38047c6ae99SBarry Smith 
381073dac72SJed Brown     Level: beginner
38247c6ae99SBarry Smith 
383e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith @*/
3867087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
38747c6ae99SBarry Smith {
38847c6ae99SBarry Smith   PetscErrorCode ierr;
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith   PetscFunctionBegin;
391171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
39247c6ae99SBarry Smith   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
39347c6ae99SBarry Smith   PetscFunctionReturn(0);
39447c6ae99SBarry Smith }
39547c6ae99SBarry Smith 
39647c6ae99SBarry Smith #undef __FUNCT__
3971411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping"
3981411c6eeSJed Brown /*@
3991411c6eeSJed Brown    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
4001411c6eeSJed Brown 
4011411c6eeSJed Brown    Collective on DM
4021411c6eeSJed Brown 
4031411c6eeSJed Brown    Input Parameter:
4041411c6eeSJed Brown .  dm - the DM that provides the mapping
4051411c6eeSJed Brown 
4061411c6eeSJed Brown    Output Parameter:
4071411c6eeSJed Brown .  ltog - the mapping
4081411c6eeSJed Brown 
4091411c6eeSJed Brown    Level: intermediate
4101411c6eeSJed Brown 
4111411c6eeSJed Brown    Notes:
4121411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMapping() or
4131411c6eeSJed Brown    MatSetLocalToGlobalMapping().
4141411c6eeSJed Brown 
4151411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
4161411c6eeSJed Brown @*/
4177087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
4181411c6eeSJed Brown {
4191411c6eeSJed Brown   PetscErrorCode ierr;
4201411c6eeSJed Brown 
4211411c6eeSJed Brown   PetscFunctionBegin;
4221411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4231411c6eeSJed Brown   PetscValidPointer(ltog,2);
4241411c6eeSJed Brown   if (!dm->ltogmap) {
4251411c6eeSJed Brown     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
4261411c6eeSJed Brown     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
4271411c6eeSJed Brown   }
4281411c6eeSJed Brown   *ltog = dm->ltogmap;
4291411c6eeSJed Brown   PetscFunctionReturn(0);
4301411c6eeSJed Brown }
4311411c6eeSJed Brown 
4321411c6eeSJed Brown #undef __FUNCT__
4331411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
4341411c6eeSJed Brown /*@
4351411c6eeSJed Brown    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
4361411c6eeSJed Brown 
4371411c6eeSJed Brown    Collective on DM
4381411c6eeSJed Brown 
4391411c6eeSJed Brown    Input Parameter:
4401411c6eeSJed Brown .  da - the distributed array that provides the mapping
4411411c6eeSJed Brown 
4421411c6eeSJed Brown    Output Parameter:
4431411c6eeSJed Brown .  ltog - the block mapping
4441411c6eeSJed Brown 
4451411c6eeSJed Brown    Level: intermediate
4461411c6eeSJed Brown 
4471411c6eeSJed Brown    Notes:
4481411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
4491411c6eeSJed Brown    MatSetLocalToGlobalMappingBlock().
4501411c6eeSJed Brown 
4511411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
4521411c6eeSJed Brown @*/
4537087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
4541411c6eeSJed Brown {
4551411c6eeSJed Brown   PetscErrorCode ierr;
4561411c6eeSJed Brown 
4571411c6eeSJed Brown   PetscFunctionBegin;
4581411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4591411c6eeSJed Brown   PetscValidPointer(ltog,2);
4601411c6eeSJed Brown   if (!dm->ltogmapb) {
4611411c6eeSJed Brown     PetscInt bs;
4621411c6eeSJed Brown     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
4631411c6eeSJed Brown     if (bs > 1) {
4641411c6eeSJed Brown       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
4651411c6eeSJed Brown       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
4661411c6eeSJed Brown     } else {
4671411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
4681411c6eeSJed Brown       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
4691411c6eeSJed Brown     }
4701411c6eeSJed Brown   }
4711411c6eeSJed Brown   *ltog = dm->ltogmapb;
4721411c6eeSJed Brown   PetscFunctionReturn(0);
4731411c6eeSJed Brown }
4741411c6eeSJed Brown 
4751411c6eeSJed Brown #undef __FUNCT__
4761411c6eeSJed Brown #define __FUNCT__ "DMGetBlockSize"
4771411c6eeSJed Brown /*@
4781411c6eeSJed Brown    DMGetBlockSize - Gets the inherent block size associated with a DM
4791411c6eeSJed Brown 
4801411c6eeSJed Brown    Not Collective
4811411c6eeSJed Brown 
4821411c6eeSJed Brown    Input Parameter:
4831411c6eeSJed Brown .  dm - the DM with block structure
4841411c6eeSJed Brown 
4851411c6eeSJed Brown    Output Parameter:
4861411c6eeSJed Brown .  bs - the block size, 1 implies no exploitable block structure
4871411c6eeSJed Brown 
4881411c6eeSJed Brown    Level: intermediate
4891411c6eeSJed Brown 
4901411c6eeSJed Brown .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
4911411c6eeSJed Brown @*/
4927087cfbeSBarry Smith PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
4931411c6eeSJed Brown {
4941411c6eeSJed Brown   PetscFunctionBegin;
4951411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4961411c6eeSJed Brown   PetscValidPointer(bs,2);
4971411c6eeSJed 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");
4981411c6eeSJed Brown   *bs = dm->bs;
4991411c6eeSJed Brown   PetscFunctionReturn(0);
5001411c6eeSJed Brown }
5011411c6eeSJed Brown 
5021411c6eeSJed Brown #undef __FUNCT__
503e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation"
50447c6ae99SBarry Smith /*@
505e727c939SJed Brown     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith     Collective on DM
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith     Input Parameter:
51047c6ae99SBarry Smith +   dm1 - the DM object
51147c6ae99SBarry Smith -   dm2 - the second, finer DM object
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith     Output Parameter:
51447c6ae99SBarry Smith +  mat - the interpolation
51547c6ae99SBarry Smith -  vec - the scaling (optional)
51647c6ae99SBarry Smith 
51747c6ae99SBarry Smith     Level: developer
51847c6ae99SBarry Smith 
51985afcc9aSBarry 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
52085afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
521d52bd9f3SBarry Smith 
522d52bd9f3SBarry Smith         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
523d52bd9f3SBarry Smith         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
52485afcc9aSBarry Smith 
52585afcc9aSBarry Smith 
526e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
52747c6ae99SBarry Smith 
52847c6ae99SBarry Smith @*/
529e727c939SJed Brown PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
53047c6ae99SBarry Smith {
53147c6ae99SBarry Smith   PetscErrorCode ierr;
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith   PetscFunctionBegin;
534171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
535171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
53625296bd5SBarry Smith   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
53747c6ae99SBarry Smith   PetscFunctionReturn(0);
53847c6ae99SBarry Smith }
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith #undef __FUNCT__
541e727c939SJed Brown #define __FUNCT__ "DMCreateInjection"
54247c6ae99SBarry Smith /*@
543e727c939SJed Brown     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith     Collective on DM
54647c6ae99SBarry Smith 
54747c6ae99SBarry Smith     Input Parameter:
54847c6ae99SBarry Smith +   dm1 - the DM object
54947c6ae99SBarry Smith -   dm2 - the second, finer DM object
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith     Output Parameter:
55247c6ae99SBarry Smith .   ctx - the injection
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith     Level: developer
55547c6ae99SBarry Smith 
55685afcc9aSBarry 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
55785afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
55885afcc9aSBarry Smith 
559e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
56047c6ae99SBarry Smith 
56147c6ae99SBarry Smith @*/
562e727c939SJed Brown PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
56347c6ae99SBarry Smith {
56447c6ae99SBarry Smith   PetscErrorCode ierr;
56547c6ae99SBarry Smith 
56647c6ae99SBarry Smith   PetscFunctionBegin;
567171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
568171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
56947c6ae99SBarry Smith   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
57047c6ae99SBarry Smith   PetscFunctionReturn(0);
57147c6ae99SBarry Smith }
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith #undef __FUNCT__
574e727c939SJed Brown #define __FUNCT__ "DMCreateColoring"
575d1e2c406SBarry Smith /*@C
576e727c939SJed Brown     DMCreateColoring - Gets coloring for a DMDA or DMComposite
57747c6ae99SBarry Smith 
57847c6ae99SBarry Smith     Collective on DM
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith     Input Parameter:
58147c6ae99SBarry Smith +   dm - the DM object
58247c6ae99SBarry Smith .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
58347c6ae99SBarry Smith -   matype - either MATAIJ or MATBAIJ
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith     Output Parameter:
58647c6ae99SBarry Smith .   coloring - the coloring
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith     Level: developer
58947c6ae99SBarry Smith 
590e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
59147c6ae99SBarry Smith 
592aab9d709SJed Brown @*/
593e727c939SJed Brown PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
59447c6ae99SBarry Smith {
59547c6ae99SBarry Smith   PetscErrorCode ierr;
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith   PetscFunctionBegin;
598171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
59947c6ae99SBarry Smith   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
60047c6ae99SBarry Smith   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
60147c6ae99SBarry Smith   PetscFunctionReturn(0);
60247c6ae99SBarry Smith }
60347c6ae99SBarry Smith 
60447c6ae99SBarry Smith #undef __FUNCT__
605950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix"
60647c6ae99SBarry Smith /*@C
607950540a4SJed Brown     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
60847c6ae99SBarry Smith 
60947c6ae99SBarry Smith     Collective on DM
61047c6ae99SBarry Smith 
61147c6ae99SBarry Smith     Input Parameter:
61247c6ae99SBarry Smith +   dm - the DM object
61347c6ae99SBarry Smith -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
61494013140SBarry Smith             any type which inherits from one of these (such as MATAIJ)
61547c6ae99SBarry Smith 
61647c6ae99SBarry Smith     Output Parameter:
61747c6ae99SBarry Smith .   mat - the empty Jacobian
61847c6ae99SBarry Smith 
619073dac72SJed Brown     Level: beginner
62047c6ae99SBarry Smith 
62194013140SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
62294013140SBarry Smith        do not need to do it yourself.
62394013140SBarry Smith 
62494013140SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
625aa219208SBarry Smith        the nonzero pattern call DMDASetMatPreallocateOnly()
62694013140SBarry Smith 
62794013140SBarry 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
62894013140SBarry Smith        internally by PETSc.
62994013140SBarry Smith 
63094013140SBarry Smith        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
631aa219208SBarry Smith        the indices for the global numbering for DMDAs which is complicated.
63294013140SBarry Smith 
633e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
63447c6ae99SBarry Smith 
635aab9d709SJed Brown @*/
636950540a4SJed Brown PetscErrorCode  DMCreateMatrix(DM dm,const MatType mtype,Mat *mat)
63747c6ae99SBarry Smith {
63847c6ae99SBarry Smith   PetscErrorCode ierr;
63947c6ae99SBarry Smith 
64047c6ae99SBarry Smith   PetscFunctionBegin;
641171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
642235683edSBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
643235683edSBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
644235683edSBarry Smith #endif
645c7b7c8a4SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
646c7b7c8a4SJed Brown   PetscValidPointer(mat,3);
647073dac72SJed Brown   if (dm->mattype) {
64825296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
649073dac72SJed Brown   } else {
65025296bd5SBarry Smith     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
651c7b7c8a4SJed Brown   }
65247c6ae99SBarry Smith   PetscFunctionReturn(0);
65347c6ae99SBarry Smith }
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith #undef __FUNCT__
656732e2eb9SMatthew G Knepley #define __FUNCT__ "DMSetMatrixPreallocateOnly"
657732e2eb9SMatthew G Knepley /*@
658950540a4SJed Brown   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
659732e2eb9SMatthew G Knepley     preallocated but the nonzero structure and zero values will not be set.
660732e2eb9SMatthew G Knepley 
661732e2eb9SMatthew G Knepley   Logically Collective on DMDA
662732e2eb9SMatthew G Knepley 
663732e2eb9SMatthew G Knepley   Input Parameter:
664732e2eb9SMatthew G Knepley + dm - the DM
665732e2eb9SMatthew G Knepley - only - PETSC_TRUE if only want preallocation
666732e2eb9SMatthew G Knepley 
667732e2eb9SMatthew G Knepley   Level: developer
668950540a4SJed Brown .seealso DMCreateMatrix()
669732e2eb9SMatthew G Knepley @*/
670732e2eb9SMatthew G Knepley PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
671732e2eb9SMatthew G Knepley {
672732e2eb9SMatthew G Knepley   PetscFunctionBegin;
673732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
674732e2eb9SMatthew G Knepley   dm->prealloc_only = only;
675732e2eb9SMatthew G Knepley   PetscFunctionReturn(0);
676732e2eb9SMatthew G Knepley }
677732e2eb9SMatthew G Knepley 
678732e2eb9SMatthew G Knepley #undef __FUNCT__
679a89ea682SMatthew G Knepley #define __FUNCT__ "DMGetWorkArray"
680a89ea682SMatthew G Knepley /*@C
681a89ea682SMatthew G Knepley   DMGetWorkArray - Gets a work array guaranteed to be at least the input size
682a89ea682SMatthew G Knepley 
683a89ea682SMatthew G Knepley   Not Collective
684a89ea682SMatthew G Knepley 
685a89ea682SMatthew G Knepley   Input Parameters:
686a89ea682SMatthew G Knepley + dm - the DM object
687a89ea682SMatthew G Knepley - size - The minium size
688a89ea682SMatthew G Knepley 
689a89ea682SMatthew G Knepley   Output Parameter:
690a89ea682SMatthew G Knepley . array - the work array
691a89ea682SMatthew G Knepley 
692a89ea682SMatthew G Knepley   Level: developer
693a89ea682SMatthew G Knepley 
694a89ea682SMatthew G Knepley .seealso DMDestroy(), DMCreate()
695a89ea682SMatthew G Knepley @*/
696a89ea682SMatthew G Knepley PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array)
697a89ea682SMatthew G Knepley {
698a89ea682SMatthew G Knepley   PetscErrorCode ierr;
699a89ea682SMatthew G Knepley 
700a89ea682SMatthew G Knepley   PetscFunctionBegin;
701a89ea682SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
702a89ea682SMatthew G Knepley   PetscValidPointer(array,3);
703a89ea682SMatthew G Knepley   if (size > dm->workSize) {
704a89ea682SMatthew G Knepley     dm->workSize = size;
705a89ea682SMatthew G Knepley     ierr = PetscFree(dm->workArray);CHKERRQ(ierr);
706a89ea682SMatthew G Knepley     ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr);
707a89ea682SMatthew G Knepley   }
708a89ea682SMatthew G Knepley   *array = dm->workArray;
709a89ea682SMatthew G Knepley   PetscFunctionReturn(0);
710a89ea682SMatthew G Knepley }
711a89ea682SMatthew G Knepley 
7124d343eeaSMatthew G Knepley #undef __FUNCT__
713e7c4fc90SDmitry Karpeev #define __FUNCT__ "DMCreateDecompositionDM"
714e7c4fc90SDmitry Karpeev /*@C
715e7c4fc90SDmitry Karpeev   DMCreateDecompositionDM - creates a DM that defines a decomposition of the original DM.
716e7c4fc90SDmitry Karpeev 
717e7c4fc90SDmitry Karpeev   Not Collective
718e7c4fc90SDmitry Karpeev 
719e7c4fc90SDmitry Karpeev   Input Parameters:
720e7c4fc90SDmitry Karpeev + dm   - the DM object
721e7c4fc90SDmitry Karpeev - name - the name of the decomposition
722e7c4fc90SDmitry Karpeev 
723e7c4fc90SDmitry Karpeev   Output Parameter:
724e7c4fc90SDmitry Karpeev . ddm  - the decomposition DM (PETSC_NULL, if no such decomposition is known)
725e7c4fc90SDmitry Karpeev 
726e7c4fc90SDmitry Karpeev   Level: advanced
727e7c4fc90SDmitry Karpeev 
728e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMCreate(), DMCreateDecomposition()
729e7c4fc90SDmitry Karpeev @*/
730e7c4fc90SDmitry Karpeev PetscErrorCode DMCreateDecompositionDM(DM dm, const char* name, DM *ddm)
731e7c4fc90SDmitry Karpeev {
732e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
733e7c4fc90SDmitry Karpeev 
734e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
735e7c4fc90SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
736e7c4fc90SDmitry Karpeev   PetscValidCharPointer(name,2);
737e7c4fc90SDmitry Karpeev   PetscValidPointer(ddm,3);
738e7c4fc90SDmitry Karpeev   if(!dm->ops->createdecompositiondm) {
739e7c4fc90SDmitry Karpeev     *ddm = PETSC_NULL;
740e7c4fc90SDmitry Karpeev   }
741e7c4fc90SDmitry Karpeev   else {
742e7c4fc90SDmitry Karpeev     ierr = (*dm->ops->createdecompositiondm)(dm,name,ddm); CHKERRQ(ierr);
743e7c4fc90SDmitry Karpeev   }
744e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
745e7c4fc90SDmitry Karpeev }
746e7c4fc90SDmitry Karpeev 
747e7c4fc90SDmitry Karpeev #undef __FUNCT__
7484d343eeaSMatthew G Knepley #define __FUNCT__ "DMCreateFieldIS"
7494f3b5142SJed Brown /*@C
7504d343eeaSMatthew G Knepley   DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field
7514d343eeaSMatthew G Knepley 
7524d343eeaSMatthew G Knepley   Not collective
7534d343eeaSMatthew G Knepley 
7544d343eeaSMatthew G Knepley   Input Parameter:
7554d343eeaSMatthew G Knepley . dm - the DM object
7564d343eeaSMatthew G Knepley 
7574d343eeaSMatthew G Knepley   Output Parameters:
75821c9b008SJed Brown + numFields - The number of fields (or PETSC_NULL if not requested)
75921c9b008SJed Brown . names     - The name for each field (or PETSC_NULL if not requested)
76021c9b008SJed Brown - fields    - The global indices for each field (or PETSC_NULL if not requested)
7614d343eeaSMatthew G Knepley 
7624d343eeaSMatthew G Knepley   Level: intermediate
7634d343eeaSMatthew G Knepley 
76421c9b008SJed Brown   Notes:
76521c9b008SJed Brown   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
76621c9b008SJed Brown   PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with
76721c9b008SJed Brown   PetscFree().
76821c9b008SJed Brown 
7694d343eeaSMatthew G Knepley .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
7704d343eeaSMatthew G Knepley @*/
77121c9b008SJed Brown PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***names, IS **fields)
7724d343eeaSMatthew G Knepley {
7734d343eeaSMatthew G Knepley   PetscErrorCode ierr;
7744d343eeaSMatthew G Knepley 
7754d343eeaSMatthew G Knepley   PetscFunctionBegin;
7764d343eeaSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7774d343eeaSMatthew G Knepley   if (numFields) {PetscValidPointer(numFields,2);}
7784d343eeaSMatthew G Knepley   if (names) {PetscValidPointer(names,3);}
7794d343eeaSMatthew G Knepley   if (fields) {PetscValidPointer(fields,4);}
7804d343eeaSMatthew G Knepley   ierr = (*dm->ops->createfieldis)(dm, numFields, names, fields);CHKERRQ(ierr);
7814d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
7824d343eeaSMatthew G Knepley }
7834d343eeaSMatthew G Knepley 
7845fe1f584SPeter Brune 
785a89ea682SMatthew G Knepley #undef __FUNCT__
786e7c4fc90SDmitry Karpeev #define __FUNCT__ "DMCreateDecomposition"
787e7c4fc90SDmitry Karpeev /*@C
788e7c4fc90SDmitry Karpeev   DMCreateDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems:
789e7c4fc90SDmitry Karpeev                           each IS contains the global indices of the dofs of the corresponding subproblem.
790e7c4fc90SDmitry Karpeev                           The optional list of DMs define the DM for each subproblem.
791e7c4fc90SDmitry Karpeev                           Generalizes DMCreateFieldIS().
792e7c4fc90SDmitry Karpeev 
793e7c4fc90SDmitry Karpeev   Not collective
794e7c4fc90SDmitry Karpeev 
795e7c4fc90SDmitry Karpeev   Input Parameter:
796e7c4fc90SDmitry Karpeev . dm - the DM object
797e7c4fc90SDmitry Karpeev 
798e7c4fc90SDmitry Karpeev   Output Parameters:
799e7c4fc90SDmitry Karpeev + len       - The number of subproblems in the decomposition (or PETSC_NULL if not requested)
800e7c4fc90SDmitry Karpeev . namelist  - The name for each subproblem (or PETSC_NULL if not requested)
801e7c4fc90SDmitry Karpeev . islist    - The global indices for each subproblem (or PETSC_NULL if not requested)
802e7c4fc90SDmitry Karpeev - dmlist    - The DMs for each subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined)
803e7c4fc90SDmitry Karpeev 
804e7c4fc90SDmitry Karpeev   Level: intermediate
805e7c4fc90SDmitry Karpeev 
806e7c4fc90SDmitry Karpeev   Notes:
807e7c4fc90SDmitry Karpeev   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
808e7c4fc90SDmitry Karpeev   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
809e7c4fc90SDmitry Karpeev   and all of the arrays should be freed with PetscFree().
810e7c4fc90SDmitry Karpeev 
811e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
812e7c4fc90SDmitry Karpeev @*/
813e7c4fc90SDmitry Karpeev PetscErrorCode DMCreateDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
814e7c4fc90SDmitry Karpeev {
815e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
816e7c4fc90SDmitry Karpeev 
817e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
818e7c4fc90SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
819e7c4fc90SDmitry Karpeev   if (len) {PetscValidPointer(len,2);}
820e7c4fc90SDmitry Karpeev   if (namelist) {PetscValidPointer(namelist,3);}
821e7c4fc90SDmitry Karpeev   if (islist) {PetscValidPointer(islist,4);}
822e7c4fc90SDmitry Karpeev   if (dmlist) {PetscValidPointer(dmlist,5);}
823e7c4fc90SDmitry Karpeev   if(!dm->ops->createdecomposition) {
824e7c4fc90SDmitry Karpeev     ierr = (*dm->ops->createfieldis)(dm, len, namelist, islist);CHKERRQ(ierr);
825e7c4fc90SDmitry Karpeev     /* By default there are no DMs associated with subproblems. */
826e7c4fc90SDmitry Karpeev     if(dmlist) *dmlist = PETSC_NULL;
827e7c4fc90SDmitry Karpeev   }
828e7c4fc90SDmitry Karpeev   else {
829e7c4fc90SDmitry Karpeev     ierr = (*dm->ops->createdecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr);
830e7c4fc90SDmitry Karpeev   }
831e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
832e7c4fc90SDmitry Karpeev }
833e7c4fc90SDmitry Karpeev 
834e7c4fc90SDmitry Karpeev 
835e7c4fc90SDmitry Karpeev #undef __FUNCT__
83647c6ae99SBarry Smith #define __FUNCT__ "DMRefine"
83747c6ae99SBarry Smith /*@
83847c6ae99SBarry Smith   DMRefine - Refines a DM object
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith   Collective on DM
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith   Input Parameter:
84347c6ae99SBarry Smith + dm   - the DM object
84447c6ae99SBarry Smith - comm - the communicator to contain the new DM object (or PETSC_NULL)
84547c6ae99SBarry Smith 
84647c6ae99SBarry Smith   Output Parameter:
847*ae0a1c52SMatthew G Knepley . dmf - the refined DM, or PETSC_NULL
848*ae0a1c52SMatthew G Knepley 
849*ae0a1c52SMatthew G Knepley   Note: If no refinement was done, the return value is PETSC_NULL
85047c6ae99SBarry Smith 
85147c6ae99SBarry Smith   Level: developer
85247c6ae99SBarry Smith 
853e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
85447c6ae99SBarry Smith @*/
8557087cfbeSBarry Smith PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
85647c6ae99SBarry Smith {
85747c6ae99SBarry Smith   PetscErrorCode ierr;
85847c6ae99SBarry Smith 
85947c6ae99SBarry Smith   PetscFunctionBegin;
860732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
86147c6ae99SBarry Smith   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
8624057135bSMatthew G Knepley   if (*dmf) {
86343842a1eSJed Brown     (*dmf)->ops->creatematrix = dm->ops->creatematrix;
864644e2e5bSBarry Smith     (*dmf)->ops->initialguess = dm->ops->initialguess;
865644e2e5bSBarry Smith     (*dmf)->ops->function     = dm->ops->function;
866644e2e5bSBarry Smith     (*dmf)->ops->functionj    = dm->ops->functionj;
867644e2e5bSBarry Smith     if (dm->ops->jacobian != DMComputeJacobianDefault) {
868644e2e5bSBarry Smith       (*dmf)->ops->jacobian     = dm->ops->jacobian;
869644e2e5bSBarry Smith     }
8708cd211a4SJed Brown     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
871644e2e5bSBarry Smith     (*dmf)->ctx     = dm->ctx;
872656b349aSBarry Smith     (*dmf)->levelup = dm->levelup + 1;
8734057135bSMatthew G Knepley   }
87447c6ae99SBarry Smith   PetscFunctionReturn(0);
87547c6ae99SBarry Smith }
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith #undef __FUNCT__
878eb3f98d2SBarry Smith #define __FUNCT__ "DMGetRefineLevel"
879eb3f98d2SBarry Smith /*@
880eb3f98d2SBarry Smith     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
881eb3f98d2SBarry Smith 
882eb3f98d2SBarry Smith     Not Collective
883eb3f98d2SBarry Smith 
884eb3f98d2SBarry Smith     Input Parameter:
885eb3f98d2SBarry Smith .   dm - the DM object
886eb3f98d2SBarry Smith 
887eb3f98d2SBarry Smith     Output Parameter:
888eb3f98d2SBarry Smith .   level - number of refinements
889eb3f98d2SBarry Smith 
890eb3f98d2SBarry Smith     Level: developer
891eb3f98d2SBarry Smith 
8926a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
893eb3f98d2SBarry Smith 
894eb3f98d2SBarry Smith @*/
895eb3f98d2SBarry Smith PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
896eb3f98d2SBarry Smith {
897eb3f98d2SBarry Smith   PetscFunctionBegin;
898eb3f98d2SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
899eb3f98d2SBarry Smith   *level = dm->levelup;
900eb3f98d2SBarry Smith   PetscFunctionReturn(0);
901eb3f98d2SBarry Smith }
902eb3f98d2SBarry Smith 
903eb3f98d2SBarry Smith #undef __FUNCT__
90447c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin"
90547c6ae99SBarry Smith /*@
90647c6ae99SBarry Smith     DMGlobalToLocalBegin - Begins updating local vectors from global vector
90747c6ae99SBarry Smith 
90847c6ae99SBarry Smith     Neighbor-wise Collective on DM
90947c6ae99SBarry Smith 
91047c6ae99SBarry Smith     Input Parameters:
91147c6ae99SBarry Smith +   dm - the DM object
91247c6ae99SBarry Smith .   g - the global vector
91347c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
91447c6ae99SBarry Smith -   l - the local vector
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith     Level: beginner
91847c6ae99SBarry Smith 
919e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
92047c6ae99SBarry Smith 
92147c6ae99SBarry Smith @*/
9227087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
92347c6ae99SBarry Smith {
92447c6ae99SBarry Smith   PetscErrorCode ierr;
92547c6ae99SBarry Smith 
92647c6ae99SBarry Smith   PetscFunctionBegin;
927171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
928843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
92947c6ae99SBarry Smith   PetscFunctionReturn(0);
93047c6ae99SBarry Smith }
93147c6ae99SBarry Smith 
93247c6ae99SBarry Smith #undef __FUNCT__
93347c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd"
93447c6ae99SBarry Smith /*@
93547c6ae99SBarry Smith     DMGlobalToLocalEnd - Ends updating local vectors from global vector
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith     Neighbor-wise Collective on DM
93847c6ae99SBarry Smith 
93947c6ae99SBarry Smith     Input Parameters:
94047c6ae99SBarry Smith +   dm - the DM object
94147c6ae99SBarry Smith .   g - the global vector
94247c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
94347c6ae99SBarry Smith -   l - the local vector
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith 
94647c6ae99SBarry Smith     Level: beginner
94747c6ae99SBarry Smith 
948e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith @*/
9517087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
95247c6ae99SBarry Smith {
95347c6ae99SBarry Smith   PetscErrorCode ierr;
95447c6ae99SBarry Smith 
95547c6ae99SBarry Smith   PetscFunctionBegin;
956171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
957843c4018SMatthew G Knepley   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
95847c6ae99SBarry Smith   PetscFunctionReturn(0);
95947c6ae99SBarry Smith }
96047c6ae99SBarry Smith 
96147c6ae99SBarry Smith #undef __FUNCT__
9629a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalBegin"
96347c6ae99SBarry Smith /*@
9649a42bb27SBarry Smith     DMLocalToGlobalBegin - updates global vectors from local vectors
9659a42bb27SBarry Smith 
9669a42bb27SBarry Smith     Neighbor-wise Collective on DM
9679a42bb27SBarry Smith 
9689a42bb27SBarry Smith     Input Parameters:
9699a42bb27SBarry Smith +   dm - the DM object
970f6813fd5SJed Brown .   l - the local vector
9719a42bb27SBarry 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
9729a42bb27SBarry Smith            base point.
973f6813fd5SJed Brown - - the global vector
9749a42bb27SBarry Smith 
9759a42bb27SBarry 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
9769a42bb27SBarry 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
9779a42bb27SBarry Smith            global array to the final global array with VecAXPY().
9789a42bb27SBarry Smith 
9799a42bb27SBarry Smith     Level: beginner
9809a42bb27SBarry Smith 
981e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
9829a42bb27SBarry Smith 
9839a42bb27SBarry Smith @*/
9847087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
9859a42bb27SBarry Smith {
9869a42bb27SBarry Smith   PetscErrorCode ierr;
9879a42bb27SBarry Smith 
9889a42bb27SBarry Smith   PetscFunctionBegin;
989171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
990843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
9919a42bb27SBarry Smith   PetscFunctionReturn(0);
9929a42bb27SBarry Smith }
9939a42bb27SBarry Smith 
9949a42bb27SBarry Smith #undef __FUNCT__
9959a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalEnd"
9969a42bb27SBarry Smith /*@
9979a42bb27SBarry Smith     DMLocalToGlobalEnd - updates global vectors from local vectors
99847c6ae99SBarry Smith 
99947c6ae99SBarry Smith     Neighbor-wise Collective on DM
100047c6ae99SBarry Smith 
100147c6ae99SBarry Smith     Input Parameters:
100247c6ae99SBarry Smith +   dm - the DM object
1003f6813fd5SJed Brown .   l - the local vector
100447c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
1005f6813fd5SJed Brown -   g - the global vector
100647c6ae99SBarry Smith 
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith     Level: beginner
100947c6ae99SBarry Smith 
1010e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith @*/
10137087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
101447c6ae99SBarry Smith {
101547c6ae99SBarry Smith   PetscErrorCode ierr;
101647c6ae99SBarry Smith 
101747c6ae99SBarry Smith   PetscFunctionBegin;
1018171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1019843c4018SMatthew G Knepley   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
102047c6ae99SBarry Smith   PetscFunctionReturn(0);
102147c6ae99SBarry Smith }
102247c6ae99SBarry Smith 
102347c6ae99SBarry Smith #undef __FUNCT__
102447c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobianDefault"
102547c6ae99SBarry Smith /*@
102647c6ae99SBarry Smith     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
102747c6ae99SBarry Smith 
102847c6ae99SBarry Smith     Collective on DM
102947c6ae99SBarry Smith 
103047c6ae99SBarry Smith     Input Parameter:
103147c6ae99SBarry Smith +   dm - the DM object
103247c6ae99SBarry Smith .   x - location to compute Jacobian at; may be ignored for linear problems
103347c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
103447c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
103547c6ae99SBarry Smith 
103647c6ae99SBarry Smith     Level: developer
103747c6ae99SBarry Smith 
1038e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
103947c6ae99SBarry Smith          DMSetFunction()
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith @*/
10427087cfbeSBarry Smith PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
104347c6ae99SBarry Smith {
104447c6ae99SBarry Smith   PetscErrorCode ierr;
1045171400e9SBarry Smith 
104647c6ae99SBarry Smith   PetscFunctionBegin;
1047171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
104847c6ae99SBarry Smith   *stflag = SAME_NONZERO_PATTERN;
104947c6ae99SBarry Smith   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
105047c6ae99SBarry Smith   if (A != B) {
105147c6ae99SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105247c6ae99SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105347c6ae99SBarry Smith   }
105447c6ae99SBarry Smith   PetscFunctionReturn(0);
105547c6ae99SBarry Smith }
105647c6ae99SBarry Smith 
105747c6ae99SBarry Smith #undef __FUNCT__
105847c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen"
105947c6ae99SBarry Smith /*@
106047c6ae99SBarry Smith     DMCoarsen - Coarsens a DM object
106147c6ae99SBarry Smith 
106247c6ae99SBarry Smith     Collective on DM
106347c6ae99SBarry Smith 
106447c6ae99SBarry Smith     Input Parameter:
106547c6ae99SBarry Smith +   dm - the DM object
106647c6ae99SBarry Smith -   comm - the communicator to contain the new DM object (or PETSC_NULL)
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith     Output Parameter:
106947c6ae99SBarry Smith .   dmc - the coarsened DM
107047c6ae99SBarry Smith 
107147c6ae99SBarry Smith     Level: developer
107247c6ae99SBarry Smith 
1073e727c939SJed Brown .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
107447c6ae99SBarry Smith 
107547c6ae99SBarry Smith @*/
10767087cfbeSBarry Smith PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
107747c6ae99SBarry Smith {
107847c6ae99SBarry Smith   PetscErrorCode ierr;
1079b17ce1afSJed Brown   DMCoarsenHookLink link;
108047c6ae99SBarry Smith 
108147c6ae99SBarry Smith   PetscFunctionBegin;
1082171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
108347c6ae99SBarry Smith   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
108443842a1eSJed Brown   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
108547c6ae99SBarry Smith   (*dmc)->ops->initialguess = dm->ops->initialguess;
108647c6ae99SBarry Smith   (*dmc)->ops->function     = dm->ops->function;
108747c6ae99SBarry Smith   (*dmc)->ops->functionj    = dm->ops->functionj;
108847c6ae99SBarry Smith   if (dm->ops->jacobian != DMComputeJacobianDefault) {
108947c6ae99SBarry Smith     (*dmc)->ops->jacobian     = dm->ops->jacobian;
109047c6ae99SBarry Smith   }
10918cd211a4SJed Brown   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1092644e2e5bSBarry Smith   (*dmc)->ctx       = dm->ctx;
1093656b349aSBarry Smith   (*dmc)->leveldown = dm->leveldown + 1;
1094b17ce1afSJed Brown   for (link=dm->coarsenhook; link; link=link->next) {
1095b17ce1afSJed Brown     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1096b17ce1afSJed Brown   }
1097b17ce1afSJed Brown   PetscFunctionReturn(0);
1098b17ce1afSJed Brown }
1099b17ce1afSJed Brown 
1100b17ce1afSJed Brown #undef __FUNCT__
1101b17ce1afSJed Brown #define __FUNCT__ "DMCoarsenHookAdd"
1102b17ce1afSJed Brown /*@
1103b17ce1afSJed Brown    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1104b17ce1afSJed Brown 
1105b17ce1afSJed Brown    Logically Collective
1106b17ce1afSJed Brown 
1107b17ce1afSJed Brown    Input Arguments:
1108b17ce1afSJed Brown +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1109b17ce1afSJed Brown .  coarsenhook - function to run when setting up a coarser level
1110b17ce1afSJed Brown .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
1111b17ce1afSJed Brown -  ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL)
1112b17ce1afSJed Brown 
1113b17ce1afSJed Brown    Calling sequence of coarsenhook:
1114b17ce1afSJed Brown $    coarsenhook(DM fine,DM coarse,void *ctx);
1115b17ce1afSJed Brown 
1116b17ce1afSJed Brown +  fine - fine level DM
1117b17ce1afSJed Brown .  coarse - coarse level DM to restrict problem to
1118b17ce1afSJed Brown -  ctx - optional user-defined function context
1119b17ce1afSJed Brown 
1120b17ce1afSJed Brown    Calling sequence for restricthook:
1121b17ce1afSJed Brown $    restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx)
1122b17ce1afSJed Brown 
1123b17ce1afSJed Brown +  fine - fine level DM
1124b17ce1afSJed Brown .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1125b17ce1afSJed Brown .  inject - matrix restricting by applying the transpose of injection
1126b17ce1afSJed Brown .  coarse - coarse level DM to update
1127b17ce1afSJed Brown -  ctx - optional user-defined function context
1128b17ce1afSJed Brown 
1129b17ce1afSJed Brown    Level: advanced
1130b17ce1afSJed Brown 
1131b17ce1afSJed Brown    Notes:
1132b17ce1afSJed Brown    This function is only needed if auxiliary data needs to be set up on coarse grids.
1133b17ce1afSJed Brown 
1134b17ce1afSJed Brown    If this function is called multiple times, the hooks will be run in the order they are added.
1135b17ce1afSJed Brown 
1136b17ce1afSJed Brown    The
1137b17ce1afSJed Brown 
1138b17ce1afSJed Brown    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1139b17ce1afSJed Brown    extract the finest level information from its context (instead of from the SNES).
1140b17ce1afSJed Brown 
1141b17ce1afSJed Brown .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1142b17ce1afSJed Brown @*/
1143b17ce1afSJed Brown PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1144b17ce1afSJed Brown {
1145b17ce1afSJed Brown   PetscErrorCode ierr;
1146b17ce1afSJed Brown   DMCoarsenHookLink link,*p;
1147b17ce1afSJed Brown 
1148b17ce1afSJed Brown   PetscFunctionBegin;
1149b17ce1afSJed Brown   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
11506bfea28cSJed Brown   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1151b17ce1afSJed Brown   ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1152b17ce1afSJed Brown   link->coarsenhook = coarsenhook;
1153b17ce1afSJed Brown   link->restricthook = restricthook;
1154b17ce1afSJed Brown   link->ctx = ctx;
11556cab3a1bSJed Brown   link->next = PETSC_NULL;
1156b17ce1afSJed Brown   *p = link;
1157b17ce1afSJed Brown   PetscFunctionReturn(0);
1158b17ce1afSJed Brown }
1159b17ce1afSJed Brown 
1160b17ce1afSJed Brown #undef __FUNCT__
1161b17ce1afSJed Brown #define __FUNCT__ "DMRestrict"
1162b17ce1afSJed Brown /*@
1163b17ce1afSJed Brown    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1164b17ce1afSJed Brown 
1165b17ce1afSJed Brown    Collective if any hooks are
1166b17ce1afSJed Brown 
1167b17ce1afSJed Brown    Input Arguments:
1168b17ce1afSJed Brown +  fine - finer DM to use as a base
1169b17ce1afSJed Brown .  restrct - restriction matrix, apply using MatRestrict()
1170b17ce1afSJed Brown .  inject - injection matrix, also use MatRestrict()
1171b17ce1afSJed Brown -  coarse - coarer DM to update
1172b17ce1afSJed Brown 
1173b17ce1afSJed Brown    Level: developer
1174b17ce1afSJed Brown 
1175b17ce1afSJed Brown .seealso: DMCoarsenHookAdd(), MatRestrict()
1176b17ce1afSJed Brown @*/
1177b17ce1afSJed Brown PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1178b17ce1afSJed Brown {
1179b17ce1afSJed Brown   PetscErrorCode ierr;
1180b17ce1afSJed Brown   DMCoarsenHookLink link;
1181b17ce1afSJed Brown 
1182b17ce1afSJed Brown   PetscFunctionBegin;
1183b17ce1afSJed Brown   for (link=fine->coarsenhook; link; link=link->next) {
1184b17ce1afSJed Brown     if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);}
1185b17ce1afSJed Brown   }
118647c6ae99SBarry Smith   PetscFunctionReturn(0);
118747c6ae99SBarry Smith }
118847c6ae99SBarry Smith 
118947c6ae99SBarry Smith #undef __FUNCT__
11905fe1f584SPeter Brune #define __FUNCT__ "DMGetCoarsenLevel"
11915fe1f584SPeter Brune /*@
11926a7d9d85SPeter Brune     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
11935fe1f584SPeter Brune 
11945fe1f584SPeter Brune     Not Collective
11955fe1f584SPeter Brune 
11965fe1f584SPeter Brune     Input Parameter:
11975fe1f584SPeter Brune .   dm - the DM object
11985fe1f584SPeter Brune 
11995fe1f584SPeter Brune     Output Parameter:
12006a7d9d85SPeter Brune .   level - number of coarsenings
12015fe1f584SPeter Brune 
12025fe1f584SPeter Brune     Level: developer
12035fe1f584SPeter Brune 
12046a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
12055fe1f584SPeter Brune 
12065fe1f584SPeter Brune @*/
12075fe1f584SPeter Brune PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
12085fe1f584SPeter Brune {
12095fe1f584SPeter Brune   PetscFunctionBegin;
12105fe1f584SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12115fe1f584SPeter Brune   *level = dm->leveldown;
12125fe1f584SPeter Brune   PetscFunctionReturn(0);
12135fe1f584SPeter Brune }
12145fe1f584SPeter Brune 
12155fe1f584SPeter Brune 
12165fe1f584SPeter Brune 
12175fe1f584SPeter Brune #undef __FUNCT__
121847c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy"
121947c6ae99SBarry Smith /*@C
122047c6ae99SBarry Smith     DMRefineHierarchy - Refines a DM object, all levels at once
122147c6ae99SBarry Smith 
122247c6ae99SBarry Smith     Collective on DM
122347c6ae99SBarry Smith 
122447c6ae99SBarry Smith     Input Parameter:
122547c6ae99SBarry Smith +   dm - the DM object
122647c6ae99SBarry Smith -   nlevels - the number of levels of refinement
122747c6ae99SBarry Smith 
122847c6ae99SBarry Smith     Output Parameter:
122947c6ae99SBarry Smith .   dmf - the refined DM hierarchy
123047c6ae99SBarry Smith 
123147c6ae99SBarry Smith     Level: developer
123247c6ae99SBarry Smith 
1233e727c939SJed Brown .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
123447c6ae99SBarry Smith 
123547c6ae99SBarry Smith @*/
12367087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
123747c6ae99SBarry Smith {
123847c6ae99SBarry Smith   PetscErrorCode ierr;
123947c6ae99SBarry Smith 
124047c6ae99SBarry Smith   PetscFunctionBegin;
1241171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
124247c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
124347c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
124447c6ae99SBarry Smith   if (dm->ops->refinehierarchy) {
124547c6ae99SBarry Smith     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
124647c6ae99SBarry Smith   } else if (dm->ops->refine) {
124747c6ae99SBarry Smith     PetscInt i;
124847c6ae99SBarry Smith 
124947c6ae99SBarry Smith     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
125047c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
125147c6ae99SBarry Smith       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
125247c6ae99SBarry Smith     }
125347c6ae99SBarry Smith   } else {
125447c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
125547c6ae99SBarry Smith   }
125647c6ae99SBarry Smith   PetscFunctionReturn(0);
125747c6ae99SBarry Smith }
125847c6ae99SBarry Smith 
125947c6ae99SBarry Smith #undef __FUNCT__
126047c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy"
126147c6ae99SBarry Smith /*@C
126247c6ae99SBarry Smith     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith     Collective on DM
126547c6ae99SBarry Smith 
126647c6ae99SBarry Smith     Input Parameter:
126747c6ae99SBarry Smith +   dm - the DM object
126847c6ae99SBarry Smith -   nlevels - the number of levels of coarsening
126947c6ae99SBarry Smith 
127047c6ae99SBarry Smith     Output Parameter:
127147c6ae99SBarry Smith .   dmc - the coarsened DM hierarchy
127247c6ae99SBarry Smith 
127347c6ae99SBarry Smith     Level: developer
127447c6ae99SBarry Smith 
1275e727c939SJed Brown .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
127647c6ae99SBarry Smith 
127747c6ae99SBarry Smith @*/
12787087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
127947c6ae99SBarry Smith {
128047c6ae99SBarry Smith   PetscErrorCode ierr;
128147c6ae99SBarry Smith 
128247c6ae99SBarry Smith   PetscFunctionBegin;
1283171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
128447c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
128547c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
128647c6ae99SBarry Smith   PetscValidPointer(dmc,3);
128747c6ae99SBarry Smith   if (dm->ops->coarsenhierarchy) {
128847c6ae99SBarry Smith     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
128947c6ae99SBarry Smith   } else if (dm->ops->coarsen) {
129047c6ae99SBarry Smith     PetscInt i;
129147c6ae99SBarry Smith 
129247c6ae99SBarry Smith     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
129347c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
129447c6ae99SBarry Smith       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
129547c6ae99SBarry Smith     }
129647c6ae99SBarry Smith   } else {
129747c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
129847c6ae99SBarry Smith   }
129947c6ae99SBarry Smith   PetscFunctionReturn(0);
130047c6ae99SBarry Smith }
130147c6ae99SBarry Smith 
130247c6ae99SBarry Smith #undef __FUNCT__
1303e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates"
130447c6ae99SBarry Smith /*@
1305e727c939SJed Brown    DMCreateAggregates - Gets the aggregates that map between
130647c6ae99SBarry Smith    grids associated with two DMs.
130747c6ae99SBarry Smith 
130847c6ae99SBarry Smith    Collective on DM
130947c6ae99SBarry Smith 
131047c6ae99SBarry Smith    Input Parameters:
131147c6ae99SBarry Smith +  dmc - the coarse grid DM
131247c6ae99SBarry Smith -  dmf - the fine grid DM
131347c6ae99SBarry Smith 
131447c6ae99SBarry Smith    Output Parameters:
131547c6ae99SBarry Smith .  rest - the restriction matrix (transpose of the projection matrix)
131647c6ae99SBarry Smith 
131747c6ae99SBarry Smith    Level: intermediate
131847c6ae99SBarry Smith 
131947c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid
132047c6ae99SBarry Smith 
1321e727c939SJed Brown .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
132247c6ae99SBarry Smith @*/
1323e727c939SJed Brown PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
132447c6ae99SBarry Smith {
132547c6ae99SBarry Smith   PetscErrorCode ierr;
132647c6ae99SBarry Smith 
132747c6ae99SBarry Smith   PetscFunctionBegin;
1328171400e9SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1329171400e9SBarry Smith   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
133047c6ae99SBarry Smith   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
133147c6ae99SBarry Smith   PetscFunctionReturn(0);
133247c6ae99SBarry Smith }
133347c6ae99SBarry Smith 
133447c6ae99SBarry Smith #undef __FUNCT__
13351a266240SBarry Smith #define __FUNCT__ "DMSetApplicationContextDestroy"
13361a266240SBarry Smith /*@C
13371a266240SBarry Smith     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
13381a266240SBarry Smith 
13391a266240SBarry Smith     Not Collective
13401a266240SBarry Smith 
13411a266240SBarry Smith     Input Parameters:
13421a266240SBarry Smith +   dm - the DM object
13431a266240SBarry Smith -   destroy - the destroy function
13441a266240SBarry Smith 
13451a266240SBarry Smith     Level: intermediate
13461a266240SBarry Smith 
1347e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
13481a266240SBarry Smith 
1349f07f9ceaSJed Brown @*/
13501a266240SBarry Smith PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
13511a266240SBarry Smith {
13521a266240SBarry Smith   PetscFunctionBegin;
1353171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13541a266240SBarry Smith   dm->ctxdestroy = destroy;
13551a266240SBarry Smith   PetscFunctionReturn(0);
13561a266240SBarry Smith }
13571a266240SBarry Smith 
13581a266240SBarry Smith #undef __FUNCT__
13591b2093e4SBarry Smith #define __FUNCT__ "DMSetApplicationContext"
1360b07ff414SBarry Smith /*@
13611b2093e4SBarry Smith     DMSetApplicationContext - Set a user context into a DM object
136247c6ae99SBarry Smith 
136347c6ae99SBarry Smith     Not Collective
136447c6ae99SBarry Smith 
136547c6ae99SBarry Smith     Input Parameters:
136647c6ae99SBarry Smith +   dm - the DM object
136747c6ae99SBarry Smith -   ctx - the user context
136847c6ae99SBarry Smith 
136947c6ae99SBarry Smith     Level: intermediate
137047c6ae99SBarry Smith 
1371e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith @*/
13741b2093e4SBarry Smith PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
137547c6ae99SBarry Smith {
137647c6ae99SBarry Smith   PetscFunctionBegin;
1377171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
137847c6ae99SBarry Smith   dm->ctx = ctx;
137947c6ae99SBarry Smith   PetscFunctionReturn(0);
138047c6ae99SBarry Smith }
138147c6ae99SBarry Smith 
138247c6ae99SBarry Smith #undef __FUNCT__
13831b2093e4SBarry Smith #define __FUNCT__ "DMGetApplicationContext"
138447c6ae99SBarry Smith /*@
13851b2093e4SBarry Smith     DMGetApplicationContext - Gets a user context from a DM object
138647c6ae99SBarry Smith 
138747c6ae99SBarry Smith     Not Collective
138847c6ae99SBarry Smith 
138947c6ae99SBarry Smith     Input Parameter:
139047c6ae99SBarry Smith .   dm - the DM object
139147c6ae99SBarry Smith 
139247c6ae99SBarry Smith     Output Parameter:
139347c6ae99SBarry Smith .   ctx - the user context
139447c6ae99SBarry Smith 
139547c6ae99SBarry Smith     Level: intermediate
139647c6ae99SBarry Smith 
1397e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
139847c6ae99SBarry Smith 
139947c6ae99SBarry Smith @*/
14001b2093e4SBarry Smith PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
140147c6ae99SBarry Smith {
140247c6ae99SBarry Smith   PetscFunctionBegin;
1403171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
14041b2093e4SBarry Smith   *(void**)ctx = dm->ctx;
140547c6ae99SBarry Smith   PetscFunctionReturn(0);
140647c6ae99SBarry Smith }
140747c6ae99SBarry Smith 
140847c6ae99SBarry Smith #undef __FUNCT__
140947c6ae99SBarry Smith #define __FUNCT__ "DMSetInitialGuess"
14107e833e3aSBarry Smith /*@C
141147c6ae99SBarry Smith     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
141247c6ae99SBarry Smith 
141347c6ae99SBarry Smith     Logically Collective on DM
141447c6ae99SBarry Smith 
141547c6ae99SBarry Smith     Input Parameter:
141647c6ae99SBarry Smith +   dm - the DM object to destroy
141747c6ae99SBarry Smith -   f - the function to compute the initial guess
141847c6ae99SBarry Smith 
141947c6ae99SBarry Smith     Level: intermediate
142047c6ae99SBarry Smith 
1421e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
142247c6ae99SBarry Smith 
1423f07f9ceaSJed Brown @*/
14247087cfbeSBarry Smith PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
142547c6ae99SBarry Smith {
142647c6ae99SBarry Smith   PetscFunctionBegin;
1427171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
142847c6ae99SBarry Smith   dm->ops->initialguess = f;
142947c6ae99SBarry Smith   PetscFunctionReturn(0);
143047c6ae99SBarry Smith }
143147c6ae99SBarry Smith 
143247c6ae99SBarry Smith #undef __FUNCT__
143347c6ae99SBarry Smith #define __FUNCT__ "DMSetFunction"
14347e833e3aSBarry Smith /*@C
143547c6ae99SBarry Smith     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
143647c6ae99SBarry Smith 
143747c6ae99SBarry Smith     Logically Collective on DM
143847c6ae99SBarry Smith 
143947c6ae99SBarry Smith     Input Parameter:
144047c6ae99SBarry Smith +   dm - the DM object
144147c6ae99SBarry Smith -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
144247c6ae99SBarry Smith 
144347c6ae99SBarry Smith     Level: intermediate
144447c6ae99SBarry Smith 
144547c6ae99SBarry Smith     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
144647c6ae99SBarry Smith            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
144747c6ae99SBarry Smith 
1448e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
144947c6ae99SBarry Smith          DMSetJacobian()
145047c6ae99SBarry Smith 
1451f07f9ceaSJed Brown @*/
14527087cfbeSBarry Smith PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
145347c6ae99SBarry Smith {
145447c6ae99SBarry Smith   PetscFunctionBegin;
1455171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
145647c6ae99SBarry Smith   dm->ops->function = f;
145747c6ae99SBarry Smith   if (f) {
145847c6ae99SBarry Smith     dm->ops->functionj = f;
145947c6ae99SBarry Smith   }
146047c6ae99SBarry Smith   PetscFunctionReturn(0);
146147c6ae99SBarry Smith }
146247c6ae99SBarry Smith 
146347c6ae99SBarry Smith #undef __FUNCT__
146447c6ae99SBarry Smith #define __FUNCT__ "DMSetJacobian"
14657e833e3aSBarry Smith /*@C
146647c6ae99SBarry Smith     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
146747c6ae99SBarry Smith 
146847c6ae99SBarry Smith     Logically Collective on DM
146947c6ae99SBarry Smith 
147047c6ae99SBarry Smith     Input Parameter:
147147c6ae99SBarry Smith +   dm - the DM object to destroy
147247c6ae99SBarry Smith -   f - the function to compute the matrix entries
147347c6ae99SBarry Smith 
147447c6ae99SBarry Smith     Level: intermediate
147547c6ae99SBarry Smith 
1476e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
147747c6ae99SBarry Smith          DMSetFunction()
147847c6ae99SBarry Smith 
1479f07f9ceaSJed Brown @*/
14807087cfbeSBarry Smith PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
148147c6ae99SBarry Smith {
148247c6ae99SBarry Smith   PetscFunctionBegin;
1483171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
148447c6ae99SBarry Smith   dm->ops->jacobian = f;
148547c6ae99SBarry Smith   PetscFunctionReturn(0);
148647c6ae99SBarry Smith }
148747c6ae99SBarry Smith 
148847c6ae99SBarry Smith #undef __FUNCT__
148908da532bSDmitry Karpeev #define __FUNCT__ "DMSetVariableBounds"
149008da532bSDmitry Karpeev /*@C
149108da532bSDmitry Karpeev     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
149208da532bSDmitry Karpeev 
149308da532bSDmitry Karpeev     Logically Collective on DM
149408da532bSDmitry Karpeev 
149508da532bSDmitry Karpeev     Input Parameter:
149608da532bSDmitry Karpeev +   dm - the DM object
149708da532bSDmitry Karpeev -   f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set)
149808da532bSDmitry Karpeev 
149908da532bSDmitry Karpeev     Level: intermediate
150008da532bSDmitry Karpeev 
1501e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
150208da532bSDmitry Karpeev          DMSetJacobian()
150308da532bSDmitry Karpeev 
150408da532bSDmitry Karpeev @*/
150508da532bSDmitry Karpeev PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
150608da532bSDmitry Karpeev {
150708da532bSDmitry Karpeev   PetscFunctionBegin;
150808da532bSDmitry Karpeev   dm->ops->computevariablebounds = f;
150908da532bSDmitry Karpeev   PetscFunctionReturn(0);
151008da532bSDmitry Karpeev }
151108da532bSDmitry Karpeev 
151208da532bSDmitry Karpeev #undef __FUNCT__
151308da532bSDmitry Karpeev #define __FUNCT__ "DMHasVariableBounds"
151408da532bSDmitry Karpeev /*@
151508da532bSDmitry Karpeev     DMHasVariableBounds - does the DM object have a variable bounds function?
151608da532bSDmitry Karpeev 
151708da532bSDmitry Karpeev     Not Collective
151808da532bSDmitry Karpeev 
151908da532bSDmitry Karpeev     Input Parameter:
152008da532bSDmitry Karpeev .   dm - the DM object to destroy
152108da532bSDmitry Karpeev 
152208da532bSDmitry Karpeev     Output Parameter:
152308da532bSDmitry Karpeev .   flg - PETSC_TRUE if the variable bounds function exists
152408da532bSDmitry Karpeev 
152508da532bSDmitry Karpeev     Level: developer
152608da532bSDmitry Karpeev 
1527e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
152808da532bSDmitry Karpeev 
152908da532bSDmitry Karpeev @*/
153008da532bSDmitry Karpeev PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
153108da532bSDmitry Karpeev {
153208da532bSDmitry Karpeev   PetscFunctionBegin;
153308da532bSDmitry Karpeev   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
153408da532bSDmitry Karpeev   PetscFunctionReturn(0);
153508da532bSDmitry Karpeev }
153608da532bSDmitry Karpeev 
153708da532bSDmitry Karpeev #undef __FUNCT__
153808da532bSDmitry Karpeev #define __FUNCT__ "DMComputeVariableBounds"
153908da532bSDmitry Karpeev /*@C
154008da532bSDmitry Karpeev     DMComputeVariableBounds - compute variable bounds used by SNESVI.
154108da532bSDmitry Karpeev 
154208da532bSDmitry Karpeev     Logically Collective on DM
154308da532bSDmitry Karpeev 
154408da532bSDmitry Karpeev     Input Parameters:
154508da532bSDmitry Karpeev +   dm - the DM object to destroy
154608da532bSDmitry Karpeev -   x  - current solution at which the bounds are computed
154708da532bSDmitry Karpeev 
154808da532bSDmitry Karpeev     Output parameters:
154908da532bSDmitry Karpeev +   xl - lower bound
155008da532bSDmitry Karpeev -   xu - upper bound
155108da532bSDmitry Karpeev 
155208da532bSDmitry Karpeev     Level: intermediate
155308da532bSDmitry Karpeev 
1554e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
155508da532bSDmitry Karpeev          DMSetFunction(), DMSetVariableBounds()
155608da532bSDmitry Karpeev 
155708da532bSDmitry Karpeev @*/
155808da532bSDmitry Karpeev PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
155908da532bSDmitry Karpeev {
156008da532bSDmitry Karpeev   PetscErrorCode ierr;
156108da532bSDmitry Karpeev   PetscFunctionBegin;
156208da532bSDmitry Karpeev   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
156308da532bSDmitry Karpeev   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
156408da532bSDmitry Karpeev   if(dm->ops->computevariablebounds) {
156508da532bSDmitry Karpeev     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr);
156608da532bSDmitry Karpeev   }
156708da532bSDmitry Karpeev   else {
156808da532bSDmitry Karpeev     ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr);
156908da532bSDmitry Karpeev     ierr = VecSet(xu,SNES_VI_INF);  CHKERRQ(ierr);
157008da532bSDmitry Karpeev   }
157108da532bSDmitry Karpeev   PetscFunctionReturn(0);
157208da532bSDmitry Karpeev }
157308da532bSDmitry Karpeev 
157408da532bSDmitry Karpeev #undef __FUNCT__
157547c6ae99SBarry Smith #define __FUNCT__ "DMComputeInitialGuess"
157647c6ae99SBarry Smith /*@
157747c6ae99SBarry Smith     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
157847c6ae99SBarry Smith 
157947c6ae99SBarry Smith     Collective on DM
158047c6ae99SBarry Smith 
158147c6ae99SBarry Smith     Input Parameter:
158247c6ae99SBarry Smith +   dm - the DM object to destroy
158347c6ae99SBarry Smith -   x - the vector to hold the initial guess values
158447c6ae99SBarry Smith 
158547c6ae99SBarry Smith     Level: developer
158647c6ae99SBarry Smith 
1587e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
158847c6ae99SBarry Smith 
158947c6ae99SBarry Smith @*/
15907087cfbeSBarry Smith PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
159147c6ae99SBarry Smith {
159247c6ae99SBarry Smith   PetscErrorCode ierr;
159347c6ae99SBarry Smith   PetscFunctionBegin;
1594171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
159547c6ae99SBarry Smith   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
159647c6ae99SBarry Smith   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
159747c6ae99SBarry Smith   PetscFunctionReturn(0);
159847c6ae99SBarry Smith }
159947c6ae99SBarry Smith 
160047c6ae99SBarry Smith #undef __FUNCT__
160147c6ae99SBarry Smith #define __FUNCT__ "DMHasInitialGuess"
160247c6ae99SBarry Smith /*@
160347c6ae99SBarry Smith     DMHasInitialGuess - does the DM object have an initial guess function
160447c6ae99SBarry Smith 
160547c6ae99SBarry Smith     Not Collective
160647c6ae99SBarry Smith 
160747c6ae99SBarry Smith     Input Parameter:
160847c6ae99SBarry Smith .   dm - the DM object to destroy
160947c6ae99SBarry Smith 
161047c6ae99SBarry Smith     Output Parameter:
161147c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
161247c6ae99SBarry Smith 
161347c6ae99SBarry Smith     Level: developer
161447c6ae99SBarry Smith 
1615e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
161647c6ae99SBarry Smith 
161747c6ae99SBarry Smith @*/
16187087cfbeSBarry Smith PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
161947c6ae99SBarry Smith {
162047c6ae99SBarry Smith   PetscFunctionBegin;
1621171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
162247c6ae99SBarry Smith   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
162347c6ae99SBarry Smith   PetscFunctionReturn(0);
162447c6ae99SBarry Smith }
162547c6ae99SBarry Smith 
162647c6ae99SBarry Smith #undef __FUNCT__
162747c6ae99SBarry Smith #define __FUNCT__ "DMHasFunction"
162847c6ae99SBarry Smith /*@
162947c6ae99SBarry Smith     DMHasFunction - does the DM object have a function
163047c6ae99SBarry Smith 
163147c6ae99SBarry Smith     Not Collective
163247c6ae99SBarry Smith 
163347c6ae99SBarry Smith     Input Parameter:
163447c6ae99SBarry Smith .   dm - the DM object to destroy
163547c6ae99SBarry Smith 
163647c6ae99SBarry Smith     Output Parameter:
163747c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
163847c6ae99SBarry Smith 
163947c6ae99SBarry Smith     Level: developer
164047c6ae99SBarry Smith 
1641e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
164247c6ae99SBarry Smith 
164347c6ae99SBarry Smith @*/
16447087cfbeSBarry Smith PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
164547c6ae99SBarry Smith {
164647c6ae99SBarry Smith   PetscFunctionBegin;
1647171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
164847c6ae99SBarry Smith   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
164947c6ae99SBarry Smith   PetscFunctionReturn(0);
165047c6ae99SBarry Smith }
165147c6ae99SBarry Smith 
165247c6ae99SBarry Smith #undef __FUNCT__
165347c6ae99SBarry Smith #define __FUNCT__ "DMHasJacobian"
165447c6ae99SBarry Smith /*@
165547c6ae99SBarry Smith     DMHasJacobian - does the DM object have a matrix function
165647c6ae99SBarry Smith 
165747c6ae99SBarry Smith     Not Collective
165847c6ae99SBarry Smith 
165947c6ae99SBarry Smith     Input Parameter:
166047c6ae99SBarry Smith .   dm - the DM object to destroy
166147c6ae99SBarry Smith 
166247c6ae99SBarry Smith     Output Parameter:
166347c6ae99SBarry Smith .   flg - PETSC_TRUE if function exists
166447c6ae99SBarry Smith 
166547c6ae99SBarry Smith     Level: developer
166647c6ae99SBarry Smith 
1667e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
166847c6ae99SBarry Smith 
166947c6ae99SBarry Smith @*/
16707087cfbeSBarry Smith PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
167147c6ae99SBarry Smith {
167247c6ae99SBarry Smith   PetscFunctionBegin;
1673171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
167447c6ae99SBarry Smith   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
167547c6ae99SBarry Smith   PetscFunctionReturn(0);
167647c6ae99SBarry Smith }
167747c6ae99SBarry Smith 
167847c6ae99SBarry Smith #undef  __FUNCT__
167908da532bSDmitry Karpeev #define __FUNCT__ "DMSetVec"
1680748fac09SDmitry Karpeev /*@C
168108da532bSDmitry Karpeev     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
168208da532bSDmitry Karpeev 
168308da532bSDmitry Karpeev     Collective on DM
168408da532bSDmitry Karpeev 
168508da532bSDmitry Karpeev     Input Parameter:
168608da532bSDmitry Karpeev +   dm - the DM object
1687e88d7f4bSDmitry Karpeev -   x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems.
168808da532bSDmitry Karpeev 
168908da532bSDmitry Karpeev     Level: developer
169008da532bSDmitry Karpeev 
1691e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
169208da532bSDmitry Karpeev          DMSetFunction(), DMSetJacobian(), DMSetVariableBounds()
169308da532bSDmitry Karpeev 
169408da532bSDmitry Karpeev @*/
169508da532bSDmitry Karpeev PetscErrorCode  DMSetVec(DM dm,Vec x)
169608da532bSDmitry Karpeev {
169708da532bSDmitry Karpeev   PetscErrorCode ierr;
169808da532bSDmitry Karpeev   PetscFunctionBegin;
169908da532bSDmitry Karpeev   if (x) {
170008da532bSDmitry Karpeev     if (!dm->x) {
170108da532bSDmitry Karpeev       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
170208da532bSDmitry Karpeev     }
170308da532bSDmitry Karpeev     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
170408da532bSDmitry Karpeev   }
170508da532bSDmitry Karpeev   else if(dm->x) {
170608da532bSDmitry Karpeev     ierr = VecDestroy(&dm->x);  CHKERRQ(ierr);
170708da532bSDmitry Karpeev   }
170808da532bSDmitry Karpeev   PetscFunctionReturn(0);
170908da532bSDmitry Karpeev }
171008da532bSDmitry Karpeev 
171108da532bSDmitry Karpeev 
171208da532bSDmitry Karpeev #undef __FUNCT__
171347c6ae99SBarry Smith #define __FUNCT__ "DMComputeFunction"
171447c6ae99SBarry Smith /*@
171547c6ae99SBarry Smith     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
171647c6ae99SBarry Smith 
171747c6ae99SBarry Smith     Collective on DM
171847c6ae99SBarry Smith 
171947c6ae99SBarry Smith     Input Parameter:
172047c6ae99SBarry Smith +   dm - the DM object to destroy
172147c6ae99SBarry Smith .   x - the location where the function is evaluationed, may be ignored for linear problems
172247c6ae99SBarry Smith -   b - the vector to hold the right hand side entries
172347c6ae99SBarry Smith 
172447c6ae99SBarry Smith     Level: developer
172547c6ae99SBarry Smith 
1726e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
172747c6ae99SBarry Smith          DMSetJacobian()
172847c6ae99SBarry Smith 
172947c6ae99SBarry Smith @*/
17307087cfbeSBarry Smith PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
173147c6ae99SBarry Smith {
173247c6ae99SBarry Smith   PetscErrorCode ierr;
173347c6ae99SBarry Smith   PetscFunctionBegin;
1734171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
173547c6ae99SBarry Smith   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1736644e2e5bSBarry Smith   PetscStackPush("DM user function");
173747c6ae99SBarry Smith   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1738644e2e5bSBarry Smith   PetscStackPop;
173947c6ae99SBarry Smith   PetscFunctionReturn(0);
174047c6ae99SBarry Smith }
174147c6ae99SBarry Smith 
174247c6ae99SBarry Smith 
174308da532bSDmitry Karpeev 
174447c6ae99SBarry Smith #undef __FUNCT__
174547c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobian"
174647c6ae99SBarry Smith /*@
174747c6ae99SBarry Smith     DMComputeJacobian - compute the matrix entries for the solver
174847c6ae99SBarry Smith 
174947c6ae99SBarry Smith     Collective on DM
175047c6ae99SBarry Smith 
175147c6ae99SBarry Smith     Input Parameter:
175247c6ae99SBarry Smith +   dm - the DM object
1753cab2e9ccSBarry Smith .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
175447c6ae99SBarry Smith .   A - matrix that defines the operator for the linear solve
175547c6ae99SBarry Smith -   B - the matrix used to construct the preconditioner
175647c6ae99SBarry Smith 
175747c6ae99SBarry Smith     Level: developer
175847c6ae99SBarry Smith 
1759e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
176047c6ae99SBarry Smith          DMSetFunction()
176147c6ae99SBarry Smith 
176247c6ae99SBarry Smith @*/
17637087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
176447c6ae99SBarry Smith {
176547c6ae99SBarry Smith   PetscErrorCode ierr;
176647c6ae99SBarry Smith 
176747c6ae99SBarry Smith   PetscFunctionBegin;
1768171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
176947c6ae99SBarry Smith   if (!dm->ops->jacobian) {
177047c6ae99SBarry Smith     ISColoring     coloring;
177147c6ae99SBarry Smith     MatFDColoring  fd;
177247c6ae99SBarry Smith 
1773171400e9SBarry Smith     ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
177447c6ae99SBarry Smith     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1775fcfd50ebSBarry Smith     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
177647c6ae99SBarry Smith     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
17770bdded8aSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
17780bdded8aSJed Brown     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
177971cd77b2SBarry Smith 
178047c6ae99SBarry Smith     dm->fd = fd;
178147c6ae99SBarry Smith     dm->ops->jacobian = DMComputeJacobianDefault;
17822533e041SBarry Smith 
178371cd77b2SBarry Smith     /* don't know why this is needed */
178471cd77b2SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
178547c6ae99SBarry Smith   }
178647c6ae99SBarry Smith   if (!x) x = dm->x;
178747c6ae99SBarry Smith   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1788cab2e9ccSBarry Smith 
178971cd77b2SBarry 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 */
1790649052a6SBarry Smith   if (x) {
1791cab2e9ccSBarry Smith     if (!dm->x) {
179271cd77b2SBarry Smith       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1793cab2e9ccSBarry Smith     }
1794cab2e9ccSBarry Smith     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1795649052a6SBarry Smith   }
1796a8248277SBarry Smith   if (A != B) {
1797a8248277SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1798a8248277SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1799a8248277SBarry Smith   }
180047c6ae99SBarry Smith   PetscFunctionReturn(0);
180147c6ae99SBarry Smith }
1802264ace61SBarry Smith 
1803cab2e9ccSBarry Smith 
1804264ace61SBarry Smith PetscFList DMList                       = PETSC_NULL;
1805264ace61SBarry Smith PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1806264ace61SBarry Smith 
1807264ace61SBarry Smith #undef __FUNCT__
1808264ace61SBarry Smith #define __FUNCT__ "DMSetType"
1809264ace61SBarry Smith /*@C
1810264ace61SBarry Smith   DMSetType - Builds a DM, for a particular DM implementation.
1811264ace61SBarry Smith 
1812264ace61SBarry Smith   Collective on DM
1813264ace61SBarry Smith 
1814264ace61SBarry Smith   Input Parameters:
1815264ace61SBarry Smith + dm     - The DM object
1816264ace61SBarry Smith - method - The name of the DM type
1817264ace61SBarry Smith 
1818264ace61SBarry Smith   Options Database Key:
1819264ace61SBarry Smith . -dm_type <type> - Sets the DM type; use -help for a list of available types
1820264ace61SBarry Smith 
1821264ace61SBarry Smith   Notes:
1822e1589f56SBarry Smith   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1823264ace61SBarry Smith 
1824264ace61SBarry Smith   Level: intermediate
1825264ace61SBarry Smith 
1826264ace61SBarry Smith .keywords: DM, set, type
1827264ace61SBarry Smith .seealso: DMGetType(), DMCreate()
1828264ace61SBarry Smith @*/
18297087cfbeSBarry Smith PetscErrorCode  DMSetType(DM dm, const DMType method)
1830264ace61SBarry Smith {
1831264ace61SBarry Smith   PetscErrorCode (*r)(DM);
1832264ace61SBarry Smith   PetscBool      match;
1833264ace61SBarry Smith   PetscErrorCode ierr;
1834264ace61SBarry Smith 
1835264ace61SBarry Smith   PetscFunctionBegin;
1836264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1837264ace61SBarry Smith   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1838264ace61SBarry Smith   if (match) PetscFunctionReturn(0);
1839264ace61SBarry Smith 
1840264ace61SBarry Smith   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
18414b91b6eaSBarry Smith   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1842264ace61SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1843264ace61SBarry Smith 
1844264ace61SBarry Smith   if (dm->ops->destroy) {
1845264ace61SBarry Smith     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1846b5c23020SJed Brown     dm->ops->destroy = PETSC_NULL;
1847264ace61SBarry Smith   }
1848264ace61SBarry Smith   ierr = (*r)(dm);CHKERRQ(ierr);
1849264ace61SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1850264ace61SBarry Smith   PetscFunctionReturn(0);
1851264ace61SBarry Smith }
1852264ace61SBarry Smith 
1853264ace61SBarry Smith #undef __FUNCT__
1854264ace61SBarry Smith #define __FUNCT__ "DMGetType"
1855264ace61SBarry Smith /*@C
1856264ace61SBarry Smith   DMGetType - Gets the DM type name (as a string) from the DM.
1857264ace61SBarry Smith 
1858264ace61SBarry Smith   Not Collective
1859264ace61SBarry Smith 
1860264ace61SBarry Smith   Input Parameter:
1861264ace61SBarry Smith . dm  - The DM
1862264ace61SBarry Smith 
1863264ace61SBarry Smith   Output Parameter:
1864264ace61SBarry Smith . type - The DM type name
1865264ace61SBarry Smith 
1866264ace61SBarry Smith   Level: intermediate
1867264ace61SBarry Smith 
1868264ace61SBarry Smith .keywords: DM, get, type, name
1869264ace61SBarry Smith .seealso: DMSetType(), DMCreate()
1870264ace61SBarry Smith @*/
18717087cfbeSBarry Smith PetscErrorCode  DMGetType(DM dm, const DMType *type)
1872264ace61SBarry Smith {
1873264ace61SBarry Smith   PetscErrorCode ierr;
1874264ace61SBarry Smith 
1875264ace61SBarry Smith   PetscFunctionBegin;
1876264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1877264ace61SBarry Smith   PetscValidCharPointer(type,2);
1878264ace61SBarry Smith   if (!DMRegisterAllCalled) {
1879264ace61SBarry Smith     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1880264ace61SBarry Smith   }
1881264ace61SBarry Smith   *type = ((PetscObject)dm)->type_name;
1882264ace61SBarry Smith   PetscFunctionReturn(0);
1883264ace61SBarry Smith }
1884264ace61SBarry Smith 
188567a56275SMatthew G Knepley #undef __FUNCT__
188667a56275SMatthew G Knepley #define __FUNCT__ "DMConvert"
188767a56275SMatthew G Knepley /*@C
188867a56275SMatthew G Knepley   DMConvert - Converts a DM to another DM, either of the same or different type.
188967a56275SMatthew G Knepley 
189067a56275SMatthew G Knepley   Collective on DM
189167a56275SMatthew G Knepley 
189267a56275SMatthew G Knepley   Input Parameters:
189367a56275SMatthew G Knepley + dm - the DM
189467a56275SMatthew G Knepley - newtype - new DM type (use "same" for the same type)
189567a56275SMatthew G Knepley 
189667a56275SMatthew G Knepley   Output Parameter:
189767a56275SMatthew G Knepley . M - pointer to new DM
189867a56275SMatthew G Knepley 
189967a56275SMatthew G Knepley   Notes:
190067a56275SMatthew G Knepley   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
190167a56275SMatthew G Knepley   the MPI communicator of the generated DM is always the same as the communicator
190267a56275SMatthew G Knepley   of the input DM.
190367a56275SMatthew G Knepley 
190467a56275SMatthew G Knepley   Level: intermediate
190567a56275SMatthew G Knepley 
190667a56275SMatthew G Knepley .seealso: DMCreate()
190767a56275SMatthew G Knepley @*/
190867a56275SMatthew G Knepley PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
190967a56275SMatthew G Knepley {
191067a56275SMatthew G Knepley   DM             B;
191167a56275SMatthew G Knepley   char           convname[256];
191267a56275SMatthew G Knepley   PetscBool      sametype, issame;
191367a56275SMatthew G Knepley   PetscErrorCode ierr;
191467a56275SMatthew G Knepley 
191567a56275SMatthew G Knepley   PetscFunctionBegin;
191667a56275SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
191767a56275SMatthew G Knepley   PetscValidType(dm,1);
191867a56275SMatthew G Knepley   PetscValidPointer(M,3);
191967a56275SMatthew G Knepley   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
192067a56275SMatthew G Knepley   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
192167a56275SMatthew G Knepley   {
192267a56275SMatthew G Knepley     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
192367a56275SMatthew G Knepley 
192467a56275SMatthew G Knepley     /*
192567a56275SMatthew G Knepley        Order of precedence:
192667a56275SMatthew G Knepley        1) See if a specialized converter is known to the current DM.
192767a56275SMatthew G Knepley        2) See if a specialized converter is known to the desired DM class.
192867a56275SMatthew G Knepley        3) See if a good general converter is registered for the desired class
192967a56275SMatthew G Knepley        4) See if a good general converter is known for the current matrix.
193067a56275SMatthew G Knepley        5) Use a really basic converter.
193167a56275SMatthew G Knepley     */
193267a56275SMatthew G Knepley 
193367a56275SMatthew G Knepley     /* 1) See if a specialized converter is known to the current DM and the desired class */
193467a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
193567a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
193667a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
193767a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
193867a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
193967a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
194067a56275SMatthew G Knepley     if (conv) goto foundconv;
194167a56275SMatthew G Knepley 
194267a56275SMatthew G Knepley     /* 2)  See if a specialized converter is known to the desired DM class. */
194367a56275SMatthew G Knepley     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
194467a56275SMatthew G Knepley     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
194567a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
194667a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
194767a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
194867a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
194967a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
195067a56275SMatthew G Knepley     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
195167a56275SMatthew G Knepley     if (conv) {
1952fcfd50ebSBarry Smith       ierr = DMDestroy(&B);CHKERRQ(ierr);
195367a56275SMatthew G Knepley       goto foundconv;
195467a56275SMatthew G Knepley     }
195567a56275SMatthew G Knepley 
195667a56275SMatthew G Knepley #if 0
195767a56275SMatthew G Knepley     /* 3) See if a good general converter is registered for the desired class */
195867a56275SMatthew G Knepley     conv = B->ops->convertfrom;
1959fcfd50ebSBarry Smith     ierr = DMDestroy(&B);CHKERRQ(ierr);
196067a56275SMatthew G Knepley     if (conv) goto foundconv;
196167a56275SMatthew G Knepley 
196267a56275SMatthew G Knepley     /* 4) See if a good general converter is known for the current matrix */
196367a56275SMatthew G Knepley     if (dm->ops->convert) {
196467a56275SMatthew G Knepley       conv = dm->ops->convert;
196567a56275SMatthew G Knepley     }
196667a56275SMatthew G Knepley     if (conv) goto foundconv;
196767a56275SMatthew G Knepley #endif
196867a56275SMatthew G Knepley 
196967a56275SMatthew G Knepley     /* 5) Use a really basic converter. */
197067a56275SMatthew G Knepley     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
197167a56275SMatthew G Knepley 
197267a56275SMatthew G Knepley     foundconv:
197367a56275SMatthew G Knepley     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
197467a56275SMatthew G Knepley     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
197567a56275SMatthew G Knepley     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
197667a56275SMatthew G Knepley   }
197767a56275SMatthew G Knepley   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
197867a56275SMatthew G Knepley   PetscFunctionReturn(0);
197967a56275SMatthew G Knepley }
1980264ace61SBarry Smith 
1981264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
1982264ace61SBarry Smith 
1983264ace61SBarry Smith #undef __FUNCT__
1984264ace61SBarry Smith #define __FUNCT__ "DMRegister"
1985264ace61SBarry Smith /*@C
1986264ace61SBarry Smith   DMRegister - See DMRegisterDynamic()
1987264ace61SBarry Smith 
1988264ace61SBarry Smith   Level: advanced
1989264ace61SBarry Smith @*/
19907087cfbeSBarry Smith PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1991264ace61SBarry Smith {
1992264ace61SBarry Smith   char fullname[PETSC_MAX_PATH_LEN];
1993264ace61SBarry Smith   PetscErrorCode ierr;
1994264ace61SBarry Smith 
1995264ace61SBarry Smith   PetscFunctionBegin;
1996264ace61SBarry Smith   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1997264ace61SBarry Smith   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1998264ace61SBarry Smith   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1999264ace61SBarry Smith   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
2000264ace61SBarry Smith   PetscFunctionReturn(0);
2001264ace61SBarry Smith }
2002264ace61SBarry Smith 
2003264ace61SBarry Smith 
2004264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
2005264ace61SBarry Smith #undef __FUNCT__
2006264ace61SBarry Smith #define __FUNCT__ "DMRegisterDestroy"
2007264ace61SBarry Smith /*@C
2008264ace61SBarry Smith    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
2009264ace61SBarry Smith 
2010264ace61SBarry Smith    Not Collective
2011264ace61SBarry Smith 
2012264ace61SBarry Smith    Level: advanced
2013264ace61SBarry Smith 
2014264ace61SBarry Smith .keywords: DM, register, destroy
2015264ace61SBarry Smith .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
2016264ace61SBarry Smith @*/
20177087cfbeSBarry Smith PetscErrorCode  DMRegisterDestroy(void)
2018264ace61SBarry Smith {
2019264ace61SBarry Smith   PetscErrorCode ierr;
2020264ace61SBarry Smith 
2021264ace61SBarry Smith   PetscFunctionBegin;
2022264ace61SBarry Smith   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
2023264ace61SBarry Smith   DMRegisterAllCalled = PETSC_FALSE;
2024264ace61SBarry Smith   PetscFunctionReturn(0);
2025264ace61SBarry Smith }
202623f975d1SBarry Smith 
202723f975d1SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
2028c6db04a5SJed Brown #include <mex.h>
202923f975d1SBarry Smith 
20303014e516SBarry Smith typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
203123f975d1SBarry Smith 
203223f975d1SBarry Smith #undef __FUNCT__
203323f975d1SBarry Smith #define __FUNCT__ "DMComputeFunction_Matlab"
203423f975d1SBarry Smith /*
203523f975d1SBarry Smith    DMComputeFunction_Matlab - Calls the function that has been set with
203623f975d1SBarry Smith                          DMSetFunctionMatlab().
203723f975d1SBarry Smith 
203823f975d1SBarry Smith    For linear problems x is null
203923f975d1SBarry Smith 
204023f975d1SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
204123f975d1SBarry Smith */
20427087cfbeSBarry Smith PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
204323f975d1SBarry Smith {
204423f975d1SBarry Smith   PetscErrorCode    ierr;
204523f975d1SBarry Smith   DMMatlabContext   *sctx;
204623f975d1SBarry Smith   int               nlhs = 1,nrhs = 4;
204723f975d1SBarry Smith   mxArray	    *plhs[1],*prhs[4];
204823f975d1SBarry Smith   long long int     lx = 0,ly = 0,ls = 0;
204923f975d1SBarry Smith 
205023f975d1SBarry Smith   PetscFunctionBegin;
205123f975d1SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
205223f975d1SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
205323f975d1SBarry Smith   PetscCheckSameComm(dm,1,y,3);
205423f975d1SBarry Smith 
205523f975d1SBarry Smith   /* call Matlab function in ctx with arguments x and y */
20561b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
205723f975d1SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
205823f975d1SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
20593014e516SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
206023f975d1SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
206123f975d1SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
206223f975d1SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)ly);
206323f975d1SBarry Smith   prhs[3] =  mxCreateString(sctx->funcname);
2064b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
206523f975d1SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
206623f975d1SBarry Smith   mxDestroyArray(prhs[0]);
206723f975d1SBarry Smith   mxDestroyArray(prhs[1]);
206823f975d1SBarry Smith   mxDestroyArray(prhs[2]);
206923f975d1SBarry Smith   mxDestroyArray(prhs[3]);
207023f975d1SBarry Smith   mxDestroyArray(plhs[0]);
207123f975d1SBarry Smith   PetscFunctionReturn(0);
207223f975d1SBarry Smith }
207323f975d1SBarry Smith 
207423f975d1SBarry Smith 
207523f975d1SBarry Smith #undef __FUNCT__
207623f975d1SBarry Smith #define __FUNCT__ "DMSetFunctionMatlab"
207723f975d1SBarry Smith /*
207823f975d1SBarry Smith    DMSetFunctionMatlab - Sets the function evaluation routine
207923f975d1SBarry Smith 
208023f975d1SBarry Smith */
20817087cfbeSBarry Smith PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
208223f975d1SBarry Smith {
208323f975d1SBarry Smith   PetscErrorCode    ierr;
208423f975d1SBarry Smith   DMMatlabContext   *sctx;
208523f975d1SBarry Smith 
208623f975d1SBarry Smith   PetscFunctionBegin;
2087171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
208823f975d1SBarry Smith   /* currently sctx is memory bleed */
20891b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
20903014e516SBarry Smith   if (!sctx) {
209123f975d1SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
20923014e516SBarry Smith   }
209323f975d1SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
20941b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
209523f975d1SBarry Smith   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
209623f975d1SBarry Smith   PetscFunctionReturn(0);
209723f975d1SBarry Smith }
20983014e516SBarry Smith 
20993014e516SBarry Smith #undef __FUNCT__
21003014e516SBarry Smith #define __FUNCT__ "DMComputeJacobian_Matlab"
21013014e516SBarry Smith /*
21023014e516SBarry Smith    DMComputeJacobian_Matlab - Calls the function that has been set with
21033014e516SBarry Smith                          DMSetJacobianMatlab().
21043014e516SBarry Smith 
21053014e516SBarry Smith    For linear problems x is null
21063014e516SBarry Smith 
21073014e516SBarry Smith .seealso: DMSetFunction(), DMGetFunction()
21083014e516SBarry Smith */
21097087cfbeSBarry Smith PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
21103014e516SBarry Smith {
21113014e516SBarry Smith   PetscErrorCode    ierr;
21123014e516SBarry Smith   DMMatlabContext   *sctx;
21133014e516SBarry Smith   int               nlhs = 2,nrhs = 5;
21143014e516SBarry Smith   mxArray	    *plhs[2],*prhs[5];
21153014e516SBarry Smith   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
21163014e516SBarry Smith 
21173014e516SBarry Smith   PetscFunctionBegin;
21183014e516SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
21193014e516SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
21203014e516SBarry Smith 
2121e3c5b3baSBarry Smith   /* call MATLAB function in ctx with arguments x, A, and B */
21221b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
21233014e516SBarry Smith   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
21243014e516SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
21253014e516SBarry Smith   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
21263014e516SBarry Smith   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
21273014e516SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
21283014e516SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
21293014e516SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lA);
21303014e516SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lB);
21313014e516SBarry Smith   prhs[4] =  mxCreateString(sctx->jacname);
2132b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
2133c980e822SBarry Smith   *str    =  (MatStructure) mxGetScalar(plhs[0]);
2134c088a8dcSBarry Smith   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
21353014e516SBarry Smith   mxDestroyArray(prhs[0]);
21363014e516SBarry Smith   mxDestroyArray(prhs[1]);
21373014e516SBarry Smith   mxDestroyArray(prhs[2]);
21383014e516SBarry Smith   mxDestroyArray(prhs[3]);
21393014e516SBarry Smith   mxDestroyArray(prhs[4]);
21403014e516SBarry Smith   mxDestroyArray(plhs[0]);
21413014e516SBarry Smith   mxDestroyArray(plhs[1]);
21423014e516SBarry Smith   PetscFunctionReturn(0);
21433014e516SBarry Smith }
21443014e516SBarry Smith 
21453014e516SBarry Smith 
21463014e516SBarry Smith #undef __FUNCT__
21473014e516SBarry Smith #define __FUNCT__ "DMSetJacobianMatlab"
21483014e516SBarry Smith /*
21493014e516SBarry Smith    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
21503014e516SBarry Smith 
21513014e516SBarry Smith */
21527087cfbeSBarry Smith PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
21533014e516SBarry Smith {
21543014e516SBarry Smith   PetscErrorCode    ierr;
21553014e516SBarry Smith   DMMatlabContext   *sctx;
21563014e516SBarry Smith 
21573014e516SBarry Smith   PetscFunctionBegin;
2158171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
21593014e516SBarry Smith   /* currently sctx is memory bleed */
21601b2093e4SBarry Smith   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
21613014e516SBarry Smith   if (!sctx) {
21623014e516SBarry Smith     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
21633014e516SBarry Smith   }
21643014e516SBarry Smith   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
21651b2093e4SBarry Smith   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
21663014e516SBarry Smith   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
21673014e516SBarry Smith   PetscFunctionReturn(0);
21683014e516SBarry Smith }
216923f975d1SBarry Smith #endif
2170b859378eSBarry Smith 
2171b859378eSBarry Smith #undef __FUNCT__
2172b859378eSBarry Smith #define __FUNCT__ "DMLoad"
2173b859378eSBarry Smith /*@C
2174b859378eSBarry Smith   DMLoad - Loads a DM that has been stored in binary or HDF5 format
2175b859378eSBarry Smith   with DMView().
2176b859378eSBarry Smith 
2177b859378eSBarry Smith   Collective on PetscViewer
2178b859378eSBarry Smith 
2179b859378eSBarry Smith   Input Parameters:
2180b859378eSBarry Smith + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2181b859378eSBarry Smith            some related function before a call to DMLoad().
2182b859378eSBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2183b859378eSBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2184b859378eSBarry Smith 
2185b859378eSBarry Smith    Level: intermediate
2186b859378eSBarry Smith 
2187b859378eSBarry Smith   Notes:
2188b859378eSBarry Smith   Defaults to the DM DA.
2189b859378eSBarry Smith 
2190b859378eSBarry Smith   Notes for advanced users:
2191b859378eSBarry Smith   Most users should not need to know the details of the binary storage
2192b859378eSBarry Smith   format, since DMLoad() and DMView() completely hide these details.
2193b859378eSBarry Smith   But for anyone who's interested, the standard binary matrix storage
2194b859378eSBarry Smith   format is
2195b859378eSBarry Smith .vb
2196b859378eSBarry Smith      has not yet been determined
2197b859378eSBarry Smith .ve
2198b859378eSBarry Smith 
2199b859378eSBarry Smith    In addition, PETSc automatically does the byte swapping for
2200b859378eSBarry Smith machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
2201b859378eSBarry Smith linux, Windows and the paragon; thus if you write your own binary
2202b859378eSBarry Smith read/write routines you have to swap the bytes; see PetscBinaryRead()
2203b859378eSBarry Smith and PetscBinaryWrite() to see how this may be done.
2204b859378eSBarry Smith 
2205b859378eSBarry Smith   Concepts: vector^loading from file
2206b859378eSBarry Smith 
2207b859378eSBarry Smith .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2208b859378eSBarry Smith @*/
2209b859378eSBarry Smith PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2210b859378eSBarry Smith {
2211b859378eSBarry Smith   PetscErrorCode ierr;
2212b859378eSBarry Smith 
2213b859378eSBarry Smith   PetscFunctionBegin;
2214b859378eSBarry Smith   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2215b859378eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2216b859378eSBarry Smith 
2217b859378eSBarry Smith   if (!((PetscObject)newdm)->type_name) {
2218b859378eSBarry Smith     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
2219b859378eSBarry Smith   }
2220b859378eSBarry Smith   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2221b859378eSBarry Smith   PetscFunctionReturn(0);
2222b859378eSBarry Smith }
2223b859378eSBarry Smith 
22247da65231SMatthew G Knepley /******************************** FEM Support **********************************/
22257da65231SMatthew G Knepley 
22267da65231SMatthew G Knepley #undef __FUNCT__
22277da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellVector"
22287da65231SMatthew G Knepley PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) {
22291d47ebbbSSatish Balay   PetscInt       f;
22301b30c384SMatthew G Knepley   PetscErrorCode ierr;
22311b30c384SMatthew G Knepley 
22327da65231SMatthew G Knepley   PetscFunctionBegin;
223374778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
22341d47ebbbSSatish Balay   for(f = 0; f < len; ++f) {
223574778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
22367da65231SMatthew G Knepley   }
22377da65231SMatthew G Knepley   PetscFunctionReturn(0);
22387da65231SMatthew G Knepley }
22397da65231SMatthew G Knepley 
22407da65231SMatthew G Knepley #undef __FUNCT__
22417da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellMatrix"
22427da65231SMatthew G Knepley PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) {
22431b30c384SMatthew G Knepley   PetscInt       f, g;
22447da65231SMatthew G Knepley   PetscErrorCode ierr;
22457da65231SMatthew G Knepley 
22467da65231SMatthew G Knepley   PetscFunctionBegin;
224774778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
22481d47ebbbSSatish Balay   for(f = 0; f < rows; ++f) {
224974778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
22501d47ebbbSSatish Balay     for(g = 0; g < cols; ++g) {
225174778d6cSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
22527da65231SMatthew G Knepley     }
225374778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
22547da65231SMatthew G Knepley   }
22557da65231SMatthew G Knepley   PetscFunctionReturn(0);
22567da65231SMatthew G Knepley }
2257e7c4fc90SDmitry Karpeev 
2258