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: #define DAGetElements(da,a,b)      DMGetElements((DM)da,a,b)
 89: #define DARestoreElements(da,a,b)  DMRestoreElements((DM)da,a,b)


 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: #define MATSEQUSFFT        "sequsfft"

102: EXTERN PetscErrorCode     DACreate1d(MPI_Comm,DAPeriodicType,PetscInt,PetscInt,PetscInt,const PetscInt[],DA *);
103: EXTERN PetscErrorCode     DACreate2d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DA *);
104: EXTERN PetscErrorCode     DACreate3d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DA*);
105: EXTERN PetscErrorCode     DACreate(MPI_Comm,PetscInt,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DA*);
106: EXTERN PetscErrorCode     DADestroy(DA);
107: EXTERN PetscErrorCode     DAView(DA,PetscViewer);

109: EXTERN PetscErrorCode     DAGlobalToLocalBegin(DA,Vec,InsertMode,Vec);
110: EXTERN PetscErrorCode     DAGlobalToLocalEnd(DA,Vec,InsertMode,Vec);
111: EXTERN PetscErrorCode     DAGlobalToNaturalBegin(DA,Vec,InsertMode,Vec);
112: EXTERN PetscErrorCode     DAGlobalToNaturalEnd(DA,Vec,InsertMode,Vec);
113: EXTERN PetscErrorCode     DANaturalToGlobalBegin(DA,Vec,InsertMode,Vec);
114: EXTERN PetscErrorCode     DANaturalToGlobalEnd(DA,Vec,InsertMode,Vec);
115: EXTERN PetscErrorCode     DALocalToLocalBegin(DA,Vec,InsertMode,Vec);
116: EXTERN PetscErrorCode     DALocalToLocalEnd(DA,Vec,InsertMode,Vec);
117: EXTERN PetscErrorCode     DALocalToGlobal(DA,Vec,InsertMode,Vec);
118: EXTERN PetscErrorCode     DALocalToGlobalBegin(DA,Vec,Vec);
119: EXTERN PetscErrorCode     DALocalToGlobalEnd(DA,Vec,Vec);
120: EXTERN PetscErrorCode     DAGetOwnershipRanges(DA,const PetscInt*[],const PetscInt*[],const PetscInt*[]);
121: EXTERN PetscErrorCode     DACreateGlobalVector(DA,Vec *);
122: EXTERN PetscErrorCode     DACreateLocalVector(DA,Vec *);
123: EXTERN PetscErrorCode     DACreateNaturalVector(DA,Vec *);
124: #define  DAGetLocalVector(da,v)      DMGetLocalVector((DM)da,v)
125: #define  DARestoreLocalVector(da,v)  DMRestoreLocalVector((DM)da,v)
126: #define  DAGetGlobalVector(da,v)     DMGetGlobalVector((DM)da,v)
127: #define  DARestoreGlobalVector(da,v) DMRestoreGlobalVector((DM)da,v)
128: EXTERN PetscErrorCode     DALoad(PetscViewer,PetscInt,PetscInt,PetscInt,DA *);
129: EXTERN PetscErrorCode     DAGetCorners(DA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
130: EXTERN PetscErrorCode     DAGetGhostCorners(DA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
131: EXTERN PetscErrorCode     DAGetInfo(DA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DAPeriodicType*,DAStencilType*);
132: EXTERN PetscErrorCode     DAGetProcessorSubset(DA,DADirection,PetscInt,MPI_Comm*);
133: EXTERN PetscErrorCode     DARefine(DA,MPI_Comm,DA*);
134: EXTERN PetscErrorCode     DACoarsen(DA,MPI_Comm,DA*);

136: EXTERN PetscErrorCode     DAGlobalToNaturalAllCreate(DA,VecScatter*);
137: EXTERN PetscErrorCode     DANaturalAllToGlobalCreate(DA,VecScatter*);

139: EXTERN PetscErrorCode     DAGetGlobalIndices(DA,PetscInt*,PetscInt**);
140: EXTERN PetscErrorCode     DAGetISLocalToGlobalMapping(DA,ISLocalToGlobalMapping*);
141: EXTERN PetscErrorCode     DAGetISLocalToGlobalMappingBlck(DA,ISLocalToGlobalMapping*);

143: EXTERN PetscErrorCode     DAGetScatter(DA,VecScatter*,VecScatter*,VecScatter*);
144: EXTERN PetscErrorCode     DAGetNeighbors(DA,const PetscMPIInt**);

146: EXTERN PetscErrorCode     DAGetAO(DA,AO*);
147: EXTERN PetscErrorCode     DASetCoordinates(DA,Vec);
148: EXTERN PetscErrorCode     DAGetCoordinates(DA,Vec *);
149: EXTERN PetscErrorCode     DAGetGhostedCoordinates(DA,Vec *);
150: EXTERN PetscErrorCode     DAGetCoordinateDA(DA,DA *);
151: EXTERN PetscErrorCode     DASetUniformCoordinates(DA,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
152: EXTERN PetscErrorCode     DASetFieldName(DA,PetscInt,const char[]);
153: EXTERN PetscErrorCode     DAGetFieldName(DA,PetscInt,char **);

155: EXTERN PetscErrorCode     DAVecGetArray(DA,Vec,void *);
156: EXTERN PetscErrorCode     DAVecRestoreArray(DA,Vec,void *);

158: EXTERN PetscErrorCode     DAVecGetArrayDOF(DA,Vec,void *);
159: EXTERN PetscErrorCode     DAVecRestoreArrayDOF(DA,Vec,void *);

161: EXTERN PetscErrorCode     DASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);

163: /*S
164:      SDA - This provides a simplified interface to the DA distributed
165:            array object in PETSc. This is intended for people who are
166:            NOT using PETSc vectors or objects but just want to distribute
167:            simple rectangular arrays amoung a number of procesors and have
168:            PETSc handle moving the ghost-values when needed.

170:           In certain applications this can serve as a replacement for 
171:           BlockComm (which is apparently being phased out?).


174:    Level: beginner

176:   Concepts: simplified distributed array

178: .seealso:  SDACreate1d(), SDACreate2d(), SDACreate3d(), SDADestroy(), DA, SDALocalToLocalBegin(),
179:            SDALocalToLocalEnd(), SDAGetCorners(), SDAGetGhostCorners()
180: S*/
181: typedef struct _n_SDA* SDA;

183: EXTERN PetscErrorCode     SDACreate3d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],SDA*);
184: EXTERN PetscErrorCode     SDACreate2d(MPI_Comm,DAPeriodicType,DAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],SDA*);
185: EXTERN PetscErrorCode     SDACreate1d(MPI_Comm,DAPeriodicType,PetscInt,PetscInt,PetscInt,const PetscInt[],SDA*);
186: EXTERN PetscErrorCode     SDADestroy(SDA);
187: EXTERN PetscErrorCode     SDALocalToLocalBegin(SDA,PetscScalar*,InsertMode,PetscScalar*);
188: EXTERN PetscErrorCode     SDALocalToLocalEnd(SDA,PetscScalar*,InsertMode,PetscScalar*);
189: EXTERN PetscErrorCode     SDAGetCorners(SDA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
190: EXTERN PetscErrorCode     SDAGetGhostCorners(SDA,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
191: EXTERN PetscErrorCode     SDAArrayView(SDA,PetscScalar*,PetscViewer);

193: EXTERN PetscErrorCode     MatRegisterDAAD(void);
194: EXTERN PetscErrorCode     MatCreateDAAD(DA,Mat*);
195: EXTERN PetscErrorCode    MatCreateSeqUSFFT(DA, DA,Mat*);

197: /*S
198:      DALocalInfo - C struct that contains information about a structured grid and a processors logical
199:               location in it.

201:    Level: beginner

203:   Concepts: distributed array

205:   Developer note: Then entries in this struct are int instead of PetscInt so that the elements may
206:                   be extracted in Fortran as if from an integer array

208: .seealso:  DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DA, DAGetLocalInfo(), DAGetInfo()
209: S*/
210: typedef struct {
211:   PetscInt       dim,dof,sw;
212:   PetscInt       mx,my,mz;    /* global number of grid points in each direction */
213:   PetscInt       xs,ys,zs;    /* starting pointd of this processor, excluding ghosts */
214:   PetscInt       xm,ym,zm;    /* number of grid points on this processor, excluding ghosts */
215:   PetscInt       gxs,gys,gzs;    /* starting point of this processor including ghosts */
216:   PetscInt       gxm,gym,gzm;    /* number of grid points on this processor including ghosts */
217:   DAPeriodicType pt;
218:   DAStencilType  st;
219:   DA             da;
220: } DALocalInfo;

222: /*MC
223:       DAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DA

225:    Synopsis:
226:    void  DAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j);
227:    
228:    Level: intermediate

230: .seealso: DAForEachPointEnd2d(), DAVecGetArray()
231: M*/
232: #define DAForEachPointBegin2d(info,i,j) {\
233:   PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\
234:   for (j=_yints; j<_yinte; j++) {\
235:     for (i=_xints; i<_xinte; i++) {\

237: /*MC
238:       DAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DA

240:    Synopsis:
241:    void  DAForEachPointEnd2d;
242:    
243:    Level: intermediate

245: .seealso: DAForEachPointBegin2d(), DAVecGetArray()
246: M*/
247: #define DAForEachPointEnd2d }}}

249: /*MC
250:       DACoor2d - Structure for holding 2d (x and y) coordinates.

252:     Level: intermediate

254:     Sample Usage:
255:       DACoor2d **coors;
256:       Vec      vcoors;
257:       DA       cda;     

259:       DAGetCoordinates(da,&vcoors); 
260:       DAGetCoordinateDA(da,&cda);
261:       DAVecGetArray(cda,vcoors,&coors);
262:       DAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
263:       for (i=mstart; i<mstart+m; i++) {
264:         for (j=nstart; j<nstart+n; j++) {
265:           x = coors[j][i].x;
266:           y = coors[j][i].y;
267:           ......
268:         }
269:       }
270:       DAVecRestoreArray(dac,vcoors,&coors);

272: .seealso: DACoor3d, DAForEachPointBegin(), DAGetCoordinateDA(), DAGetCoordinates(), DAGetGhostCoordinates()
273: M*/
274: typedef struct {PetscScalar x,y;} DACoor2d;

276: /*MC
277:       DACoor3d - Structure for holding 3d (x, y and z) coordinates.

279:     Level: intermediate

281:     Sample Usage:
282:       DACoor3d ***coors;
283:       Vec      vcoors;
284:       DA       cda;     

286:       DAGetCoordinates(da,&vcoors); 
287:       DAGetCoordinateDA(da,&cda);
288:       DAVecGetArray(cda,vcoors,&coors);
289:       DAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
290:       for (i=mstart; i<mstart+m; i++) {
291:         for (j=nstart; j<nstart+n; j++) {
292:           for (k=pstart; k<pstart+p; k++) {
293:             x = coors[k][j][i].x;
294:             y = coors[k][j][i].y;
295:             z = coors[k][j][i].z;
296:           ......
297:         }
298:       }
299:       DAVecRestoreArray(dac,vcoors,&coors);

301: .seealso: DACoor2d, DAForEachPointBegin(), DAGetCoordinateDA(), DAGetCoordinates(), DAGetGhostCoordinates()
302: M*/
303: typedef struct {PetscScalar x,y,z;} DACoor3d;
304: 
305: EXTERN PetscErrorCode   DAGetLocalInfo(DA,DALocalInfo*);
306: typedef PetscErrorCode (*DALocalFunction1)(DALocalInfo*,void*,void*,void*);
307: EXTERN PetscErrorCode   DAFormFunctionLocal(DA, DALocalFunction1, Vec, Vec, void *);
308: EXTERN PetscErrorCode   DAFormFunctionLocalGhost(DA, DALocalFunction1, Vec, Vec, void *);
309: EXTERN PetscErrorCode   DAFormJacobianLocal(DA, DALocalFunction1, Vec, Mat, void *);
310: EXTERN PetscErrorCode   DAFormFunction1(DA,Vec,Vec,void*);
311: EXTERN PetscErrorCode   DAFormFunction(DA,PetscErrorCode (*)(void),Vec,Vec,void*);
312: EXTERN PetscErrorCode   DAFormFunctioni1(DA,PetscInt,Vec,PetscScalar*,void*);
313: EXTERN PetscErrorCode   DAFormFunctionib1(DA,PetscInt,Vec,PetscScalar*,void*);
314: EXTERN PetscErrorCode   DAComputeJacobian1WithAdic(DA,Vec,Mat,void*);
315: EXTERN PetscErrorCode   DAComputeJacobian1WithAdifor(DA,Vec,Mat,void*);
316: EXTERN PetscErrorCode   DAMultiplyByJacobian1WithAdic(DA,Vec,Vec,Vec,void*);
317: EXTERN PetscErrorCode   DAMultiplyByJacobian1WithAdifor(DA,Vec,Vec,Vec,void*);
318: EXTERN PetscErrorCode   DAMultiplyByJacobian1WithAD(DA,Vec,Vec,Vec,void*);
319: EXTERN PetscErrorCode   DAComputeJacobian1(DA,Vec,Mat,void*);
320: EXTERN PetscErrorCode   DAGetLocalFunction(DA,DALocalFunction1*);
321: EXTERN PetscErrorCode   DASetLocalFunction(DA,DALocalFunction1);
322: EXTERN PetscErrorCode   DASetLocalFunctioni(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
323: EXTERN PetscErrorCode   DASetLocalFunctionib(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
324: EXTERN PetscErrorCode   DAGetLocalJacobian(DA,DALocalFunction1*);
325: EXTERN PetscErrorCode   DASetLocalJacobian(DA,DALocalFunction1);
326: EXTERN PetscErrorCode   DASetLocalAdicFunction_Private(DA,DALocalFunction1);

328: /*MC
329:        DASetLocalAdicFunction - Caches in a DA a local function computed by ADIC/ADIFOR

331:    Collective on DA

333:    Synopsis:
334:    PetscErrorCode DASetLocalAdicFunction(DA da,DALocalFunction1 ad_lf)
335:    
336:    Input Parameter:
337: +  da - initial distributed array
338: -  ad_lf - the local function as computed by ADIC/ADIFOR

340:    Level: intermediate

342: .keywords:  distributed array, refine

344: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
345:           DASetLocalJacobian()
346: M*/
347: #if defined(PETSC_HAVE_ADIC)
348: #  define DASetLocalAdicFunction(a,d) DASetLocalAdicFunction_Private(a,(DALocalFunction1)d)
349: #else
350: #  define DASetLocalAdicFunction(a,d) DASetLocalAdicFunction_Private(a,0)
351: #endif

353: EXTERN PetscErrorCode   DASetLocalAdicMFFunction_Private(DA,DALocalFunction1);
354: #if defined(PETSC_HAVE_ADIC)
355: #  define DASetLocalAdicMFFunction(a,d) DASetLocalAdicMFFunction_Private(a,(DALocalFunction1)d)
356: #else
357: #  define DASetLocalAdicMFFunction(a,d) DASetLocalAdicMFFunction_Private(a,0)
358: #endif
359: EXTERN PetscErrorCode   DASetLocalAdicFunctioni_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
360: #if defined(PETSC_HAVE_ADIC)
361: #  define DASetLocalAdicFunctioni(a,d) DASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
362: #else
363: #  define DASetLocalAdicFunctioni(a,d) DASetLocalAdicFunctioni_Private(a,0)
364: #endif
365: EXTERN PetscErrorCode   DASetLocalAdicMFFunctioni_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
366: #if defined(PETSC_HAVE_ADIC)
367: #  define DASetLocalAdicMFFunctioni(a,d) DASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
368: #else
369: #  define DASetLocalAdicMFFunctioni(a,d) DASetLocalAdicMFFunctioni_Private(a,0)
370: #endif

372: EXTERN PetscErrorCode   DASetLocalAdicFunctionib_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
373: #if defined(PETSC_HAVE_ADIC)
374: #  define DASetLocalAdicFunctionib(a,d) DASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
375: #else
376: #  define DASetLocalAdicFunctionib(a,d) DASetLocalAdicFunctionib_Private(a,0)
377: #endif
378: EXTERN PetscErrorCode   DASetLocalAdicMFFunctionib_Private(DA,PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*));
379: #if defined(PETSC_HAVE_ADIC)
380: #  define DASetLocalAdicMFFunctionib(a,d) DASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DALocalInfo*,MatStencil*,void*,void*,void*))d)
381: #else
382: #  define DASetLocalAdicMFFunctionib(a,d) DASetLocalAdicMFFunctionib_Private(a,0)
383: #endif

385: EXTERN PetscErrorCode   DAFormFunctioniTest1(DA,void*);

387:  #include petscmat.h

389: /*S
390:      DM - Abstract PETSc object that manages an abstract grid object
391:           
392:    Level: intermediate

394:   Concepts: grids, grid refinement

396:    Notes: The DA object and the DMComposite object are examples of DMs

398:           Though the DA objects require the petscsnes.h include files the DM library is
399:     NOT dependent on the SNES or KSP library. In fact, the KSP and SNES libraries depend on
400:     DM. (This is not great design, but not trivial to fix).

402: .seealso:  DMCompositeCreate(), DA, DMComposite
403: S*/
404: typedef struct _p_DM* DM;

406: EXTERN PetscErrorCode   DMView(DM,PetscViewer);
407: EXTERN PetscErrorCode   DMDestroy(DM);
408: EXTERN PetscErrorCode   DMCreateGlobalVector(DM,Vec*);
409: EXTERN PetscErrorCode   DMCreateLocalVector(DM,Vec*);
410: EXTERN PetscErrorCode   DMGetLocalVector(DM,Vec *);
411: EXTERN PetscErrorCode   DMRestoreLocalVector(DM,Vec *);
412: EXTERN PetscErrorCode   DMGetGlobalVector(DM,Vec *);
413: EXTERN PetscErrorCode   DMRestoreGlobalVector(DM,Vec *);
414: EXTERN PetscErrorCode   DMGetColoring(DM,ISColoringType,ISColoring*);
415: EXTERN PetscErrorCode   DMGetMatrix(DM, const MatType,Mat*);
416: EXTERN PetscErrorCode   DMGetInterpolation(DM,DM,Mat*,Vec*);
417: EXTERN PetscErrorCode   DMGetInjection(DM,DM,VecScatter*);
418: EXTERN PetscErrorCode   DMRefine(DM,MPI_Comm,DM*);
419: EXTERN PetscErrorCode   DMCoarsen(DM,MPI_Comm,DM*);
420: EXTERN PetscErrorCode   DMRefineHierarchy(DM,PetscInt,DM**);
421: EXTERN PetscErrorCode   DMCoarsenHierarchy(DM,PetscInt,DM**);
422: EXTERN PetscErrorCode   DMGetInterpolationScale(DM,DM,Mat,Vec*);
423: EXTERN PetscErrorCode   DMGetAggregates(DM,DM,Mat*);
424: EXTERN PetscErrorCode   DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
425: EXTERN PetscErrorCode   DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
426: EXTERN PetscErrorCode   DMLocalToGlobal(DM,Vec,InsertMode,Vec);
427: EXTERN PetscErrorCode   DMGetElements(DM,PetscInt *,const PetscInt*[]);
428: EXTERN PetscErrorCode   DMRestoreElements(DM,PetscInt *,const PetscInt*[]);
429: EXTERN PetscErrorCode   DMFinalizePackage(void);

431: EXTERN PetscErrorCode   DAGetColoring(DA,ISColoringType,ISColoring *);
432: EXTERN PetscErrorCode   DAGetMatrix(DA, const MatType,Mat *);
433: EXTERN PetscErrorCode   DASetGetMatrix(DA,PetscErrorCode (*)(DA, const MatType,Mat *));
434: EXTERN PetscErrorCode   DAGetInterpolation(DA,DA,Mat*,Vec*);
435: EXTERN PetscErrorCode   DAGetAggregates(DA,DA,Mat*);
436: EXTERN PetscErrorCode   DAGetInjection(DA,DA,VecScatter*);
437: EXTERN PetscErrorCode   DASetBlockFills(DA,PetscInt*,PetscInt*);
438: EXTERN PetscErrorCode   DASetMatPreallocateOnly(DA,PetscTruth);
439: EXTERN PetscErrorCode   DASetRefinementFactor(DA,PetscInt,PetscInt,PetscInt);
440: EXTERN PetscErrorCode   DAGetRefinementFactor(DA,PetscInt*,PetscInt*,PetscInt*);

442: EXTERN PetscErrorCode   DAGetAdicArray(DA,PetscTruth,void**,void**,PetscInt*);
443: EXTERN PetscErrorCode   DARestoreAdicArray(DA,PetscTruth,void**,void**,PetscInt*);
444: EXTERN PetscErrorCode   DAGetAdicMFArray(DA,PetscTruth,void**,void**,PetscInt*);
445: EXTERN PetscErrorCode   DAGetAdicMFArray4(DA,PetscTruth,void**,void**,PetscInt*);
446: EXTERN PetscErrorCode   DAGetAdicMFArray9(DA,PetscTruth,void**,void**,PetscInt*);
447: EXTERN PetscErrorCode   DAGetAdicMFArrayb(DA,PetscTruth,void**,void**,PetscInt*);
448: EXTERN PetscErrorCode   DARestoreAdicMFArray(DA,PetscTruth,void**,void**,PetscInt*);
449: EXTERN PetscErrorCode   DAGetArray(DA,PetscTruth,void**);
450: EXTERN PetscErrorCode   DARestoreArray(DA,PetscTruth,void**);
451: EXTERN PetscErrorCode   ad_DAGetArray(DA,PetscTruth,void**);
452: EXTERN PetscErrorCode   ad_DARestoreArray(DA,PetscTruth,void**);
453: EXTERN PetscErrorCode   admf_DAGetArray(DA,PetscTruth,void**);
454: EXTERN PetscErrorCode   admf_DARestoreArray(DA,PetscTruth,void**);

456:  #include petscpf.h
457: EXTERN PetscErrorCode   DACreatePF(DA,PF*);

459: /*S
460:      DMComposite - Abstract PETSc object that manages treating several distinct vectors as if they
461:         were one.   The DMComposite routines allow one to manage a nonlinear solver that works on a
462:         vector that consists of several distinct parts. This is mostly used for LNKS solvers, 
463:         that is design optimization problems that are written as a nonlinear system

465:    Level: beginner

467:   Concepts: multi-component, LNKS solvers

469: .seealso:  DMCompositeCreate(), DMCompositeDestroy(), DM
470: S*/
471: typedef struct _p_DMComposite* DMComposite;

473: EXTERN PetscErrorCode   DMCompositeCreate(MPI_Comm,DMComposite*);
474: EXTERN PetscErrorCode   DMCompositeDestroy(DMComposite);
475: EXTERN PetscErrorCode   DMCompositeAddArray(DMComposite,PetscMPIInt,PetscInt);
476: EXTERN PetscErrorCode   DMCompositeAddDM(DMComposite,DM);
477: EXTERN PetscErrorCode   DMCompositeSetCoupling(DMComposite,PetscErrorCode (*)(DMComposite,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt));
478: EXTERN PetscErrorCode   DMCompositeSetContext(DMComposite,void*);
479: EXTERN PetscErrorCode   DMCompositeGetContext(DMComposite,void**);
480: EXTERN PetscErrorCode   DMCompositeAddVecScatter(DMComposite,VecScatter);
481: EXTERN PetscErrorCode   DMCompositeScatter(DMComposite,Vec,...);
482: EXTERN PetscErrorCode   DMCompositeGather(DMComposite,Vec,...);
483: EXTERN PetscErrorCode   DMCompositeGetAccess(DMComposite,Vec,...);
484: EXTERN PetscErrorCode   DMCompositeGetNumberDM(DMComposite,PetscInt*);
485: EXTERN PetscErrorCode   DMCompositeRestoreAccess(DMComposite,Vec,...);
486: EXTERN PetscErrorCode   DMCompositeGetLocalVectors(DMComposite,...);
487: EXTERN PetscErrorCode   DMCompositeGetEntries(DMComposite,...);
488: EXTERN PetscErrorCode   DMCompositeRestoreLocalVectors(DMComposite,...);
489: EXTERN PetscErrorCode   DMCompositeCreateGlobalVector(DMComposite,Vec*);
490: EXTERN PetscErrorCode   DMCompositeCreateLocalVector(DMComposite,Vec*);
491: EXTERN PetscErrorCode   DMCompositeGetLocalISs(DMComposite,IS*[]);
492: EXTERN PetscErrorCode   DMCompositeGetGlobalISs(DMComposite,IS*[]);
493: EXTERN PetscErrorCode   DMCompositeRefine(DMComposite,MPI_Comm,DMComposite*);
494: EXTERN PetscErrorCode   DMCompositeGetInterpolation(DMComposite,DMComposite,Mat*,Vec*);
495: EXTERN PetscErrorCode   DMCompositeGetMatrix(DMComposite,const MatType,Mat*);
496: EXTERN PetscErrorCode   DMCompositeGetColoring(DMComposite,ISColoringType,ISColoring*);
497: EXTERN PetscErrorCode   DMCompositeGlobalToLocalBegin(DMComposite,Vec,InsertMode,Vec);
498: EXTERN PetscErrorCode   DMCompositeGlobalToLocalEnd(DMComposite,Vec,InsertMode,Vec);

500: /*S
501:      Slice - Abstract PETSc object that manages distributed field data for a simple unstructured matrix

503:    Level: beginner

505:   Concepts: distributed array

507: .seealso:  DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), VecScatter, DACreate(), DMCompositeCreate(), DMComposite
508: S*/
509: typedef struct _p_Sliced* Sliced;

511: EXTERN PetscErrorCode   SlicedView(Sliced,PetscViewer);
512: EXTERN PetscErrorCode   SlicedCreate(MPI_Comm,Sliced*);
513: EXTERN PetscErrorCode   SlicedDestroy(Sliced);
514: EXTERN PetscErrorCode   SlicedCreateGlobalVector(Sliced,Vec*);
515: EXTERN PetscErrorCode   SlicedCreateLocalVector(Sliced,Vec*);
516: EXTERN PetscErrorCode   SlicedGetMatrix(Sliced, const MatType,Mat*);
517: EXTERN PetscErrorCode   SlicedGetGlobalIndices(Sliced,PetscInt*[]);
518: EXTERN PetscErrorCode   SlicedSetPreallocation(Sliced,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
519: EXTERN PetscErrorCode   SlicedSetGhosts(Sliced,PetscInt,PetscInt,PetscInt,const PetscInt[]);
520: EXTERN PetscErrorCode   SlicedGlobalToLocalBegin(Sliced,Vec,InsertMode,Vec);
521: EXTERN PetscErrorCode   SlicedGlobalToLocalEnd(Sliced,Vec,InsertMode,Vec);


524: typedef struct NLF_DAAD* NLF;

526: #include <petscbag.h>

528: EXTERN PetscErrorCode  PetscViewerBinaryMatlabOpen(MPI_Comm, const char [], PetscViewer*);
529: EXTERN PetscErrorCode  PetscViewerBinaryMatlabDestroy(PetscViewer);
530: EXTERN PetscErrorCode  PetscViewerBinaryMatlabOutputBag(PetscViewer, const char [], PetscBag);
531: EXTERN PetscErrorCode  PetscViewerBinaryMatlabOutputVec(PetscViewer, const char [], Vec);
532: EXTERN PetscErrorCode  PetscViewerBinaryMatlabOutputVecDA(PetscViewer, const char [], Vec, DA);


535: /*S
536:   ADDA - Abstract PETSc object that manages distributed field data for a single structured grid
537:          These are for any number of dimensions.

539:   Level: advanced. 

541:   Concepts: distributed array
542: .seealso: DA, DACreate(), ADDACreate()
543: S*/
544: typedef struct _p_ADDA* ADDA;


548: PetscErrorCode  ADDACreate(MPI_Comm,PetscInt,PetscInt*,PetscInt*,PetscInt,PetscTruth*,ADDA*);
549: PetscErrorCode  ADDADestroy(ADDA);

551: /* DM interface functions */
552: PetscErrorCode  ADDAView(ADDA,PetscViewer);
553: PetscErrorCode  ADDACreateGlobalVector(ADDA,Vec*);
554: PetscErrorCode  ADDAGetColoring(ADDA,ISColoringType,ISColoring*);
555: PetscErrorCode  ADDAGetMatrix(ADDA,const MatType, Mat*);
556: PetscErrorCode  ADDAGetInterpolation(ADDA,ADDA,Mat*,Vec*);
557: PetscErrorCode  ADDARefine(ADDA, MPI_Comm,ADDA *);
558: PetscErrorCode  ADDACoarsen(ADDA, MPI_Comm, ADDA*);
559: PetscErrorCode  ADDAGetInjection(ADDA, ADDA, VecScatter*);
560: PetscErrorCode  ADDAGetAggregates(ADDA, ADDA, Mat *);

562: /* functions only supported by ADDA */
563: PetscErrorCode  ADDASetRefinement(ADDA, PetscInt *,PetscInt);
564: PetscErrorCode  ADDAGetCorners(ADDA, PetscInt **, PetscInt **);
565: PetscErrorCode  ADDAGetGhostCorners(ADDA, PetscInt **, PetscInt **);
566: PetscErrorCode  ADDAGetMatrixNS(ADDA, ADDA, const MatType , Mat *);

568: /* functions to set values in vectors and matrices */
569: struct _ADDAIdx_s {
570:   PetscInt     *x;               /* the coordinates, user has to make sure it is the correct size! */
571:   PetscInt     d;                /* indexes the dof */
572: };
573: typedef struct _ADDAIdx_s ADDAIdx;

575: PetscErrorCode  ADDAMatSetValues(Mat, ADDA, PetscInt, const ADDAIdx[], ADDA, PetscInt,
576:                                                   const ADDAIdx[], const PetscScalar[], InsertMode);

578: PetscTruth ADDAHCiterStartup(const PetscInt, const PetscInt *const, const PetscInt *const, PetscInt *const);
579: PetscTruth ADDAHCiter(const PetscInt, const PetscInt *const, const PetscInt *const, PetscInt *const);

582: #endif