xref: /petsc/src/dm/impls/swarm/data_bucket.h (revision e600fa544e2bb197ca2af9b6e65ea465976dec56)
1 #if !defined(__DMSWARM_DATA_BUCKET_H__)
2 #define __DMSWARM_DATA_BUCKET_H__
3 
4 #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
5 
6 #define DMSWARM_DATA_BUCKET_BUFFER_DEFAULT -1
7 #define DMSWARM_DATAFIELD_POINT_ACCESS_GUARD
8 
9 /* Logging flag */
10 #define DMSWARM_DATA_BUCKET_LOG
11 
12 typedef enum { DATABUCKET_VIEW_STDOUT=0, DATABUCKET_VIEW_ASCII, DATABUCKET_VIEW_BINARY, DATABUCKET_VIEW_HDF5 } DMSwarmDataBucketViewType;
13 
14 struct _p_DMSwarmDataField {
15         char          *registration_function;
16         PetscInt      L,bs;
17         PetscBool     active;
18         size_t        atomic_size;
19         char          *name; /* what are they called */
20         void          *data; /* the data - an array of structs */
21   PetscDataType petsc_type;
22 };
23 
24 struct _p_DMSwarmDataBucket {
25         PetscInt  L;             /* number in use */
26         PetscInt  buffer;        /* memory buffer used for re-allocation */
27         PetscInt  allocated;     /* number allocated, this will equal datafield->L */
28         PetscBool finalised;     /* DEPRECATED */
29         PetscInt  nfields;       /* how many fields of this type */
30         DMSwarmDataField *field; /* the data */
31 };
32 
33 #define DMSWARM_DATAFIELD_point_access(data,index,atomic_size) (void*)((char*)(data) + (index)*(atomic_size))
34 #define DMSWARM_DATAFIELD_point_access_offset(data,index,atomic_size,offset) (void*)((char*)(data) + (index)*(atomic_size) + (offset))
35 
36 PETSC_INTERN PetscErrorCode DMSwarmDataFieldStringInList(const char [],const PetscInt ,const DMSwarmDataField [],PetscBool *);
37 PETSC_INTERN PetscErrorCode DMSwarmDataFieldStringFindInList(const char [],const PetscInt,const DMSwarmDataField [],PetscInt *);
38 
39 PETSC_INTERN PetscErrorCode DMSwarmDataFieldCreate(const char [],const char [],const size_t,const PetscInt,DMSwarmDataField *);
40 PETSC_INTERN PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *);
41 PETSC_INTERN PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *);
42 PETSC_INTERN PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *);
43 PETSC_INTERN PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket,PetscBool *);
44 PETSC_INTERN PetscErrorCode DMSwarmDataBucketRegisterField(DMSwarmDataBucket,const char [],const char [],size_t,DMSwarmDataField *);
45 
46 PETSC_INTERN PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField,PetscInt *);
47 PETSC_INTERN PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField,PetscInt);
48 PETSC_INTERN PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField,const PetscInt);
49 PETSC_INTERN PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField,const PetscInt,const PetscInt);
50 PETSC_INTERN PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField);
51 PETSC_INTERN PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField,const PetscInt,void **);
52 PETSC_INTERN PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField,const size_t,const PetscInt,void **);
53 PETSC_INTERN PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField);
54 PETSC_INTERN PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField,const size_t);
55 PETSC_INTERN PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField,size_t *);
56 
57 PETSC_INTERN PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField,void **);
58 PETSC_INTERN PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField,void **);
59 
60 PETSC_INTERN PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField,const PetscInt,const void *);
61 PETSC_INTERN PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt,const DMSwarmDataField,const PetscInt,const DMSwarmDataField);
62 PETSC_INTERN PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField,const PetscInt);
63 
64 PETSC_INTERN PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket,const char [],DMSwarmDataField *);
65 PETSC_INTERN PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket,const char [],PetscBool *);
66 PETSC_INTERN PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket);
67 PETSC_INTERN PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket,const PetscInt,const PetscInt);
68 PETSC_INTERN PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket,const PetscInt,const PetscInt);
69 PETSC_INTERN PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket,PetscInt *,PetscInt *,PetscInt *);
70 PETSC_INTERN PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm,DMSwarmDataBucket,PetscInt *,PetscInt *,PetscInt *);
71 PETSC_INTERN PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket,PetscInt *,DMSwarmDataField *[]);
72 
73 PETSC_INTERN PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket,const PetscInt,const DMSwarmDataBucket,const PetscInt);
74 PETSC_INTERN PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket,const PetscInt,const PetscInt [],DMSwarmDataBucket *);
75 PETSC_INTERN PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket,const PetscInt);
76 
77 PETSC_INTERN PetscErrorCode DMSwarmDataBucketView(MPI_Comm,DMSwarmDataBucket,const char [],DMSwarmDataBucketViewType);
78 
79 PETSC_INTERN PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket);
80 PETSC_INTERN PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket);
81 PETSC_INTERN PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket,const PetscInt);
82 
83 PETSC_INTERN PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket,DMSwarmDataBucket *);
84 PETSC_INTERN PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket,DMSwarmDataBucket);
85 
86 /* helpers for parallel send/recv */
87 PETSC_INTERN PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket,size_t *,void **);
88 PETSC_INTERN PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket,void **);
89 PETSC_INTERN PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket,const PetscInt,void *);
90 PETSC_INTERN PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket,const PetscInt,void *);
91 
92 #endif
93