Actual source code: petscda.h
1: /*
2: Regular array object, for easy parallelism of simple grid
3: problems on regular distributed arrays.
4: */
7: #include petscvec.h
8: #include petscao.h
11: /*S
12: DA - Abstract PETSc object that manages distributed field data for a single structured grid
14: Level: beginner
16: Concepts: distributed array
18: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), VecScatter, DACreate(), DM, DMComposite
19: S*/
20: typedef struct _p_DA* DA;
22: /*E
23: DAStencilType - Determines if the stencil extends only along the coordinate directions, or also
24: to the northeast, northwest etc
26: Level: beginner
28: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DA, DACreate()
29: E*/
30: typedef enum { DA_STENCIL_STAR,DA_STENCIL_BOX } DAStencilType;
32: /*MC
33: DA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
34: (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k)
36: Level: beginner
38: .seealso: DA_STENCIL_BOX, DAStencilType
39: M*/
41: /*MC
42: DA_STENCIL_Box - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
43: be in the stencil.
45: Level: beginner
47: .seealso: DA_STENCIL_STAR, DAStencilType
48: M*/
50: /*E
51: DAPeriodicType - Is the domain periodic in one or more directions
53: Level: beginner
55: DA_XYZGHOSTED means that ghost points are put around all the physical boundaries
56: in the local representation of the Vec (i.e. DACreate/GetLocalVector().
58: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DA, DACreate()
59: E*/
60: typedef enum { DA_NONPERIODIC,DA_XPERIODIC,DA_YPERIODIC,DA_XYPERIODIC,
61: DA_XYZPERIODIC,DA_XZPERIODIC,DA_YZPERIODIC,DA_ZPERIODIC,DA_XYZGHOSTED} DAPeriodicType;
64: /*E
65: DAInterpolationType - Defines the type of interpolation that will be returned by
66: DAGetInterpolation.
68: Level: beginner
70: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DA, DAGetInterpolation(), DASetInterpolationType(), DACreate()
71: E*/
72: typedef enum { DA_Q0, DA_Q1 } DAInterpolationType;
74: EXTERN PetscErrorCode DASetInterpolationType(DA,DAInterpolationType);
76: /*E
77: DAElementType - Defines the type of elements that will be returned by
78: DAGetElements.
80: Level: beginner
82: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DA, DAGetInterpolation(), DASetInterpolationType(),
83: DASetElementType(), DAGetElements(), DARestoreElements(), DACreate()
84: E*/
85: typedef enum { DA_ELEMENT_P1, DA_ELEMENT_Q1 } DAElementType;
87: EXTERN PetscErrorCode DASetElementType(DA,DAElementType);
88: EXTERN PetscErrorCode DAGetElements(DA,PetscInt *,const PetscInt*[]);
89: EXTERN PetscErrorCode DARestoreElements(DA,PetscInt *,const PetscInt*[]);
92: #define DAXPeriodic(pt) ((pt)==DA_XPERIODIC||(pt)==DA_XYPERIODIC||(pt)==DA_XZPERIODIC||(pt)==DA_XYZPERIODIC)
93: #define DAYPeriodic(pt) ((pt)==DA_YPERIODIC||(pt)==DA_XYPERIODIC||(pt)==DA_YZPERIODIC||(pt)==DA_XYZPERIODIC)
94: #define DAZPeriodic(pt) ((pt)==DA_ZPERIODIC||(pt)==DA_XZPERIODIC||(pt)==DA_YZPERIODIC||(pt)==DA_XYZPERIODIC)
96: typedef enum { DA_X,DA_Y,DA_Z } DADirection;
100: EXTERN PetscErrorCode DACreate1d(MPI_Comm,DAPeriodicType,PetscInt,PetscInt,PetscInt,PetscInt*,DA *);
101: EXTERN PetscErrorCode DACreate2d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,DA *);
102: EXTERN PetscErrorCode DACreate3d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscInt*,DA*);
103: EXTERN PetscErrorCode DACreate(MPI_Comm,PetscInt,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscInt*,DA*);
104: EXTERN PetscErrorCode DADestroy(DA);
105: EXTERN PetscErrorCode DAView(DA,PetscViewer);
107: EXTERN PetscErrorCode DAGlobalToLocalBegin(DA,Vec,InsertMode,Vec);
108: EXTERN PetscErrorCode DAGlobalToLocalEnd(DA,Vec,InsertMode,Vec);
109: EXTERN PetscErrorCode DAGlobalToNaturalBegin(DA,Vec,InsertMode,Vec);
110: EXTERN PetscErrorCode DAGlobalToNaturalEnd(DA,Vec,InsertMode,Vec);
111: EXTERN PetscErrorCode DANaturalToGlobalBegin(DA,Vec,InsertMode,Vec);
112: EXTERN PetscErrorCode DANaturalToGlobalEnd(DA,Vec,InsertMode,Vec);
113: EXTERN PetscErrorCode DALocalToLocalBegin(DA,Vec,InsertMode,Vec);
114: EXTERN PetscErrorCode DALocalToLocalEnd(DA,Vec,InsertMode,Vec);
115: EXTERN PetscErrorCode DALocalToGlobal(DA,Vec,InsertMode,Vec);
116: EXTERN PetscErrorCode DALocalToGlobalBegin(DA,Vec,Vec);
117: EXTERN PetscErrorCode DALocalToGlobalEnd(DA,Vec,Vec);
118: EXTERN PetscErrorCode DAGetOwnershipRange(DA,PetscInt **,PetscInt **,PetscInt **);
119: EXTERN PetscErrorCode DACreateGlobalVector(DA,Vec *);
120: EXTERN PetscErrorCode DACreateNaturalVector(DA,Vec *);
121: EXTERN PetscErrorCode DACreateLocalVector(DA,Vec *);
122: EXTERN PetscErrorCode DAGetLocalVector(DA,Vec *);
123: EXTERN PetscErrorCode DARestoreLocalVector(DA,Vec *);
124: EXTERN PetscErrorCode DAGetGlobalVector(DA,Vec *);
125: EXTERN PetscErrorCode DARestoreGlobalVector(DA,Vec *);
126: EXTERN PetscErrorCode DALoad(PetscViewer,PetscInt,PetscInt,PetscInt,DA *);
127: EXTERN PetscErrorCode DAGetCorners(DA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
128: EXTERN PetscErrorCode DAGetGhostCorners(DA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
129: EXTERN PetscErrorCode DAGetInfo(DA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DAPeriodicType*,DAStencilType*);
130: EXTERN PetscErrorCode DAGetProcessorSubset(DA,DADirection,PetscInt,MPI_Comm*);
131: EXTERN PetscErrorCode DARefine(DA,MPI_Comm,DA*);
132: EXTERN PetscErrorCode DACoarsen(DA,MPI_Comm,DA*);
134: EXTERN PetscErrorCode DAGlobalToNaturalAllCreate(DA,VecScatter*);
135: EXTERN PetscErrorCode DANaturalAllToGlobalCreate(DA,VecScatter*);
137: EXTERN PetscErrorCode DAGetGlobalIndices(DA,PetscInt*,PetscInt**);
138: EXTERN PetscErrorCode DAGetISLocalToGlobalMapping(DA,ISLocalToGlobalMapping*);
139: EXTERN PetscErrorCode DAGetISLocalToGlobalMappingBlck(DA,ISLocalToGlobalMapping*);
141: EXTERN PetscErrorCode DAGetScatter(DA,VecScatter*,VecScatter*,VecScatter*);
143: EXTERN PetscErrorCode DAGetAO(DA,AO*);
144: EXTERN PetscErrorCode DASetCoordinates(DA,Vec);
145: EXTERN PetscErrorCode DAGetCoordinates(DA,Vec *);
146: EXTERN PetscErrorCode DAGetGhostedCoordinates(DA,Vec *);
147: EXTERN PetscErrorCode DAGetCoordinateDA(DA,DA *);
148: EXTERN PetscErrorCode DASetUniformCoordinates(DA,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
149: EXTERN PetscErrorCode DASetFieldName(DA,PetscInt,const char[]);
150: EXTERN PetscErrorCode DAGetFieldName(DA,PetscInt,char **);
152: EXTERN PetscErrorCode DAVecGetArray(DA,Vec,void *);
153: EXTERN PetscErrorCode DAVecRestoreArray(DA,Vec,void *);
155: EXTERN PetscErrorCode DAVecGetArrayDOF(DA,Vec,void *);
156: EXTERN PetscErrorCode DAVecRestoreArrayDOF(DA,Vec,void *);
158: EXTERN PetscErrorCode DASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
160: /*S
161: SDA - This provides a simplified interface to the DA distributed
162: array object in PETSc. This is intended for people who are
163: NOT using PETSc vectors or objects but just want to distribute
164: simple rectangular arrays amoung a number of procesors and have
165: PETSc handle moving the ghost-values when needed.
167: In certain applications this can serve as a replacement for
168: BlockComm (which is apparently being phased out?).
171: Level: beginner
173: Concepts: simplified distributed array
175: .seealso: SDACreate1d(), SDACreate2d(), SDACreate3d(), SDADestroy(), DA, SDALocalToLocalBegin(),
176: SDALocalToLocalEnd(), SDAGetCorners(), SDAGetGhostCorners()
177: S*/
178: typedef struct _n_SDA* SDA;
180: EXTERN PetscErrorCode SDACreate3d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscInt*,SDA*);
181: EXTERN PetscErrorCode SDACreate2d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,SDA*);
182: EXTERN PetscErrorCode SDACreate1d(MPI_Comm,DAPeriodicType,PetscInt,PetscInt,PetscInt,PetscInt*,SDA*);
183: EXTERN PetscErrorCode SDADestroy(SDA);
184: EXTERN PetscErrorCode SDALocalToLocalBegin(SDA,PetscScalar*,InsertMode,PetscScalar*);
185: EXTERN PetscErrorCode SDALocalToLocalEnd(SDA,PetscScalar*,InsertMode,PetscScalar*);
186: EXTERN PetscErrorCode SDAGetCorners(SDA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
187: EXTERN PetscErrorCode SDAGetGhostCorners(SDA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
188: EXTERN PetscErrorCode SDAArrayView(SDA,PetscScalar*,PetscViewer);
190: EXTERN PetscErrorCode MatRegisterDAAD(void);
191: EXTERN PetscErrorCode MatCreateDAAD(DA,Mat*);
193: /*S
194: DALocalInfo - C struct that contains information about a structured grid and a processors logical
195: location in it.
197: Level: beginner
199: Concepts: distributed array
201: Developer note: Then entries in this struct are int instead of PetscInt so that the elements may
202: be extracted in Fortran as if from an integer array
204: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DA, DAGetLocalInfo(), DAGetInfo()
205: S*/
206: typedef struct {
207: PetscInt dim,dof,sw;
208: PetscInt mx,my,mz; /* global number of grid points in each direction */
209: PetscInt xs,ys,zs; /* starting pointd of this processor, excluding ghosts */
210: PetscInt xm,ym,zm; /* number of grid points on this processor, excluding ghosts */
211: PetscInt gxs,gys,gzs; /* starting point of this processor including ghosts */
212: PetscInt gxm,gym,gzm; /* number of grid points on this processor including ghosts */
213: DAPeriodicType pt;
214: DAStencilType st;
215: DA da;
216: } DALocalInfo;
218: /*MC
219: DAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DA
221: Synopsis:
222: void DAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j);
223:
224: Level: intermediate
226: .seealso: DAForEachPointEnd2d(), DAVecGetArray()
227: M*/
228: #define DAForEachPointBegin2d(info,i,j) {\
229: PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\
230: for (j=_yints; j<_yinte; j++) {\
231: for (i=_xints; i<_xinte; i++) {\
233: /*MC
234: DAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DA
236: Synopsis:
237: void DAForEachPointEnd2d;
238:
239: Level: intermediate
241: .seealso: DAForEachPointBegin2d(), DAVecGetArray()
242: M*/
243: #define DAForEachPointEnd2d }}}
245: /*MC
246: DACoor2d - Structure for holding 2d (x and y) coordinates.
248: Level: intermediate
250: Sample Usage:
251: DACoor2d **coors;
252: Vec vcoors;
253: DA cda;
255: DAGetCoordinates(da,&vcoors);
256: DAGetCoordinateDA(da,&cda);
257: DAVecGetArray(cda,vcoors,&coors);
258: DAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
259: for (i=mstart; i<mstart+m; i++) {
260: for (j=nstart; j<nstart+n; j++) {
261: x = coors[j][i].x;
262: y = coors[j][i].y;
263: ......
264: }
265: }
266: DAVecRestoreArray(dac,vcoors,&coors);
268: .seealso: DACoor3d, DAForEachPointBegin(), DAGetCoordinateDA(), DAGetCoordinates(), DAGetGhostCoordinates()
269: M*/
270: typedef struct {PetscScalar x,y;} DACoor2d;
272: /*MC
273: DACoor3d - Structure for holding 3d (x, y and z) coordinates.
275: Level: intermediate
277: Sample Usage:
278: DACoor3d **coors;
279: Vec vcoors;
280: DA cda;
282: DAGetCoordinates(da,&vcoors);
283: DAGetCoordinateDA(da,&cda);
284: DAVecGetArray(cda,vcoors,&coors);
285: DAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
286: for (i=mstart; i<mstart+m; i++) {
287: for (j=nstart; j<nstart+n; j++) {
288: for (k=pstart; k<pstart+p; k++) {
289: x = coors[k][j][i].x;
290: y = coors[k][j][i].y;
291: z = coors[k][j][i].z;
292: ......
293: }
294: }
295: DAVecRestoreArray(dac,vcoors,&coors);
297: .seealso: DACoor2d, DAForEachPointBegin(), DAGetCoordinateDA(), DAGetCoordinates(), DAGetGhostCoordinates()
298: M*/
299: typedef struct {PetscScalar x,y,z;} DACoor3d;
300:
301: EXTERN PetscErrorCode DAGetLocalInfo(DA,DALocalInfo*);
302: typedef PetscErrorCode (*DALocalFunction1)(DALocalInfo*,void*,void*,void*);
303: EXTERN PetscErrorCode DAFormFunctionLocal(DA, DALocalFunction1, Vec, Vec, void *);
304: EXTERN PetscErrorCode DAFormFunctionLocalGhost(DA, DALocalFunction1, Vec, Vec, void *);
305: EXTERN PetscErrorCode DAFormJacobianLocal(DA, DALocalFunction1, Vec, Mat, void *);
306: EXTERN PetscErrorCode DAFormFunction1(DA,Vec,Vec,void*);
307: EXTERN PetscErrorCode DAFormFunction(DA,PetscErrorCode (*)(void),Vec,Vec,void*);
308: EXTERN PetscErrorCode DAFormFunctioni1(DA,PetscInt,Vec,PetscScalar*,void*);
309: EXTERN PetscErrorCode DAFormFunctionib1(DA,PetscInt,Vec,PetscScalar*,void*);
310: EXTERN PetscErrorCode DAComputeJacobian1WithAdic(DA,Vec,Mat,void*);
311: EXTERN PetscErrorCode DAComputeJacobian1WithAdifor(DA,Vec,Mat,void*);
312: EXTERN PetscErrorCode DAMultiplyByJacobian1WithAdic(DA,Vec,Vec,Vec,void*);
313: EXTERN PetscErrorCode DAMultiplyByJacobian1WithAdifor(DA,Vec,Vec,Vec,void*);
314: EXTERN PetscErrorCode DAMultiplyByJacobian1WithAD(DA,Vec,Vec,Vec,void*);
315: EXTERN PetscErrorCode DAComputeJacobian1(DA,Vec,Mat,void*);
316: EXTERN PetscErrorCode DAGetLocalFunction(DA,DALocalFunction1*);
317: EXTERN PetscErrorCode DASetLocalFunction(DA,DALocalFunction1);
318: EXTERN PetscErrorCode DASetLocalFunctioni(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
319: EXTERN PetscErrorCode DASetLocalFunctionib(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
320: EXTERN PetscErrorCode DAGetLocalJacobian(DA,DALocalFunction1*);
321: EXTERN PetscErrorCode DASetLocalJacobian(DA,DALocalFunction1);
322: EXTERN PetscErrorCode DASetLocalAdicFunction_Private(DA,DALocalFunction1);
324: /*MC
325: DASetLocalAdicFunction - Caches in a DA a local function computed by ADIC/ADIFOR
327: Collective on DA
329: Synopsis:
330: PetscErrorCode DASetLocalAdicFunction(DA da,DALocalFunction1 ad_lf)
331:
332: Input Parameter:
333: + da - initial distributed array
334: - ad_lf - the local function as computed by ADIC/ADIFOR
336: Level: intermediate
338: .keywords: distributed array, refine
340: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
341: DASetLocalJacobian()
342: M*/
343: #if defined(PETSC_HAVE_ADIC)
344: # define DASetLocalAdicFunction(a,d) DASetLocalAdicFunction_Private(a,(DALocalFunction1)d)
345: #else
346: # define DASetLocalAdicFunction(a,d) DASetLocalAdicFunction_Private(a,0)
347: #endif
349: EXTERN PetscErrorCode DASetLocalAdicMFFunction_Private(DA,DALocalFunction1);
350: #if defined(PETSC_HAVE_ADIC)
351: # define DASetLocalAdicMFFunction(a,d) DASetLocalAdicMFFunction_Private(a,(DALocalFunction1)d)
352: #else
353: # define DASetLocalAdicMFFunction(a,d) DASetLocalAdicMFFunction_Private(a,0)
354: #endif
355: EXTERN PetscErrorCode DASetLocalAdicFunctioni_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
356: #if defined(PETSC_HAVE_ADIC)
357: # define DASetLocalAdicFunctioni(a,d) DASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
358: #else
359: # define DASetLocalAdicFunctioni(a,d) DASetLocalAdicFunctioni_Private(a,0)
360: #endif
361: EXTERN PetscErrorCode DASetLocalAdicMFFunctioni_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
362: #if defined(PETSC_HAVE_ADIC)
363: # define DASetLocalAdicMFFunctioni(a,d) DASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
364: #else
365: # define DASetLocalAdicMFFunctioni(a,d) DASetLocalAdicMFFunctioni_Private(a,0)
366: #endif
368: EXTERN PetscErrorCode DASetLocalAdicFunctionib_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
369: #if defined(PETSC_HAVE_ADIC)
370: # define DASetLocalAdicFunctionib(a,d) DASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
371: #else
372: # define DASetLocalAdicFunctionib(a,d) DASetLocalAdicFunctionib_Private(a,0)
373: #endif
374: EXTERN PetscErrorCode DASetLocalAdicMFFunctionib_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
375: #if defined(PETSC_HAVE_ADIC)
376: # define DASetLocalAdicMFFunctionib(a,d) DASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
377: #else
378: # define DASetLocalAdicMFFunctionib(a,d) DASetLocalAdicMFFunctionib_Private(a,0)
379: #endif
381: EXTERN PetscErrorCode DAFormFunctioniTest1(DA,void*);
383: #include petscmat.h
384: EXTERN PetscErrorCode DAGetColoring(DA,ISColoringType,ISColoring *);
385: EXTERN PetscErrorCode DAGetMatrix(DA, MatType,Mat *);
386: EXTERN PetscErrorCode DASetGetMatrix(DA,PetscErrorCode (*)(DA, MatType,Mat *));
387: EXTERN PetscErrorCode DAGetInterpolation(DA,DA,Mat*,Vec*);
388: EXTERN PetscErrorCode DAGetAggregates(DA,DA,Mat*);
389: EXTERN PetscErrorCode DAGetInjection(DA,DA,VecScatter*);
390: EXTERN PetscErrorCode DASetBlockFills(DA,PetscInt*,PetscInt*);
391: EXTERN PetscErrorCode DASetMatPreallocateOnly(DA,PetscTruth);
392: EXTERN PetscErrorCode DASetRefinementFactor(DA,PetscInt,PetscInt,PetscInt);
393: EXTERN PetscErrorCode DAGetRefinementFactor(DA,PetscInt*,PetscInt*,PetscInt*);
395: EXTERN PetscErrorCode DAGetAdicArray(DA,PetscTruth,void**,void**,PetscInt*);
396: EXTERN PetscErrorCode DARestoreAdicArray(DA,PetscTruth,void**,void**,PetscInt*);
397: EXTERN PetscErrorCode DAGetAdicMFArray(DA,PetscTruth,void**,void**,PetscInt*);
398: EXTERN PetscErrorCode DAGetAdicMFArray4(DA,PetscTruth,void**,void**,PetscInt*);
399: EXTERN PetscErrorCode DAGetAdicMFArray9(DA,PetscTruth,void**,void**,PetscInt*);
400: EXTERN PetscErrorCode DAGetAdicMFArrayb(DA,PetscTruth,void**,void**,PetscInt*);
401: EXTERN PetscErrorCode DARestoreAdicMFArray(DA,PetscTruth,void**,void**,PetscInt*);
402: EXTERN PetscErrorCode DAGetArray(DA,PetscTruth,void**);
403: EXTERN PetscErrorCode DARestoreArray(DA,PetscTruth,void**);
404: EXTERN PetscErrorCode ad_DAGetArray(DA,PetscTruth,void**);
405: EXTERN PetscErrorCode ad_DARestoreArray(DA,PetscTruth,void**);
406: EXTERN PetscErrorCode admf_DAGetArray(DA,PetscTruth,void**);
407: EXTERN PetscErrorCode admf_DARestoreArray(DA,PetscTruth,void**);
409: #include petscpf.h
410: EXTERN PetscErrorCode DACreatePF(DA,PF*);
412: /*S
413: DMComposite - Abstract PETSc object that manages treating several distinct vectors as if they
414: were one. The DMComposite routines allow one to manage a nonlinear solver that works on a
415: vector that consists of several distinct parts. This is mostly used for LNKS solvers,
416: that is design optimization problems that are written as a nonlinear system
418: Level: beginner
420: Concepts: multi-component, LNKS solvers
422: .seealso: DMCompositeCreate(), DMCompositeDestroy(), DM
423: S*/
424: typedef struct _p_DMComposite* DMComposite;
426: EXTERN PetscErrorCode DMCompositeCreate(MPI_Comm,DMComposite*);
427: EXTERN PetscErrorCode DMCompositeDestroy(DMComposite);
428: EXTERN PetscErrorCode DMCompositeAddArray(DMComposite,PetscMPIInt,PetscInt);
429: EXTERN PetscErrorCode DMCompositeAddDA(DMComposite,DA);
430: EXTERN PetscErrorCode DMCompositeAddVecScatter(DMComposite,VecScatter);
431: EXTERN PetscErrorCode DMCompositeScatter(DMComposite,Vec,...);
432: EXTERN PetscErrorCode DMCompositeGather(DMComposite,Vec,...);
433: EXTERN PetscErrorCode DMCompositeGetAccess(DMComposite,Vec,...);
434: EXTERN PetscErrorCode DMCompositeRestoreAccess(DMComposite,Vec,...);
435: EXTERN PetscErrorCode DMCompositeGetLocalVectors(DMComposite,...);
436: EXTERN PetscErrorCode DMCompositeGetEntries(DMComposite,...);
437: EXTERN PetscErrorCode DMCompositeRestoreLocalVectors(DMComposite,...);
438: EXTERN PetscErrorCode DMCompositeCreateGlobalVector(DMComposite,Vec*);
439: EXTERN PetscErrorCode DMCompositeGetGlobalIndices(DMComposite,...);
440: EXTERN PetscErrorCode DMCompositeRefine(DMComposite,MPI_Comm,DMComposite*);
441: EXTERN PetscErrorCode DMCompositeGetInterpolation(DMComposite,DMComposite,Mat*,Vec*);
442: EXTERN PetscErrorCode DMCompositeGetMatrix(DMComposite,MatType,Mat*);
443: EXTERN PetscErrorCode DMCompositeGetColoring(DMComposite,ISColoringType,ISColoring*);
445: /*S
446: Slice - Abstract PETSc object that manages distributed field data for a simple unstructured matrix
448: Level: beginner
450: Concepts: distributed array
452: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), VecScatter, DACreate(), DMCompositeCreate(), DMComposite
453: S*/
454: typedef struct _p_Sliced* Sliced;
456: EXTERN PetscErrorCode SlicedView(Sliced,PetscViewer);
457: EXTERN PetscErrorCode SlicedCreate(MPI_Comm,Sliced*);
458: EXTERN PetscErrorCode SlicedDestroy(Sliced);
459: EXTERN PetscErrorCode SlicedCreateGlobalVector(Sliced,Vec*);
460: EXTERN PetscErrorCode SlicedGetMatrix(Sliced, MatType,Mat*);
461: EXTERN PetscErrorCode SlicedGetGlobalIndices(Sliced,PetscInt*[]);
462: EXTERN PetscErrorCode SlicedSetPreallocation(Sliced,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
463: EXTERN PetscErrorCode SlicedSetGhosts(Sliced,PetscInt,PetscInt,PetscInt,const PetscInt[]);
465: /*S
466: DM - Abstract PETSc object that manages an abstract grid object
467:
468: Level: intermediate
470: Concepts: grids, grid refinement
472: Notes: The DA object and the DMComposite object are examples of DMs
474: Though the DA objects require the petscsnes.h include files the DM library is
475: NOT dependent on the SNES or KSP library. In fact, the KSP and SNES libraries depend on
476: DM. (This is not great design, but not trivial to fix).
478: .seealso: DMCompositeCreate(), DA, DMComposite
479: S*/
480: typedef struct _p_DM* DM;
482: EXTERN PetscErrorCode DMView(DM,PetscViewer);
483: EXTERN PetscErrorCode DMDestroy(DM);
484: EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
485: EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
486: EXTERN PetscErrorCode DMGetColoring(DM,ISColoringType,ISColoring*);
487: EXTERN PetscErrorCode DMGetMatrix(DM, MatType,Mat*);
488: EXTERN PetscErrorCode DMGetInterpolation(DM,DM,Mat*,Vec*);
489: EXTERN PetscErrorCode DMGetInjection(DM,DM,VecScatter*);
490: EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
491: EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
492: EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM**);
493: EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM**);
494: EXTERN PetscErrorCode DMGetInterpolationScale(DM,DM,Mat,Vec*);
495: EXTERN PetscErrorCode DMGetAggregates(DM,DM,Mat*);
496: EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
497: EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
498: EXTERN PetscErrorCode DMLocalToGlobal(DM,Vec,InsertMode,Vec);
500: EXTERN PetscErrorCode DMFinalizePackage(void);
502: typedef struct NLF_DAAD* NLF;
504: #include <petscbag.h>
506: EXTERN PetscErrorCode PetscViewerBinaryMatlabOpen(MPI_Comm, const char [], PetscViewer*);
507: EXTERN PetscErrorCode PetscViewerBinaryMatlabDestroy(PetscViewer);
508: EXTERN PetscErrorCode PetscViewerBinaryMatlabOutputBag(PetscViewer, const char [], PetscBag);
509: EXTERN PetscErrorCode PetscViewerBinaryMatlabOutputVec(PetscViewer, const char [], Vec);
510: EXTERN PetscErrorCode PetscViewerBinaryMatlabOutputVecDA(PetscViewer, const char [], Vec, DA);
513: /*S
514: ADDA - Abstract PETSc object that manages distributed field data for a single structured grid
515: These are for any number of dimensions.
517: Level: advanced.
519: Concepts: distributed array
520: .seealso: DA, DACreate(), ADDACreate()
521: S*/
522: typedef struct _p_ADDA* ADDA;
526: PetscErrorCode ADDACreate(MPI_Comm,PetscInt,PetscInt*,PetscInt*,PetscInt,PetscTruth*,ADDA*);
527: PetscErrorCode ADDADestroy(ADDA);
529: /* DM interface functions */
530: PetscErrorCode ADDAView(ADDA,PetscViewer);
531: PetscErrorCode ADDACreateGlobalVector(ADDA,Vec*);
532: PetscErrorCode ADDAGetColoring(ADDA,ISColoringType,ISColoring*);
533: PetscErrorCode ADDAGetMatrix(ADDA,MatType, Mat*);
534: PetscErrorCode ADDAGetInterpolation(ADDA,ADDA,Mat*,Vec*);
535: PetscErrorCode ADDARefine(ADDA, MPI_Comm,ADDA *);
536: PetscErrorCode ADDACoarsen(ADDA, MPI_Comm, ADDA*);
537: PetscErrorCode ADDAGetInjection(ADDA, ADDA, VecScatter*);
538: PetscErrorCode ADDAGetAggregates(ADDA, ADDA, Mat *);
540: /* functions only supported by ADDA */
541: PetscErrorCode ADDASetRefinement(ADDA, PetscInt *,PetscInt);
542: PetscErrorCode ADDAGetCorners(ADDA, PetscInt **, PetscInt **);
543: PetscErrorCode ADDAGetGhostCorners(ADDA, PetscInt **, PetscInt **);
544: PetscErrorCode ADDAGetMatrixNS(ADDA, ADDA, MatType , Mat *);
546: /* functions to set values in vectors and matrices */
547: struct _ADDAIdx_s {
548: PetscInt *x; /* the coordinates, user has to make sure it is the correct size! */
549: PetscInt d; /* indexes the dof */
550: };
551: typedef struct _ADDAIdx_s ADDAIdx;
553: PetscErrorCode ADDAMatSetValues(Mat, ADDA, PetscInt, const ADDAIdx[], ADDA, PetscInt,
554: const ADDAIdx[], const PetscScalar[], InsertMode);
556: PetscTruth ADDAHCiterStartup(const PetscInt, const PetscInt *const, const PetscInt *const, PetscInt *const);
557: PetscTruth ADDAHCiter(const PetscInt, const PetscInt *const, const PetscInt *const, PetscInt *const);
560: #endif