152849c42SDave May 252849c42SDave May #ifndef __DATA_BUCKET_H__ 352849c42SDave May #define __DATA_BUCKET_H__ 452849c42SDave May 552849c42SDave May 652849c42SDave May #include <petsc.h> 74b46c5e1SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 852849c42SDave May 952849c42SDave May 1052849c42SDave May #define DATAFIELD_POINT_ACCESS_GUARD 1152849c42SDave May 1252849c42SDave May /* Logging flag */ 134b46c5e1SDave May #define DATA_BUCKET_LOG 1452849c42SDave May 1552849c42SDave May 1652849c42SDave May typedef enum { DATABUCKET_VIEW_STDOUT=0, DATABUCKET_VIEW_ASCII, DATABUCKET_VIEW_BINARY, DATABUCKET_VIEW_HDF5 } DataBucketViewType; 1752849c42SDave May 1852849c42SDave May 1952849c42SDave May 2052849c42SDave May struct _p_DataField { 2152849c42SDave May char *registeration_function; 22*52c3ed93SDave May PetscInt L,bs; 2352849c42SDave May PetscBool active; 2452849c42SDave May size_t atomic_size; 2552849c42SDave May char *name; /* what are they called */ 2652849c42SDave May void *data; /* the data - an array of structs */ 274b46c5e1SDave May PetscDataType petsc_type; 2852849c42SDave May }; 2952849c42SDave May 3052849c42SDave May struct _p_DataBucket { 315c18a9d6SDave May PetscInt L; /* number in use */ 325c18a9d6SDave May PetscInt buffer; /* memory buffer used for re-allocation */ 335c18a9d6SDave May PetscInt allocated; /* number allocated, this will equal datafield->L */ 3452849c42SDave May PetscBool finalised; /* DEPRECIATED */ 355c18a9d6SDave May PetscInt nfields; /* how many fields of this type */ 3652849c42SDave May DataField *field; /* the data */ 3752849c42SDave May }; 3852849c42SDave May 3952849c42SDave May #define __DATATFIELD_point_access(data,index,atomic_size) (void*)((char*)(data) + (index)*(atomic_size)) 4052849c42SDave May #define __DATATFIELD_point_access_offset(data,index,atomic_size,offset) (void*)((char*)(data) + (index)*(atomic_size) + (offset)) 4152849c42SDave May 4252849c42SDave May 4352849c42SDave May 442eac95f8SDave May PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val); 452eac95f8SDave May PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index); 4652849c42SDave May 472eac95f8SDave May PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF); 482eac95f8SDave May PetscErrorCode DataFieldDestroy(DataField *DF); 492eac95f8SDave May PetscErrorCode DataBucketCreate(DataBucket *DB); 502eac95f8SDave May PetscErrorCode DataBucketDestroy(DataBucket *DB); 512eac95f8SDave May PetscErrorCode DataBucketRegisterField( 5252849c42SDave May DataBucket db, 5352849c42SDave May const char registeration_function[], 5452849c42SDave May const char field_name[], 5552849c42SDave May size_t atomic_size,DataField *_gfield ); 5652849c42SDave May 5752849c42SDave May 587fbf63aeSDave May PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum); 59*52c3ed93SDave May PetscErrorCode DataFieldSetBlockSize(DataField df,PetscInt blocksize); 607fbf63aeSDave May PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L); 617fbf63aeSDave May PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end); 627fbf63aeSDave May PetscErrorCode DataFieldGetAccess(const DataField gfield); 637fbf63aeSDave May PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p); 647fbf63aeSDave May PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p); 657fbf63aeSDave May PetscErrorCode DataFieldRestoreAccess(DataField gfield); 667fbf63aeSDave May PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size); 677fbf63aeSDave May PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size); 6852849c42SDave May 697fbf63aeSDave May PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data); 707fbf63aeSDave May PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data); 7152849c42SDave May 727fbf63aeSDave May PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx); 73a2d1c729SDave May PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,const PetscInt pid_y,const DataField field_y); 747fbf63aeSDave May PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index); 7552849c42SDave May 762eac95f8SDave May PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield); 777fbf63aeSDave May PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found); 787fbf63aeSDave May PetscErrorCode DataBucketFinalize(DataBucket db); 797fbf63aeSDave May PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer); 807fbf63aeSDave May PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer); 817fbf63aeSDave May PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated); 827fbf63aeSDave May PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated); 837fbf63aeSDave May PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[]); 8452849c42SDave May 85a2d1c729SDave May PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,const DataBucket yb,const PetscInt pid_y); 867fbf63aeSDave May PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB); 877fbf63aeSDave May PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index); 8852849c42SDave May 897fbf63aeSDave May //PetscErrorCode DataBucketLoadFromFile(const char filename[],DataBucketViewType type,DataBucket *db); 907fbf63aeSDave May PetscErrorCode DataBucketLoadFromFile(MPI_Comm comm,const char filename[],DataBucketViewType type,DataBucket *db); 917fbf63aeSDave May //PetscErrorCode DataBucketView(DataBucket db,const char filename[],DataBucketViewType type); 927fbf63aeSDave May PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type); 9352849c42SDave May 947fbf63aeSDave May PetscErrorCode DataBucketAddPoint(DataBucket db); 957fbf63aeSDave May PetscErrorCode DataBucketRemovePoint(DataBucket db); 967fbf63aeSDave May PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index); 9752849c42SDave May 987fbf63aeSDave May PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB); 997fbf63aeSDave May PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2); 10052849c42SDave May 10152849c42SDave May /* helpers for parallel send/recv */ 1027fbf63aeSDave May PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf); 1037fbf63aeSDave May PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf); 1047fbf63aeSDave May PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf); 1057fbf63aeSDave May PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data); 10652849c42SDave May 10752849c42SDave May 10852849c42SDave May #endif 10952849c42SDave May 110