Actual source code: dmimpl.h

petsc-3.8.4 2018-03-24
Report Typos and Errors


  3: #if !defined(_DMIMPL_H)
  4: #define _DMIMPL_H

  6:  #include <petscdm.h>
  7:  #include <petsc/private/petscimpl.h>
  8:  #include <petsc/private/petscdsimpl.h>

 10: PETSC_EXTERN PetscBool DMRegisterAllCalled;
 11: PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
 12: typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt field, MatNullSpace *nullSpace);

 14: typedef struct _DMOps *DMOps;
 15: struct _DMOps {
 16:   PetscErrorCode (*view)(DM,PetscViewer);
 17:   PetscErrorCode (*load)(DM,PetscViewer);
 18:   PetscErrorCode (*clone)(DM,DM*);
 19:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,DM);
 20:   PetscErrorCode (*setup)(DM);
 21:   PetscErrorCode (*createdefaultsection)(DM);
 22:   PetscErrorCode (*createdefaultconstraints)(DM);
 23:   PetscErrorCode (*createglobalvector)(DM,Vec*);
 24:   PetscErrorCode (*createlocalvector)(DM,Vec*);
 25:   PetscErrorCode (*getlocaltoglobalmapping)(DM);
 26:   PetscErrorCode (*createfieldis)(DM,PetscInt*,char***,IS**);
 27:   PetscErrorCode (*createcoordinatedm)(DM,DM*);

 29:   PetscErrorCode (*getcoloring)(DM,ISColoringType,ISColoring*);
 30:   PetscErrorCode (*creatematrix)(DM, Mat*);
 31:   PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*);
 32:   PetscErrorCode (*createrestriction)(DM,DM,Mat*);
 33:   PetscErrorCode (*getaggregates)(DM,DM,Mat*);
 34:   PetscErrorCode (*getinjection)(DM,DM,Mat*);

 36:   PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
 37:   PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*);
 38:   PetscErrorCode (*refinehierarchy)(DM,PetscInt,DM*);
 39:   PetscErrorCode (*coarsenhierarchy)(DM,PetscInt,DM*);
 40:   PetscErrorCode (*adaptlabel)(DM,DMLabel,DM*);
 41:   PetscErrorCode (*adaptmetric)(DM,Vec,DMLabel,DM*);

 43:   PetscErrorCode (*globaltolocalbegin)(DM,Vec,InsertMode,Vec);
 44:   PetscErrorCode (*globaltolocalend)(DM,Vec,InsertMode,Vec);
 45:   PetscErrorCode (*localtoglobalbegin)(DM,Vec,InsertMode,Vec);
 46:   PetscErrorCode (*localtoglobalend)(DM,Vec,InsertMode,Vec);
 47:   PetscErrorCode (*localtolocalbegin)(DM,Vec,InsertMode,Vec);
 48:   PetscErrorCode (*localtolocalend)(DM,Vec,InsertMode,Vec);

 50:   PetscErrorCode (*destroy)(DM);

 52:   PetscErrorCode (*computevariablebounds)(DM,Vec,Vec);

 54:   PetscErrorCode (*createsubdm)(DM,PetscInt,PetscInt*,IS*,DM*);
 55:   PetscErrorCode (*createfielddecomposition)(DM,PetscInt*,char***,IS**,DM**);
 56:   PetscErrorCode (*createdomaindecomposition)(DM,PetscInt*,char***,IS**,IS**,DM**);
 57:   PetscErrorCode (*createddscatters)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

 59:   PetscErrorCode (*getdimpoints)(DM,PetscInt,PetscInt*,PetscInt*);
 60:   PetscErrorCode (*locatepoints)(DM,Vec,DMPointLocationType,PetscSF);
 61:   PetscErrorCode (*getneighbors)(DM,PetscInt*,const PetscMPIInt**);
 62:   PetscErrorCode (*getboundingbox)(DM,PetscReal*,PetscReal*);
 63:   PetscErrorCode (*getlocalboundingbox)(DM,PetscReal*,PetscReal*);
 64:   PetscErrorCode (*locatepointssubdomain)(DM,Vec,PetscMPIInt**);

 66:   PetscErrorCode (*projectfunctionlocal)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
 67:   PetscErrorCode (*projectfunctionlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
 68:   PetscErrorCode (*projectfieldlocal)(DM,PetscReal,Vec,void(**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
 69:   PetscErrorCode (*projectfieldlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Vec,void(**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
 70:   PetscErrorCode (*computel2diff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
 71:   PetscErrorCode (*computel2gradientdiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal [],const PetscReal[],PetscInt, PetscScalar *,void *),void **,Vec,const PetscReal[],PetscReal *);
 72:   PetscErrorCode (*computel2fielddiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
 73: };

 75: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
 76: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
 77: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);

 79: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
 80: struct _DMCoarsenHookLink {
 81:   PetscErrorCode (*coarsenhook)(DM,DM,void*);              /* Run once, when coarse DM is created */
 82:   PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*); /* Run each time a new problem is restricted to a coarse grid */
 83:   void *ctx;
 84:   DMCoarsenHookLink next;
 85: };

 87: typedef struct _DMRefineHookLink *DMRefineHookLink;
 88: struct _DMRefineHookLink {
 89:   PetscErrorCode (*refinehook)(DM,DM,void*);     /* Run once, when a fine DM is created */
 90:   PetscErrorCode (*interphook)(DM,Mat,DM,void*); /* Run each time a new problem is interpolated to a fine grid */
 91:   void *ctx;
 92:   DMRefineHookLink next;
 93: };

 95: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
 96: struct _DMSubDomainHookLink {
 97:   PetscErrorCode (*ddhook)(DM,DM,void*);
 98:   PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*);
 99:   void *ctx;
100:   DMSubDomainHookLink next;
101: };

103: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
104: struct _DMGlobalToLocalHookLink {
105:   PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
106:   PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
107:   void *ctx;
108:   DMGlobalToLocalHookLink next;
109: };

111: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
112: struct _DMLocalToGlobalHookLink {
113:   PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
114:   PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
115:   void *ctx;
116:   DMLocalToGlobalHookLink next;
117: };

119: typedef enum {DMVEC_STATUS_IN,DMVEC_STATUS_OUT} DMVecStatus;
120: typedef struct _DMNamedVecLink *DMNamedVecLink;
121: struct _DMNamedVecLink {
122:   Vec X;
123:   char *name;
124:   DMVecStatus status;
125:   DMNamedVecLink next;
126: };

128: typedef struct _DMWorkLink *DMWorkLink;
129: struct _DMWorkLink {
130:   size_t     bytes;
131:   void       *mem;
132:   DMWorkLink next;
133: };

135: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users  via DMGetGlobalVector(), DMGetLocalVector() */

137: struct _n_DMLabelLink {
138:   DMLabel              label;
139:   PetscBool            output;
140:   struct _n_DMLabelLink *next;
141: };
142: typedef struct _n_DMLabelLink *DMLabelLink;

144: struct _n_DMLabelLinkList {
145:   PetscInt refct;
146:   DMLabelLink next;
147: };
148: typedef struct _n_DMLabelLinkList *DMLabelLinkList;

150: typedef struct _n_Boundary *DMBoundary;

152: struct _n_Boundary {
153:   DSBoundary  dsboundary;
154:   DMLabel     label;
155:   DMBoundary  next;
156: };

158: PETSC_EXTERN PetscErrorCode DMDestroyLabelLinkList(DM);

160: struct _p_DM {
161:   PETSCHEADER(struct _DMOps);
162:   Vec                     localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
163:   Vec                     globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
164:   DMNamedVecLink          namedglobal;
165:   DMNamedVecLink          namedlocal;
166:   DMWorkLink              workin,workout;
167:   DMLabelLinkList         labels;            /* Linked list of labels */
168:   DMLabel                 depthLabel;        /* Optimized access to depth label */
169:   void                    *ctx;    /* a user context */
170:   PetscErrorCode          (*ctxdestroy)(void**);
171:   Vec                     x;       /* location at which the functions/Jacobian are computed */
172:   ISColoringType          coloringtype;
173:   MatFDColoring           fd;
174:   VecType                 vectype;  /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
175:   MatType                 mattype;  /* type of matrix created with DMCreateMatrix() */
176:   PetscInt                bs;
177:   ISLocalToGlobalMapping  ltogmap;
178:   PetscBool               prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
179:   PetscBool               structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
180:   PetscInt                levelup,leveldown;  /* if the DM has been obtained by refining (or coarsening) this indicates how many times that process has been used to generate this DM */
181:   PetscBool               setupcalled;        /* Indicates that the DM has been set up, methods that modify a DM such that a fresh setup is required should reset this flag */
182:   void                    *data;
183:   /* Hierarchy / Submeshes */
184:   DM                      coarseMesh;
185:   DM                      fineMesh;
186:   DMCoarsenHookLink       coarsenhook; /* For transfering auxiliary problem data to coarser grids */
187:   DMRefineHookLink        refinehook;
188:   DMSubDomainHookLink     subdomainhook;
189:   DMGlobalToLocalHookLink gtolhook;
190:   DMLocalToGlobalHookLink ltoghook;
191:   /* Topology */
192:   PetscInt                dim;                  /* The topological dimension */
193:   /* Flexible communication */
194:   PetscSF                 sf;                   /* SF for parallel point overlap */
195:   PetscSF                 defaultSF;            /* SF for parallel dof overlap using default section */
196:   PetscSF                 sfNatural;            /* SF mapping to the "natural" ordering */
197:   PetscBool               useNatural;           /* Create the natural SF */
198:   /* Allows a non-standard data layout */
199:   PetscSection            defaultSection;       /* Layout for local vectors */
200:   PetscSection            defaultGlobalSection; /* Layout for global vectors */
201:   PetscLayout             map;
202:   /* Constraints */
203:   PetscSection            defaultConstraintSection;
204:   Mat                     defaultConstraintMat;
205:   /* Coordinates */
206:   PetscInt                dimEmbed;             /* The dimension of the embedding space */
207:   DM                      coordinateDM;         /* Layout for coordinates (default section) */
208:   Vec                     coordinates;          /* Coordinate values in global vector */
209:   Vec                     coordinatesLocal;     /* Coordinate values in local  vector */
210:   PetscBool               periodic;             /* Is the DM periodic? */
211:   PetscReal              *L, *maxCell;          /* Size of periodic box and max cell size for determining periodicity */
212:   DMBoundaryType         *bdtype;               /* Indicates type of topological boundary */
213:   /* Null spaces -- of course I should make this have a variable number of fields */
214:   /*   I now believe this might not be the right way: see below */
215:   NullSpaceFunc           nullspaceConstructors[10];
216:   /* Fields are represented by objects */
217:   PetscDS                 prob;
218:   DMBoundary              boundary;          /* List of boundary conditions */
219:   /* Output structures */
220:   DM                      dmBC;                 /* The DM with boundary conditions in the global DM */
221:   PetscInt                outputSequenceNum;    /* The current sequence number for output */
222:   PetscReal               outputSequenceVal;    /* The current sequence value for output */

224:   PetscObject             dmksp,dmsnes,dmts;
225: };

227: PETSC_EXTERN PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocatePoints, DM_Coarsen, DM_Refine, DM_CreateInterpolation, DM_CreateRestriction;

229: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
230: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);
231: PETSC_EXTERN PetscErrorCode DMCreateSubDM_Section_Private(DM,PetscInt,PetscInt[],IS*,DM*);

233: /*

235:           Composite Vectors

237:       Single global representation
238:       Individual global representations
239:       Single local representation
240:       Individual local representations

242:       Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????

244:        DM da_u, da_v, da_p

246:        DM dm_velocities

248:        DM dm

250:        DMDACreate(,&da_u);
251:        DMDACreate(,&da_v);
252:        DMCompositeCreate(,&dm_velocities);
253:        DMCompositeAddDM(dm_velocities,(DM)du);
254:        DMCompositeAddDM(dm_velocities,(DM)dv);

256:        DMDACreate(,&da_p);
257:        DMCompositeCreate(,&dm_velocities);
258:        DMCompositeAddDM(dm,(DM)dm_velocities);
259:        DMCompositeAddDM(dm,(DM)dm_p);


262:     Access parts of composite vectors (Composite only)
263:     ---------------------------------
264:       DMCompositeGetAccess  - access the global vector as subvectors and array (for redundant arrays)
265:       ADD for local vector -

267:     Element access
268:     --------------
269:       From global vectors
270:          -DAVecGetArray   - for DMDA
271:          -VecGetArray - for DMSliced
272:          ADD for DMComposite???  maybe

274:       From individual vector
275:           -DAVecGetArray - for DMDA
276:           -VecGetArray -for sliced
277:          ADD for DMComposite??? maybe

279:       From single local vector
280:           ADD         * single local vector as arrays?

282:    Communication
283:    -------------
284:       DMGlobalToLocal - global vector to single local vector

286:       DMCompositeScatter/Gather - direct to individual local vectors and arrays   CHANGE name closer to GlobalToLocal?

288:    Obtaining vectors
289:    -----------------
290:       DMCreateGlobal/Local
291:       DMGetGlobal/Local
292:       DMCompositeGetLocalVectors   - gives individual local work vectors and arrays


295: ?????   individual global vectors   ????

297: */

299: #if defined(PETSC_HAVE_HDF5)
300: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
301: #endif

303: #endif