Actual source code: dmimpl.h
petsc-3.11.3 2019-06-26
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>
9: #include <petsc/private/isimpl.h>
11: PETSC_EXTERN PetscBool DMRegisterAllCalled;
12: PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
13: typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt field, MatNullSpace *nullSpace);
15: typedef struct _DMOps *DMOps;
16: struct _DMOps {
17: PetscErrorCode (*view)(DM,PetscViewer);
18: PetscErrorCode (*load)(DM,PetscViewer);
19: PetscErrorCode (*clone)(DM,DM*);
20: PetscErrorCode (*setfromoptions)(PetscOptionItems*,DM);
21: PetscErrorCode (*setup)(DM);
22: PetscErrorCode (*createdefaultsection)(DM);
23: PetscErrorCode (*createdefaultconstraints)(DM);
24: PetscErrorCode (*createglobalvector)(DM,Vec*);
25: PetscErrorCode (*createlocalvector)(DM,Vec*);
26: PetscErrorCode (*getlocaltoglobalmapping)(DM);
27: PetscErrorCode (*createfieldis)(DM,PetscInt*,char***,IS**);
28: PetscErrorCode (*createcoordinatedm)(DM,DM*);
29: PetscErrorCode (*createcoordinatefield)(DM,DMField*);
31: PetscErrorCode (*getcoloring)(DM,ISColoringType,ISColoring*);
32: PetscErrorCode (*creatematrix)(DM, Mat*);
33: PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*);
34: PetscErrorCode (*createrestriction)(DM,DM,Mat*);
35: PetscErrorCode (*createmassmatrix)(DM,DM,Mat*);
36: PetscErrorCode (*getaggregates)(DM,DM,Mat*);
37: PetscErrorCode (*hascreateinjection)(DM,PetscBool*);
38: PetscErrorCode (*getinjection)(DM,DM,Mat*);
40: PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
41: PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*);
42: PetscErrorCode (*refinehierarchy)(DM,PetscInt,DM*);
43: PetscErrorCode (*coarsenhierarchy)(DM,PetscInt,DM*);
44: PetscErrorCode (*adaptlabel)(DM,DMLabel,DM*);
45: PetscErrorCode (*adaptmetric)(DM,Vec,DMLabel,DM*);
47: PetscErrorCode (*globaltolocalbegin)(DM,Vec,InsertMode,Vec);
48: PetscErrorCode (*globaltolocalend)(DM,Vec,InsertMode,Vec);
49: PetscErrorCode (*localtoglobalbegin)(DM,Vec,InsertMode,Vec);
50: PetscErrorCode (*localtoglobalend)(DM,Vec,InsertMode,Vec);
51: PetscErrorCode (*localtolocalbegin)(DM,Vec,InsertMode,Vec);
52: PetscErrorCode (*localtolocalend)(DM,Vec,InsertMode,Vec);
54: PetscErrorCode (*destroy)(DM);
56: PetscErrorCode (*computevariablebounds)(DM,Vec,Vec);
58: PetscErrorCode (*createsubdm)(DM,PetscInt,const PetscInt*,IS*,DM*);
59: PetscErrorCode (*createsuperdm)(DM*,PetscInt,IS**,DM*);
60: PetscErrorCode (*createfielddecomposition)(DM,PetscInt*,char***,IS**,DM**);
61: PetscErrorCode (*createdomaindecomposition)(DM,PetscInt*,char***,IS**,IS**,DM**);
62: PetscErrorCode (*createddscatters)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
64: PetscErrorCode (*getdimpoints)(DM,PetscInt,PetscInt*,PetscInt*);
65: PetscErrorCode (*locatepoints)(DM,Vec,DMPointLocationType,PetscSF);
66: PetscErrorCode (*getneighbors)(DM,PetscInt*,const PetscMPIInt**);
67: PetscErrorCode (*getboundingbox)(DM,PetscReal*,PetscReal*);
68: PetscErrorCode (*getlocalboundingbox)(DM,PetscReal*,PetscReal*);
69: PetscErrorCode (*locatepointssubdomain)(DM,Vec,PetscMPIInt**);
71: PetscErrorCode (*projectfunctionlocal)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
72: PetscErrorCode (*projectfunctionlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
73: 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);
74: 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);
75: PetscErrorCode (*computel2diff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
76: PetscErrorCode (*computel2gradientdiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal [],const PetscReal[],PetscInt, PetscScalar *,void *),void **,Vec,const PetscReal[],PetscReal *);
77: PetscErrorCode (*computel2fielddiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
79: PetscErrorCode (*getcompatibility)(DM,DM,PetscBool*,PetscBool*);
80: };
82: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
83: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
84: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
86: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
87: struct _DMCoarsenHookLink {
88: PetscErrorCode (*coarsenhook)(DM,DM,void*); /* Run once, when coarse DM is created */
89: PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*); /* Run each time a new problem is restricted to a coarse grid */
90: void *ctx;
91: DMCoarsenHookLink next;
92: };
94: typedef struct _DMRefineHookLink *DMRefineHookLink;
95: struct _DMRefineHookLink {
96: PetscErrorCode (*refinehook)(DM,DM,void*); /* Run once, when a fine DM is created */
97: PetscErrorCode (*interphook)(DM,Mat,DM,void*); /* Run each time a new problem is interpolated to a fine grid */
98: void *ctx;
99: DMRefineHookLink next;
100: };
102: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
103: struct _DMSubDomainHookLink {
104: PetscErrorCode (*ddhook)(DM,DM,void*);
105: PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*);
106: void *ctx;
107: DMSubDomainHookLink next;
108: };
110: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
111: struct _DMGlobalToLocalHookLink {
112: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
113: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
114: void *ctx;
115: DMGlobalToLocalHookLink next;
116: };
118: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
119: struct _DMLocalToGlobalHookLink {
120: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
121: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
122: void *ctx;
123: DMLocalToGlobalHookLink next;
124: };
126: typedef enum {DMVEC_STATUS_IN,DMVEC_STATUS_OUT} DMVecStatus;
127: typedef struct _DMNamedVecLink *DMNamedVecLink;
128: struct _DMNamedVecLink {
129: Vec X;
130: char *name;
131: DMVecStatus status;
132: DMNamedVecLink next;
133: };
135: typedef struct _DMWorkLink *DMWorkLink;
136: struct _DMWorkLink {
137: size_t bytes;
138: void *mem;
139: DMWorkLink next;
140: };
142: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
144: struct _n_DMLabelLink {
145: DMLabel label;
146: PetscBool output;
147: struct _n_DMLabelLink *next;
148: };
149: typedef struct _n_DMLabelLink *DMLabelLink;
151: struct _n_DMLabelLinkList {
152: PetscInt refct;
153: DMLabelLink next;
154: };
155: typedef struct _n_DMLabelLinkList *DMLabelLinkList;
157: typedef struct _n_Boundary *DMBoundary;
159: struct _n_Boundary {
160: DSBoundary dsboundary;
161: DMLabel label;
162: DMBoundary next;
163: };
165: typedef struct _n_Field {
166: PetscObject disc; /* Field discretization, or a PetscContainer with the field name */
167: DMLabel label; /* Label defining the domain of definition of the field */
168: PetscBool adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
169: } RegionField;
171: typedef struct _n_Space {
172: PetscDS ds;
173: DMLabel label;
174: } DMSpace;
176: PETSC_EXTERN PetscErrorCode DMDestroyLabelLinkList(DM);
178: struct _p_DM {
179: PETSCHEADER(struct _DMOps);
180: Vec localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
181: Vec globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
182: DMNamedVecLink namedglobal;
183: DMNamedVecLink namedlocal;
184: DMWorkLink workin,workout;
185: DMLabelLinkList labels; /* Linked list of labels */
186: DMLabel depthLabel; /* Optimized access to depth label */
187: void *ctx; /* a user context */
188: PetscErrorCode (*ctxdestroy)(void**);
189: Vec x; /* location at which the functions/Jacobian are computed */
190: ISColoringType coloringtype;
191: MatFDColoring fd;
192: VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
193: MatType mattype; /* type of matrix created with DMCreateMatrix() */
194: PetscInt bs;
195: ISLocalToGlobalMapping ltogmap;
196: PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
197: PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
198: 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 */
199: 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 */
200: PetscBool setfromoptionscalled;
201: void *data;
202: /* Hierarchy / Submeshes */
203: DM coarseMesh;
204: DM fineMesh;
205: DMCoarsenHookLink coarsenhook; /* For transfering auxiliary problem data to coarser grids */
206: DMRefineHookLink refinehook;
207: DMSubDomainHookLink subdomainhook;
208: DMGlobalToLocalHookLink gtolhook;
209: DMLocalToGlobalHookLink ltoghook;
210: /* Topology */
211: PetscInt dim; /* The topological dimension */
212: /* Flexible communication */
213: PetscSF sfMigration; /* SF for point distribution created during distribution */
214: PetscSF sf; /* SF for parallel point overlap */
215: PetscSF defaultSF; /* SF for parallel dof overlap using default section */
216: PetscSF sfNatural; /* SF mapping to the "natural" ordering */
217: PetscBool useNatural; /* Create the natural SF */
218: /* Allows a non-standard data layout */
219: PetscBool adjacency[2]; /* [use cone() or support() first, use the transitive closure] */
220: PetscSection defaultSection; /* Layout for local vectors */
221: PetscSection defaultGlobalSection; /* Layout for global vectors */
222: PetscLayout map;
223: /* Constraints */
224: PetscSection defaultConstraintSection;
225: Mat defaultConstraintMat;
226: /* Coordinates */
227: PetscInt dimEmbed; /* The dimension of the embedding space */
228: DM coordinateDM; /* Layout for coordinates (default section) */
229: Vec coordinates; /* Coordinate values in global vector */
230: Vec coordinatesLocal; /* Coordinate values in local vector */
231: PetscBool periodic; /* Is the DM periodic? */
232: DMField coordinateField; /* Coordinates as an abstract field */
233: PetscReal *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
234: DMBoundaryType *bdtype; /* Indicates type of topological boundary */
235: /* Null spaces -- of course I should make this have a variable number of fields */
236: /* I now believe this might not be the right way: see below */
237: NullSpaceFunc nullspaceConstructors[10];
238: NullSpaceFunc nearnullspaceConstructors[10];
239: /* Fields are represented by objects */
240: PetscInt Nf; /* Number of fields defined on the total domain */
241: RegionField *fields; /* Array of discretization fields with regions of validity */
242: DMBoundary boundary; /* List of boundary conditions */
243: PetscInt Nds; /* Number of discrete systems defined on the total domain */
244: DMSpace *probs; /* Array of discrete systems */
245: /* Output structures */
246: DM dmBC; /* The DM with boundary conditions in the global DM */
247: PetscInt outputSequenceNum; /* The current sequence number for output */
248: PetscReal outputSequenceVal; /* The current sequence value for output */
250: PetscObject dmksp,dmsnes,dmts;
251: };
253: PETSC_EXTERN PetscLogEvent DM_Convert;
254: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
255: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
256: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
257: PETSC_EXTERN PetscLogEvent DM_Coarsen;
258: PETSC_EXTERN PetscLogEvent DM_Refine;
259: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
260: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
262: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
263: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);
265: PETSC_EXTERN PetscErrorCode DMView_GLVis(DM,PetscViewer,PetscErrorCode(*)(DM,PetscViewer));
267: /*
269: Composite Vectors
271: Single global representation
272: Individual global representations
273: Single local representation
274: Individual local representations
276: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
278: DM da_u, da_v, da_p
280: DM dm_velocities
282: DM dm
284: DMDACreate(,&da_u);
285: DMDACreate(,&da_v);
286: DMCompositeCreate(,&dm_velocities);
287: DMCompositeAddDM(dm_velocities,(DM)du);
288: DMCompositeAddDM(dm_velocities,(DM)dv);
290: DMDACreate(,&da_p);
291: DMCompositeCreate(,&dm_velocities);
292: DMCompositeAddDM(dm,(DM)dm_velocities);
293: DMCompositeAddDM(dm,(DM)dm_p);
296: Access parts of composite vectors (Composite only)
297: ---------------------------------
298: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
299: ADD for local vector -
301: Element access
302: --------------
303: From global vectors
304: -DAVecGetArray - for DMDA
305: -VecGetArray - for DMSliced
306: ADD for DMComposite??? maybe
308: From individual vector
309: -DAVecGetArray - for DMDA
310: -VecGetArray -for sliced
311: ADD for DMComposite??? maybe
313: From single local vector
314: ADD * single local vector as arrays?
316: Communication
317: -------------
318: DMGlobalToLocal - global vector to single local vector
320: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
322: Obtaining vectors
323: -----------------
324: DMCreateGlobal/Local
325: DMGetGlobal/Local
326: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
329: ????? individual global vectors ????
331: */
333: #if defined(PETSC_HAVE_HDF5)
334: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
335: #endif
337: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
338: {
340: #if defined(PETSC_USE_DEBUG)
341: {
342: PetscInt dof;
344: *start = *end = 0; /* Silence overzealous compiler warning */
345: if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
346: PetscSectionGetOffset(dm->defaultSection, point, start);
347: PetscSectionGetDof(dm->defaultSection, point, &dof);
348: *end = *start + dof;
349: }
350: #else
351: {
352: const PetscSection s = dm->defaultSection;
353: *start = s->atlasOff[point - s->pStart];
354: *end = *start + s->atlasDof[point - s->pStart];
355: }
356: #endif
357: return(0);
358: }
360: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
361: {
363: #if defined(PETSC_USE_DEBUG)
364: {
365: PetscInt dof;
367: *start = *end = 0; /* Silence overzealous compiler warning */
368: if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
369: PetscSectionGetFieldOffset(dm->defaultSection, point, field, start);
370: PetscSectionGetFieldDof(dm->defaultSection, point, field, &dof);
371: *end = *start + dof;
372: }
373: #else
374: {
375: const PetscSection s = dm->defaultSection->field[field];
376: *start = s->atlasOff[point - s->pStart];
377: *end = *start + s->atlasDof[point - s->pStart];
378: }
379: #endif
380: return(0);
381: }
383: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
384: {
386: #if defined(PETSC_USE_DEBUG)
387: {
389: PetscInt dof,cdof;
390: *start = *end = 0; /* Silence overzealous compiler warning */
391: if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
392: if (!dm->defaultGlobalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default global Section. It will be created automatically by DMGetDefaultGlobalSection()");
393: PetscSectionGetOffset(dm->defaultGlobalSection, point, start);
394: PetscSectionGetDof(dm->defaultGlobalSection, point, &dof);
395: PetscSectionGetConstraintDof(dm->defaultGlobalSection, point, &cdof);
396: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
397: }
398: #else
399: {
400: const PetscSection s = dm->defaultGlobalSection;
401: const PetscInt dof = s->atlasDof[point - s->pStart];
402: const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
403: *start = s->atlasOff[point - s->pStart];
404: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
405: }
406: #endif
407: return(0);
408: }
410: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
411: {
413: #if defined(PETSC_USE_DEBUG)
414: {
415: PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
417: *start = *end = 0; /* Silence overzealous compiler warning */
418: if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
419: if (!dm->defaultGlobalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default global Section. It will be crated automatically by DMGetDefaultGlobalSection()");
420: PetscSectionGetOffset(dm->defaultGlobalSection, point, start);
421: PetscSectionGetOffset(dm->defaultSection, point, &loff);
422: PetscSectionGetFieldOffset(dm->defaultSection, point, field, &lfoff);
423: PetscSectionGetFieldDof(dm->defaultSection, point, field, &fdof);
424: PetscSectionGetFieldConstraintDof(dm->defaultSection, point, field, &fcdof);
425: *start = *start < 0 ? *start - (lfoff-loff) : *start + lfoff-loff;
426: for (f = 0; f < field; ++f) {
427: PetscSectionGetFieldConstraintDof(dm->defaultSection, point, f, &ffcdof);
428: *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
429: }
430: *end = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
431: }
432: #else
433: {
434: const PetscSection s = dm->defaultSection;
435: const PetscSection fs = dm->defaultSection->field[field];
436: const PetscSection gs = dm->defaultGlobalSection;
437: const PetscInt loff = s->atlasOff[point - s->pStart];
438: const PetscInt goff = gs->atlasOff[point - s->pStart];
439: const PetscInt lfoff = fs->atlasOff[point - s->pStart];
440: const PetscInt fdof = fs->atlasDof[point - s->pStart];
441: const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
442: PetscInt ffcdof = 0, f;
444: for (f = 0; f < field; ++f) {
445: const PetscSection ffs = dm->defaultSection->field[f];
446: ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
447: }
448: *start = goff + (goff < 0 ? loff-lfoff + ffcdof : lfoff-loff - ffcdof);
449: *end = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
450: }
451: #endif
452: return(0);
453: }
455: #endif