xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision dbe06d34ff6929b64e5adad1313d8fffaf0196ff)
152849c42SDave May 
252849c42SDave May #include "data_bucket.h"
352849c42SDave May 
452849c42SDave May /* string helpers */
52eac95f8SDave May PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val)
652849c42SDave May {
75c18a9d6SDave May 	PetscInt i;
852849c42SDave May 
952849c42SDave May 	*val = PETSC_FALSE;
1052849c42SDave May 	for (i=0; i<N; i++) {
1152849c42SDave May 		if (strcmp( name, gfield[i]->name) == 0 ) {
1252849c42SDave May 			*val = PETSC_TRUE;
132eac95f8SDave May       PetscFunctionReturn(0);
1452849c42SDave May 		}
1552849c42SDave May 	}
162eac95f8SDave May   PetscFunctionReturn(0);
1752849c42SDave May }
1852849c42SDave May 
192eac95f8SDave May PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index)
2052849c42SDave May {
215c18a9d6SDave May 	PetscInt i;
2252849c42SDave May 
2352849c42SDave May 	*index = -1;
2452849c42SDave May 	for (i=0; i<N; i++) {
2552849c42SDave May 		if (strcmp( name, gfield[i]->name ) == 0) {
2652849c42SDave May 			*index = i;
272eac95f8SDave May       PetscFunctionReturn(0);
2852849c42SDave May 		}
2952849c42SDave May 	}
302eac95f8SDave May   PetscFunctionReturn(0);
3152849c42SDave May }
3252849c42SDave May 
332eac95f8SDave May #undef __FUNCT__
342eac95f8SDave May #define __FUNCT__ "DataFieldCreate"
352eac95f8SDave May PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF)
3652849c42SDave May {
3752849c42SDave May 	DataField df;
3852849c42SDave May 
3952849c42SDave May 	df = malloc( sizeof(struct _p_DataField) );
4052849c42SDave May 	memset( df, 0, sizeof(struct _p_DataField) );
4152849c42SDave May 
4252849c42SDave May 
4352849c42SDave May 	asprintf( &df->registeration_function, "%s", registeration_function );
4452849c42SDave May 	asprintf( &df->name, "%s", name );
4552849c42SDave May 	df->atomic_size = size;
4652849c42SDave May 	df->L = L;
4752849c42SDave May 
4852849c42SDave May 	df->data = malloc( size * L ); /* allocate something so we don't have to reallocate */
4952849c42SDave May 	memset( df->data, 0, size * L );
5052849c42SDave May 
5152849c42SDave May 	*DF = df;
522eac95f8SDave May   PetscFunctionReturn(0);
5352849c42SDave May }
5452849c42SDave May 
552eac95f8SDave May #undef __FUNCT__
562eac95f8SDave May #define __FUNCT__ "DataFieldDestroy"
572eac95f8SDave May PetscErrorCode DataFieldDestroy(DataField *DF)
5852849c42SDave May {
5952849c42SDave May 	DataField df = *DF;
6052849c42SDave May 
6152849c42SDave May 	free(df->registeration_function);
6252849c42SDave May 	free(df->name);
6352849c42SDave May 	free(df->data);
6452849c42SDave May 	free(df);
6552849c42SDave May 
6652849c42SDave May 	*DF = NULL;
672eac95f8SDave May   PetscFunctionReturn(0);
6852849c42SDave May }
6952849c42SDave May 
7052849c42SDave May /* data bucket */
712eac95f8SDave May #undef __FUNCT__
722eac95f8SDave May #define __FUNCT__ "DataBucketCreate"
732eac95f8SDave May PetscErrorCode DataBucketCreate(DataBucket *DB)
7452849c42SDave May {
7552849c42SDave May 	DataBucket db;
7652849c42SDave May 
7752849c42SDave May 
7852849c42SDave May 	db = malloc( sizeof(struct _p_DataBucket) );
7952849c42SDave May 	memset( db, 0, sizeof(struct _p_DataBucket) );
8052849c42SDave May 
8152849c42SDave May 	db->finalised = PETSC_FALSE;
8252849c42SDave May 
8352849c42SDave May 	/* create empty spaces for fields */
8452849c42SDave May 	db->L         = 0;
8552849c42SDave May 	db->buffer    = 1;
8652849c42SDave May 	db->allocated = 1;
8752849c42SDave May 
8852849c42SDave May 	db->nfields   = 0;
8952849c42SDave May 	db->field     = malloc(sizeof(DataField));
9052849c42SDave May 
9152849c42SDave May 	*DB = db;
922eac95f8SDave May   PetscFunctionReturn(0);
9352849c42SDave May }
9452849c42SDave May 
952eac95f8SDave May #undef __FUNCT__
962eac95f8SDave May #define __FUNCT__ "DataBucketDestroy"
972eac95f8SDave May PetscErrorCode DataBucketDestroy(DataBucket *DB)
9852849c42SDave May {
9952849c42SDave May 	DataBucket db = *DB;
1005c18a9d6SDave May 	PetscInt f;
101*dbe06d34SDave May 	PetscErrorCode ierr;
10252849c42SDave May 
10352849c42SDave May 	/* release fields */
10452849c42SDave May 	for (f=0; f<db->nfields; f++) {
105*dbe06d34SDave May 		ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
10652849c42SDave May 	}
10752849c42SDave May 
10852849c42SDave May 	/* this will catch the initially allocated objects in the event that no fields are registered */
10952849c42SDave May 	if (db->field != NULL) {
11052849c42SDave May 		free(db->field);
11152849c42SDave May 	}
11252849c42SDave May 
11352849c42SDave May 	free(db);
11452849c42SDave May 
11552849c42SDave May 	*DB = NULL;
1162eac95f8SDave May   PetscFunctionReturn(0);
11752849c42SDave May }
11852849c42SDave May 
1192eac95f8SDave May #undef __FUNCT__
1202eac95f8SDave May #define __FUNCT__ "DataBucketRegisterField"
1212eac95f8SDave May PetscErrorCode DataBucketRegisterField(
12252849c42SDave May                               DataBucket db,
12352849c42SDave May                               const char registeration_function[],
12452849c42SDave May                               const char field_name[],
12552849c42SDave May                               size_t atomic_size, DataField *_gfield)
12652849c42SDave May {
12752849c42SDave May 	PetscBool val;
12852849c42SDave May 	DataField *field,fp;
129*dbe06d34SDave May 	PetscErrorCode ierr;
13052849c42SDave May 
13152849c42SDave May 	/* check we haven't finalised the registration of fields */
13252849c42SDave May 	/*
13352849c42SDave May    if(db->finalised==PETSC_TRUE) {
13452849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
13552849c42SDave May    ERROR();
13652849c42SDave May    }
13752849c42SDave May    */
13852849c42SDave May 
13952849c42SDave May 	/* check for repeated name */
14052849c42SDave May 	StringInList( field_name, db->nfields, (const DataField*)db->field, &val );
1412eac95f8SDave May 	if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
14252849c42SDave May 
14352849c42SDave May 	/* create new space for data */
14452849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
14552849c42SDave May 	db->field     = field;
14652849c42SDave May 
14752849c42SDave May 	/* add field */
148*dbe06d34SDave May 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, db->allocated, &fp );CHKERRQ(ierr);
14952849c42SDave May 	db->field[ db->nfields ] = fp;
15052849c42SDave May 
15152849c42SDave May 	db->nfields++;
15252849c42SDave May 
15352849c42SDave May 	if (_gfield != NULL) {
15452849c42SDave May 		*_gfield = fp;
15552849c42SDave May 	}
1562eac95f8SDave May   PetscFunctionReturn(0);
15752849c42SDave May }
15852849c42SDave May 
15952849c42SDave May /*
16052849c42SDave May  #define DataBucketRegisterField(db,name,size,k) {\
16152849c42SDave May  char *location;\
16252849c42SDave May  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
16352849c42SDave May  _DataBucketRegisterField( (db), location, (name), (size), (k) );\
16452849c42SDave May  free(location);\
16552849c42SDave May  }
16652849c42SDave May  */
16752849c42SDave May 
1682eac95f8SDave May #undef __FUNCT__
1692eac95f8SDave May #define __FUNCT__ "DataBucketGetDataFieldByName"
1702eac95f8SDave May PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
17152849c42SDave May {
1725c18a9d6SDave May 	PetscInt idx;
17352849c42SDave May 	PetscBool found;
17452849c42SDave May 
17552849c42SDave May 	StringInList(name,db->nfields,(const DataField*)db->field,&found);
1762eac95f8SDave May 	if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name);
1772eac95f8SDave May 
17852849c42SDave May 	StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);
17952849c42SDave May 
18052849c42SDave May 	*gfield = db->field[idx];
1812eac95f8SDave May   PetscFunctionReturn(0);
18252849c42SDave May }
18352849c42SDave May 
1847fbf63aeSDave May #undef __FUNCT__
1857fbf63aeSDave May #define __FUNCT__ "DataBucketQueryDataFieldByName"
1867fbf63aeSDave May PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
18752849c42SDave May {
18852849c42SDave May 	*found = PETSC_FALSE;
18952849c42SDave May 	StringInList(name,db->nfields,(const DataField*)db->field,found);
1907fbf63aeSDave May   PetscFunctionReturn(0);
19152849c42SDave May }
19252849c42SDave May 
1937fbf63aeSDave May #undef __FUNCT__
1947fbf63aeSDave May #define __FUNCT__ "DataBucketFinalize"
1957fbf63aeSDave May PetscErrorCode DataBucketFinalize(DataBucket db)
19652849c42SDave May {
19752849c42SDave May 	db->finalised = PETSC_TRUE;
1987fbf63aeSDave May   PetscFunctionReturn(0);
19952849c42SDave May }
20052849c42SDave May 
2017fbf63aeSDave May #undef __FUNCT__
2027fbf63aeSDave May #define __FUNCT__ "DataFieldGetNumEntries"
2037fbf63aeSDave May PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum)
20452849c42SDave May {
20552849c42SDave May 	*sum = df->L;
2067fbf63aeSDave May   PetscFunctionReturn(0);
20752849c42SDave May }
20852849c42SDave May 
2097fbf63aeSDave May #undef __FUNCT__
2107fbf63aeSDave May #define __FUNCT__ "DataFieldSetSize"
2117fbf63aeSDave May PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L)
21252849c42SDave May {
21352849c42SDave May 	void *tmp_data;
21452849c42SDave May 
215a233d522SDave May 	if (new_L <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be <= 0");
2167fbf63aeSDave May 
2177fbf63aeSDave May 	if (new_L == df->L) PetscFunctionReturn(0);
21852849c42SDave May 
21952849c42SDave May 	if (new_L > df->L) {
22052849c42SDave May 
22152849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L) );
22252849c42SDave May 		df->data = tmp_data;
22352849c42SDave May 
22452849c42SDave May 		/* init new contents */
22552849c42SDave May 		memset( ( ((char*)df->data)+df->L*df->atomic_size), 0, (new_L-df->L)*df->atomic_size );
22652849c42SDave May 
2277fbf63aeSDave May 	} else {
22852849c42SDave May 		/* reallocate pointer list, add +1 in case new_L = 0 */
22952849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L+1) );
23052849c42SDave May 		df->data = tmp_data;
23152849c42SDave May 	}
23252849c42SDave May 
23352849c42SDave May 	df->L = new_L;
2347fbf63aeSDave May   PetscFunctionReturn(0);
23552849c42SDave May }
23652849c42SDave May 
2377fbf63aeSDave May #undef __FUNCT__
2387fbf63aeSDave May #define __FUNCT__ "DataFieldZeroBlock"
2397fbf63aeSDave May PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end)
24052849c42SDave May {
2417fbf63aeSDave May 	if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
2427fbf63aeSDave May 
2437fbf63aeSDave May 	if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
2447fbf63aeSDave May 
245a233d522SDave May 	if (end > df->L) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if end(%D) >= array size(%D)",end,df->L);
24652849c42SDave May 
24752849c42SDave May 	memset( ( ((char*)df->data)+start*df->atomic_size), 0, (end-start)*df->atomic_size );
2487fbf63aeSDave May   PetscFunctionReturn(0);
24952849c42SDave May }
25052849c42SDave May 
25152849c42SDave May /*
25252849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
25352849c42SDave May  */
2547fbf63aeSDave May #undef __FUNCT__
2557fbf63aeSDave May #define __FUNCT__ "DataBucketSetSizes"
2567fbf63aeSDave May PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
25752849c42SDave May {
2585c18a9d6SDave May 	PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
259*dbe06d34SDave May   PetscErrorCode ierr;
26052849c42SDave May 
2617fbf63aeSDave May 	if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()");
26252849c42SDave May 
26352849c42SDave May 	current_allocated = db->allocated;
26452849c42SDave May 
26552849c42SDave May 	new_used   = L;
26652849c42SDave May 	new_unused = current_allocated - new_used;
26752849c42SDave May 	new_buffer = db->buffer;
26852849c42SDave May 	if (buffer >= 0) { /* update the buffer value */
26952849c42SDave May 		new_buffer = buffer;
27052849c42SDave May 	}
27152849c42SDave May 	new_allocated = new_used + new_buffer;
27252849c42SDave May 
27352849c42SDave May 	/* action */
27452849c42SDave May 	if (new_allocated > current_allocated) {
27552849c42SDave May 		/* increase size to new_used + new_buffer */
27652849c42SDave May 		for (f=0; f<db->nfields; f++) {
277*dbe06d34SDave May 			ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
27852849c42SDave May 		}
27952849c42SDave May 
28052849c42SDave May 		db->L         = new_used;
28152849c42SDave May 		db->buffer    = new_buffer;
28252849c42SDave May 		db->allocated = new_used + new_buffer;
28352849c42SDave May 	}
28452849c42SDave May 	else {
28552849c42SDave May 		if (new_unused > 2 * new_buffer) {
28652849c42SDave May 
28752849c42SDave May 			/* shrink array to new_used + new_buffer */
28852849c42SDave May 			for (f=0; f<db->nfields; f++) {
289*dbe06d34SDave May 				ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
29052849c42SDave May 			}
29152849c42SDave May 
29252849c42SDave May 			db->L         = new_used;
29352849c42SDave May 			db->buffer    = new_buffer;
29452849c42SDave May 			db->allocated = new_used + new_buffer;
29552849c42SDave May 		}
29652849c42SDave May 		else {
29752849c42SDave May 			db->L      = new_used;
29852849c42SDave May 			db->buffer = new_buffer;
29952849c42SDave May 		}
30052849c42SDave May 	}
30152849c42SDave May 
30252849c42SDave May 	/* zero all entries from db->L to db->allocated */
30352849c42SDave May 	for (f=0; f<db->nfields; f++) {
30452849c42SDave May 		DataField field = db->field[f];
305*dbe06d34SDave May 		ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
30652849c42SDave May 	}
3077fbf63aeSDave May   PetscFunctionReturn(0);
30852849c42SDave May }
30952849c42SDave May 
3107fbf63aeSDave May #undef __FUNCT__
3117fbf63aeSDave May #define __FUNCT__ "DataBucketSetInitialSizes"
3127fbf63aeSDave May PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
31352849c42SDave May {
3145c18a9d6SDave May 	PetscInt f;
315*dbe06d34SDave May   PetscErrorCode ierr;
316*dbe06d34SDave May 
317*dbe06d34SDave May 	ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
31852849c42SDave May 
31952849c42SDave May 	for (f=0; f<db->nfields; f++) {
32052849c42SDave May 		DataField field = db->field[f];
321*dbe06d34SDave May 		ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
32252849c42SDave May 	}
3237fbf63aeSDave May   PetscFunctionReturn(0);
32452849c42SDave May }
32552849c42SDave May 
3267fbf63aeSDave May #undef __FUNCT__
3277fbf63aeSDave May #define __FUNCT__ "DataBucketGetSizes"
3287fbf63aeSDave May PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
32952849c42SDave May {
33052849c42SDave May 	if (L) { *L = db->L; }
33152849c42SDave May 	if (buffer) { *buffer = db->buffer; }
33252849c42SDave May 	if (allocated) { *allocated = db->allocated; }
3337fbf63aeSDave May   PetscFunctionReturn(0);
33452849c42SDave May }
33552849c42SDave May 
3367fbf63aeSDave May #undef __FUNCT__
3377fbf63aeSDave May #define __FUNCT__ "DataBucketGetGlobalSizes"
3387fbf63aeSDave May PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
33952849c42SDave May {
3405c18a9d6SDave May 	PetscInt _L,_buffer,_allocated;
3415c18a9d6SDave May 	PetscInt ierr;
34252849c42SDave May 
3435c18a9d6SDave May 	_L = db->L;
3445c18a9d6SDave May 	_buffer = db->buffer;
3455c18a9d6SDave May 	_allocated = db->allocated;
34652849c42SDave May 
3475c18a9d6SDave May 	if (L) {         ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm); }
3485c18a9d6SDave May 	if (buffer) {    ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm); }
3495c18a9d6SDave May 	if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm); }
3507fbf63aeSDave May   PetscFunctionReturn(0);
35152849c42SDave May }
35252849c42SDave May 
3537fbf63aeSDave May #undef __FUNCT__
3547fbf63aeSDave May #define __FUNCT__ "DataBucketGetGlobalSizes"
3557fbf63aeSDave May PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[])
35652849c42SDave May {
35752849c42SDave May 	if (L) {      *L      = db->nfields; }
35852849c42SDave May 	if (fields) { *fields = db->field; }
3597fbf63aeSDave May   PetscFunctionReturn(0);
36052849c42SDave May }
36152849c42SDave May 
3627fbf63aeSDave May #undef __FUNCT__
3637fbf63aeSDave May #define __FUNCT__ "DataFieldGetAccess"
3647fbf63aeSDave May PetscErrorCode DataFieldGetAccess(const DataField gfield)
36552849c42SDave May {
3667fbf63aeSDave May 	if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name);
3677fbf63aeSDave May 
36852849c42SDave May 	gfield->active = PETSC_TRUE;
3697fbf63aeSDave May   PetscFunctionReturn(0);
37052849c42SDave May }
37152849c42SDave May 
3727fbf63aeSDave May #undef __FUNCT__
3737fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPoint"
3747fbf63aeSDave May PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p)
37552849c42SDave May {
37652849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
37752849c42SDave May 	/* debug mode */
3785c18a9d6SDave May 	/* check poPetscInt is valid */
3797fbf63aeSDave May 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
3807fbf63aeSDave May 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
38152849c42SDave May 
3827fbf63aeSDave May 	if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess() before poPetscInt data can be retrivied",gfield->name);
38352849c42SDave May #endif
38452849c42SDave May 
38552849c42SDave May 	//*ctx_p  = (void*)( ((char*)gfield->data) + pid * gfield->atomic_size);
38652849c42SDave May 	*ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
3877fbf63aeSDave May   PetscFunctionReturn(0);
38852849c42SDave May }
38952849c42SDave May 
3907fbf63aeSDave May #undef __FUNCT__
3917fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPointOffset"
3927fbf63aeSDave May PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
39352849c42SDave May {
39452849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
39552849c42SDave May 	/* debug mode */
39652849c42SDave May 
3975c18a9d6SDave May 	/* check poPetscInt is valid */
3987fbf63aeSDave May 	/* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*//* Note compiler realizes this can never happen with an unsigned PetscInt */
3997fbf63aeSDave May 	if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
40052849c42SDave May 
4015c18a9d6SDave May 	/* check poPetscInt is valid */
4027fbf63aeSDave May 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
403a233d522SDave May 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
40452849c42SDave May 
405a233d522SDave May 	if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess() before point data can be retrivied",gfield->name);
40652849c42SDave May #endif
40752849c42SDave May 
40852849c42SDave May 	*ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
4097fbf63aeSDave May   PetscFunctionReturn(0);
41052849c42SDave May }
41152849c42SDave May 
4127fbf63aeSDave May #undef __FUNCT__
4137fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPointOffset"
4147fbf63aeSDave May PetscErrorCode DataFieldRestoreAccess(DataField gfield)
41552849c42SDave May {
4167fbf63aeSDave May 	if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name );
4177fbf63aeSDave May 
41852849c42SDave May 	gfield->active = PETSC_FALSE;
4197fbf63aeSDave May   PetscFunctionReturn(0);
42052849c42SDave May }
42152849c42SDave May 
4227fbf63aeSDave May #undef __FUNCT__
4237fbf63aeSDave May #define __FUNCT__ "DataFieldVerifyAccess"
4247fbf63aeSDave May PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size)
42552849c42SDave May {
42652849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
4277fbf63aeSDave May 	if (gfield->atomic_size != size) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" must be mapped to %zu bytes, your intended structure is %zu bytes in length.",
42852849c42SDave May            gfield->name, gfield->atomic_size, size );
42952849c42SDave May #endif
4307fbf63aeSDave May   PetscFunctionReturn(0);
43152849c42SDave May }
43252849c42SDave May 
4337fbf63aeSDave May #undef __FUNCT__
4347fbf63aeSDave May #define __FUNCT__ "DataFieldGetAtomicSize"
4357fbf63aeSDave May PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size)
43652849c42SDave May {
43752849c42SDave May   if (size) { *size = gfield->atomic_size; }
4387fbf63aeSDave May   PetscFunctionReturn(0);
43952849c42SDave May }
44052849c42SDave May 
4417fbf63aeSDave May #undef __FUNCT__
4427fbf63aeSDave May #define __FUNCT__ "DataFieldGetEntries"
4437fbf63aeSDave May PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data)
44452849c42SDave May {
44552849c42SDave May   if (data) {
44652849c42SDave May     *data = gfield->data;
44752849c42SDave May   }
4487fbf63aeSDave May   PetscFunctionReturn(0);
44952849c42SDave May }
45052849c42SDave May 
4517fbf63aeSDave May #undef __FUNCT__
4527fbf63aeSDave May #define __FUNCT__ "DataFieldRestoreEntries"
4537fbf63aeSDave May PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data)
45452849c42SDave May {
45552849c42SDave May   if (data) {
45652849c42SDave May     *data = NULL;
45752849c42SDave May   }
4587fbf63aeSDave May   PetscFunctionReturn(0);
45952849c42SDave May }
46052849c42SDave May 
46152849c42SDave May /* y = x */
4627fbf63aeSDave May #undef __FUNCT__
4637fbf63aeSDave May #define __FUNCT__ "DataBucketCopyPoint"
4647fbf63aeSDave May PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,
4655c18a9d6SDave May                          const DataBucket yb,const PetscInt pid_y)
46652849c42SDave May {
4675c18a9d6SDave May 	PetscInt f;
468*dbe06d34SDave May 	PetscErrorCode ierr;
469*dbe06d34SDave May 
47052849c42SDave May 	for (f=0; f<xb->nfields; f++) {
47152849c42SDave May 		void *dest;
47252849c42SDave May 		void *src;
47352849c42SDave May 
474*dbe06d34SDave May 		ierr = DataFieldGetAccess( xb->field[f] );CHKERRQ(ierr);
47552849c42SDave May 		if (xb != yb) { DataFieldGetAccess( yb->field[f] ); }
47652849c42SDave May 
477*dbe06d34SDave May 		ierr = DataFieldAccessPoint( xb->field[f],pid_x, &src );CHKERRQ(ierr);
478*dbe06d34SDave May 		ierr = DataFieldAccessPoint( yb->field[f],pid_y, &dest );CHKERRQ(ierr);
47952849c42SDave May 
48052849c42SDave May 		memcpy( dest, src, xb->field[f]->atomic_size );
48152849c42SDave May 
482*dbe06d34SDave May 		ierr = DataFieldRestoreAccess( xb->field[f] );CHKERRQ(ierr);
48352849c42SDave May 		if (xb != yb) { DataFieldRestoreAccess( yb->field[f] ); }
48452849c42SDave May 	}
4857fbf63aeSDave May   PetscFunctionReturn(0);
48652849c42SDave May }
48752849c42SDave May 
4887fbf63aeSDave May #undef __FUNCT__
4897fbf63aeSDave May #define __FUNCT__ "DataBucketCreateFromSubset"
4907fbf63aeSDave May PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB)
49152849c42SDave May {
4925c18a9d6SDave May 	PetscInt nfields;
49352849c42SDave May 	DataField *fields;
49452849c42SDave May 	DataBucketCreate(DB);
4955c18a9d6SDave May 	PetscInt f,L,buffer,allocated,p;
496*dbe06d34SDave May 	PetscErrorCode ierr;
49752849c42SDave May 
49852849c42SDave May 	/* copy contents of DBIn */
499*dbe06d34SDave May 	ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
500*dbe06d34SDave May 	ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
50152849c42SDave May 
50252849c42SDave May 	for (f=0; f<nfields; f++) {
503*dbe06d34SDave May 		ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
50452849c42SDave May 	}
505*dbe06d34SDave May 	ierr = DataBucketFinalize(*DB);CHKERRQ(ierr);
50652849c42SDave May 
507*dbe06d34SDave May 	ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
50852849c42SDave May 
50952849c42SDave May 	/* now copy the desired guys from DBIn => DB */
51052849c42SDave May 	for (p=0; p<N; p++) {
511*dbe06d34SDave May 		ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
51252849c42SDave May 	}
5137fbf63aeSDave May   PetscFunctionReturn(0);
51452849c42SDave May }
51552849c42SDave May 
51652849c42SDave May // insert into an exisitng location
5177fbf63aeSDave May #undef __FUNCT__
5187fbf63aeSDave May #define __FUNCT__ "DataFieldInsertPoint"
5197fbf63aeSDave May PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx)
52052849c42SDave May {
52152849c42SDave May 
52252849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
5235c18a9d6SDave May 	/* check poPetscInt is valid */
524a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
525a233d522SDave May 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
52652849c42SDave May #endif
52752849c42SDave May 
52852849c42SDave May   //	memcpy( (void*)((char*)field->data + index*field->atomic_size), ctx, field->atomic_size );
52952849c42SDave May 	memcpy( __DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size );
5307fbf63aeSDave May   PetscFunctionReturn(0);
53152849c42SDave May }
53252849c42SDave May 
53352849c42SDave May // remove data at index - replace with last point
534a233d522SDave May #undef __FUNCT__
535a233d522SDave May #define __FUNCT__ "DataBucketRemovePointAtIndex"
5367fbf63aeSDave May PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index)
53752849c42SDave May {
5385c18a9d6SDave May 	PetscInt f;
539*dbe06d34SDave May 	PetscErrorCode ierr;
54052849c42SDave May 
54152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
5425c18a9d6SDave May 	/* check poPetscInt is valid */
543a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
544a233d522SDave May 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
54552849c42SDave May #endif
54652849c42SDave May 
547a233d522SDave May 	if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
548a233d522SDave May 		SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"You should not be trying to remove point at index=%D since it's < db->L = %D", index, db->L );
54952849c42SDave May 	}
55052849c42SDave May 
551a233d522SDave May 	if (index != db->L-1) { /* not last point in list */
55252849c42SDave May 		for (f=0; f<db->nfields; f++) {
55352849c42SDave May 			DataField field = db->field[f];
55452849c42SDave May 
55552849c42SDave May 			/* copy then remove */
556*dbe06d34SDave May 			ierr = DataFieldCopyPoint( db->L-1,field, index,field );CHKERRQ(ierr);
55752849c42SDave May 
55852849c42SDave May 			//DataFieldZeroPoint(field,index);
55952849c42SDave May 		}
56052849c42SDave May 	}
56152849c42SDave May 
56252849c42SDave May 	/* decrement size */
56352849c42SDave May 	/* this will zero out an crap at the end of the list */
564*dbe06d34SDave May 	ierr = DataBucketRemovePoint(db);CHKERRQ(ierr);
5657fbf63aeSDave May   PetscFunctionReturn(0);
56652849c42SDave May }
56752849c42SDave May 
56852849c42SDave May /* copy x into y */
5697fbf63aeSDave May #undef __FUNCT__
5707fbf63aeSDave May #define __FUNCT__ "DataFieldCopyPoint"
5717fbf63aeSDave May PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,
5725c18a9d6SDave May                         const PetscInt pid_y,const DataField field_y )
57352849c42SDave May {
57452849c42SDave May 
57552849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
5765c18a9d6SDave May 	/* check poPetscInt is valid */
577a233d522SDave May 	if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
578a233d522SDave May 	if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
57952849c42SDave May 
580a233d522SDave May 	if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
581a233d522SDave May 	if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
58252849c42SDave May 
58352849c42SDave May 	if( field_y->atomic_size != field_x->atomic_size ) {
584a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
58552849c42SDave May 	}
58652849c42SDave May #endif
58752849c42SDave May 	/*
58852849c42SDave May    memcpy( (void*)((char*)field_y->data + pid_y*field_y->atomic_size),
58952849c42SDave May    (void*)((char*)field_x->data + pid_x*field_x->atomic_size),
59052849c42SDave May    field_x->atomic_size );
59152849c42SDave May    */
59252849c42SDave May 	memcpy(		__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
59352849c42SDave May          __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
59452849c42SDave May          field_y->atomic_size );
5957fbf63aeSDave May   PetscFunctionReturn(0);
59652849c42SDave May }
59752849c42SDave May 
59852849c42SDave May 
59952849c42SDave May // zero only the datafield at this point
6007fbf63aeSDave May #undef __FUNCT__
6017fbf63aeSDave May #define __FUNCT__ "DataFieldZeroPoint"
6027fbf63aeSDave May PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index)
60352849c42SDave May {
60452849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
6055c18a9d6SDave May 	/* check poPetscInt is valid */
606a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
607a233d522SDave May 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
60852849c42SDave May #endif
60952849c42SDave May 
61052849c42SDave May   //	memset( (void*)((char*)field->data + index*field->atomic_size), 0, field->atomic_size );
61152849c42SDave May 	memset( __DATATFIELD_point_access(field->data,index,field->atomic_size), 0, field->atomic_size );
6127fbf63aeSDave May   PetscFunctionReturn(0);
61352849c42SDave May }
61452849c42SDave May 
61552849c42SDave May // zero ALL data for this point
6167fbf63aeSDave May #undef __FUNCT__
6177fbf63aeSDave May #define __FUNCT__ "DataBucketZeroPoint"
6187fbf63aeSDave May PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index)
61952849c42SDave May {
6205c18a9d6SDave May 	PetscInt f;
621*dbe06d34SDave May 	PetscErrorCode ierr;
62252849c42SDave May 
6235c18a9d6SDave May 	/* check poPetscInt is valid */
624a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
625a233d522SDave May 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
62652849c42SDave May 
62752849c42SDave May 	for (f=0; f<db->nfields; f++) {
62852849c42SDave May 		DataField field = db->field[f];
62952849c42SDave May 
630*dbe06d34SDave May 		ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr);
63152849c42SDave May 	}
6327fbf63aeSDave May   PetscFunctionReturn(0);
63352849c42SDave May }
63452849c42SDave May 
63552849c42SDave May /* increment */
6367fbf63aeSDave May #undef __FUNCT__
6377fbf63aeSDave May #define __FUNCT__ "DataBucketAddPoint"
6387fbf63aeSDave May PetscErrorCode DataBucketAddPoint(DataBucket db)
63952849c42SDave May {
640*dbe06d34SDave May 	PetscErrorCode ierr;
641*dbe06d34SDave May 
642*dbe06d34SDave May 	ierr = DataBucketSetSizes(db,db->L+1,-1);CHKERRQ(ierr);
6437fbf63aeSDave May   PetscFunctionReturn(0);
64452849c42SDave May }
64552849c42SDave May 
6467fbf63aeSDave May /* decrement */
6477fbf63aeSDave May #undef __FUNCT__
6487fbf63aeSDave May #define __FUNCT__ "DataBucketRemovePoint"
6497fbf63aeSDave May PetscErrorCode DataBucketRemovePoint(DataBucket db)
6507fbf63aeSDave May {
651*dbe06d34SDave May 	PetscErrorCode ierr;
652*dbe06d34SDave May 
653*dbe06d34SDave May 	ierr = DataBucketSetSizes(db,db->L-1,-1);CHKERRQ(ierr);
6547fbf63aeSDave May   PetscFunctionReturn(0);
6557fbf63aeSDave May }
6567fbf63aeSDave May 
6577fbf63aeSDave May #undef __FUNCT__
6587fbf63aeSDave May #define __FUNCT__ "_DataFieldViewBinary"
6597fbf63aeSDave May PetscErrorCode _DataFieldViewBinary(DataField field,FILE *fp)
66052849c42SDave May {
66152849c42SDave May 	fprintf(fp,"<DataField>\n");
66252849c42SDave May 	fprintf(fp,"%d\n", field->L);
66352849c42SDave May 	fprintf(fp,"%zu\n",field->atomic_size);
66452849c42SDave May 	fprintf(fp,"%s\n", field->registeration_function);
66552849c42SDave May 	fprintf(fp,"%s\n", field->name);
66652849c42SDave May 
66752849c42SDave May 	fwrite(field->data, field->atomic_size, field->L, fp);
66852849c42SDave May   /*
66952849c42SDave May    printf("  ** wrote %zu bytes for DataField \"%s\" \n", field->atomic_size * field->L, field->name );
67052849c42SDave May    */
67152849c42SDave May 	fprintf(fp,"\n</DataField>\n");
6727fbf63aeSDave May   PetscFunctionReturn(0);
67352849c42SDave May }
67452849c42SDave May 
6757fbf63aeSDave May #undef __FUNCT__
6767fbf63aeSDave May #define __FUNCT__ "_DataBucketRegisterFieldFromFile"
6777fbf63aeSDave May PetscErrorCode _DataBucketRegisterFieldFromFile(FILE *fp,DataBucket db)
67852849c42SDave May {
67952849c42SDave May 	PetscBool val;
68052849c42SDave May 	DataField *field;
68152849c42SDave May 
68252849c42SDave May 	DataField gfield;
68352849c42SDave May 	char dummy[100];
68452849c42SDave May 	char registeration_function[5000];
68552849c42SDave May 	char field_name[5000];
6865c18a9d6SDave May 	PetscInt L;
68752849c42SDave May 	size_t atomic_size,strL;
688*dbe06d34SDave May 	PetscErrorCode ierr;
68952849c42SDave May 
69052849c42SDave May 	/* check we haven't finalised the registration of fields */
69152849c42SDave May 	/*
69252849c42SDave May    if(db->finalised==PETSC_TRUE) {
69352849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
69452849c42SDave May    ERROR();
69552849c42SDave May    }
69652849c42SDave May    */
69752849c42SDave May 
69852849c42SDave May 
69952849c42SDave May 	/* read file contents */
70052849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
70152849c42SDave May 
70252849c42SDave May 	fscanf( fp, "%d\n",&L); //printf("read(L): %d\n", L);
70352849c42SDave May 
70452849c42SDave May 	fscanf( fp, "%zu\n",&atomic_size); //printf("read(size): %zu\n",atomic_size);
70552849c42SDave May 
70652849c42SDave May 	fgets(registeration_function,4999,fp); //printf("read(reg func): %s", registeration_function );
70752849c42SDave May 	strL = strlen(registeration_function);
70852849c42SDave May 	if (strL > 1) {
70952849c42SDave May 		registeration_function[strL-1] = 0;
71052849c42SDave May 	}
71152849c42SDave May 
71252849c42SDave May 	fgets(field_name,4999,fp); //printf("read(name): %s", field_name );
71352849c42SDave May 	strL = strlen(field_name);
71452849c42SDave May 	if (strL > 1) {
71552849c42SDave May 		field_name[strL-1] = 0;
71652849c42SDave May 	}
71752849c42SDave May 
7184b46c5e1SDave May #ifdef DATA_BUCKET_LOG
719b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"  ** read L=%D; atomic_size=%zu; reg_func=\"%s\"; name=\"%s\" \n", L,atomic_size,registeration_function,field_name);
72052849c42SDave May #endif
72152849c42SDave May 
72252849c42SDave May 
72352849c42SDave May 	/* check for repeated name */
72452849c42SDave May 	StringInList( field_name, db->nfields, (const DataField*)db->field, &val );
725a233d522SDave May 	if (val == PETSC_TRUE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot add same field twice");
72652849c42SDave May 
72752849c42SDave May 	/* create new space for data */
72852849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
72952849c42SDave May 	db->field     = field;
73052849c42SDave May 
73152849c42SDave May 	/* add field */
732*dbe06d34SDave May 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield );CHKERRQ(ierr);
73352849c42SDave May 
73452849c42SDave May 	/* copy contents of file */
73552849c42SDave May 	fread(gfield->data, gfield->atomic_size, gfield->L, fp);
7364b46c5e1SDave May #ifdef DATA_BUCKET_LOG
737b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"  ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name );
73852849c42SDave May #endif
73952849c42SDave May 	/* finish reading meta data */
74052849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
74152849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
74252849c42SDave May 
74352849c42SDave May 	db->field[ db->nfields ] = gfield;
74452849c42SDave May 
74552849c42SDave May 	db->nfields++;
7467fbf63aeSDave May   PetscFunctionReturn(0);
74752849c42SDave May }
74852849c42SDave May 
7497fbf63aeSDave May #undef __FUNCT__
7507fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderWrite_v00"
7517fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderWrite_v00(FILE *fp)
75252849c42SDave May {
75352849c42SDave May 	fprintf(fp,"<DataBucketHeader>\n");
75452849c42SDave May 	fprintf(fp,"type=DataBucket\n");
75552849c42SDave May 	fprintf(fp,"format=ascii\n");
75652849c42SDave May 	fprintf(fp,"version=0.0\n");
75752849c42SDave May 	fprintf(fp,"options=\n");
75852849c42SDave May 	fprintf(fp,"</DataBucketHeader>\n");
7597fbf63aeSDave May   PetscFunctionReturn(0);
76052849c42SDave May }
7617fbf63aeSDave May 
7627fbf63aeSDave May #undef __FUNCT__
7637fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderRead_v00"
7647fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderRead_v00(FILE *fp)
76552849c42SDave May {
76652849c42SDave May 	char dummy[100];
76752849c42SDave May 	size_t strL;
76852849c42SDave May 
76952849c42SDave May 	// header open
77052849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
77152849c42SDave May 
77252849c42SDave May 	// type
77352849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
77452849c42SDave May 	strL = strlen(dummy);
77552849c42SDave May 	if (strL > 1) { dummy[strL-1] = 0; }
77652849c42SDave May 	if (strcmp(dummy,"type=DataBucket") != 0) {
777a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Data file doesn't contain a DataBucket type");
77852849c42SDave May 	}
77952849c42SDave May 
78052849c42SDave May 	// format
78152849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
78252849c42SDave May 
78352849c42SDave May 	// version
78452849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
78552849c42SDave May 	strL = strlen(dummy);
78652849c42SDave May 	if (strL > 1) { dummy[strL-1] = 0; }
78752849c42SDave May 	if (strcmp(dummy,"version=0.0") != 0) {
788a233d522SDave May 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"DataBucket file must be parsed with version=0.0 : You tried %s", dummy);
78952849c42SDave May 	}
79052849c42SDave May 
79152849c42SDave May 	// options
79252849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
79352849c42SDave May 	// header close
79452849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
7957fbf63aeSDave May   PetscFunctionReturn(0);
79652849c42SDave May }
79752849c42SDave May 
7987fbf63aeSDave May #undef __FUNCT__
7997fbf63aeSDave May #define __FUNCT__ "_DataBucketLoadFromFileBinary_SEQ"
8007fbf63aeSDave May PetscErrorCode _DataBucketLoadFromFileBinary_SEQ(const char filename[],DataBucket *_db)
80152849c42SDave May {
80252849c42SDave May 	DataBucket db;
80352849c42SDave May 	FILE *fp;
8045c18a9d6SDave May 	PetscInt L,buffer,f,nfields;
805*dbe06d34SDave May 	PetscErrorCode ierr;
80652849c42SDave May 
80752849c42SDave May 
8084b46c5e1SDave May #ifdef DATA_BUCKET_LOG
809b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
81052849c42SDave May #endif
81152849c42SDave May 
81252849c42SDave May 	/* open file */
81352849c42SDave May 	fp = fopen(filename,"rb");
81452849c42SDave May 	if (fp == NULL) {
815a233d522SDave May 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file with name %s", filename);
81652849c42SDave May 	}
81752849c42SDave May 
81852849c42SDave May 	/* read header */
819*dbe06d34SDave May 	ierr = _DataBucketViewAscii_HeaderRead_v00(fp);CHKERRQ(ierr);
82052849c42SDave May 
82152849c42SDave May 	fscanf(fp,"%d\n%d\n%d\n",&L,&buffer,&nfields);
82252849c42SDave May 
823*dbe06d34SDave May 	ierr = DataBucketCreate(&db);CHKERRQ(ierr);
82452849c42SDave May 
82552849c42SDave May 	for (f=0; f<nfields; f++) {
826*dbe06d34SDave May 		ierr = _DataBucketRegisterFieldFromFile(fp,db);CHKERRQ(ierr);
82752849c42SDave May 	}
82852849c42SDave May 	fclose(fp);
82952849c42SDave May 
830*dbe06d34SDave May 	ierr = DataBucketFinalize(db);CHKERRQ(ierr);
83152849c42SDave May 
83252849c42SDave May   /*
83352849c42SDave May    DataBucketSetSizes(db,L,buffer);
83452849c42SDave May    */
83552849c42SDave May 	db->L = L;
83652849c42SDave May 	db->buffer = buffer;
83752849c42SDave May 	db->allocated = L + buffer;
83852849c42SDave May 
83952849c42SDave May 	*_db = db;
8407fbf63aeSDave May   PetscFunctionReturn(0);
84152849c42SDave May }
84252849c42SDave May 
8437fbf63aeSDave May #undef __FUNCT__
8447fbf63aeSDave May #define __FUNCT__ "DataBucketLoadFromFile"
8457fbf63aeSDave May PetscErrorCode DataBucketLoadFromFile(MPI_Comm comm,const char filename[],DataBucketViewType type,DataBucket *db)
84652849c42SDave May {
8475c18a9d6SDave May 	PetscMPIInt nproc,rank;
848997fa542SDave May 	PetscErrorCode ierr;
84952849c42SDave May 
850997fa542SDave May 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
851997fa542SDave May 	ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
85252849c42SDave May 
8534b46c5e1SDave May #ifdef DATA_BUCKET_LOG
854b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
85552849c42SDave May #endif
85652849c42SDave May 	if (type == DATABUCKET_VIEW_STDOUT) {
85752849c42SDave May 
85852849c42SDave May 	} else if (type == DATABUCKET_VIEW_ASCII) {
859a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
86052849c42SDave May 	} else if (type == DATABUCKET_VIEW_BINARY) {
86152849c42SDave May 		if (nproc == 1) {
862*dbe06d34SDave May 			ierr = _DataBucketLoadFromFileBinary_SEQ(filename,db);CHKERRQ(ierr);
86352849c42SDave May 		} else {
86452849c42SDave May 			char *name;
86552849c42SDave May 
86652849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
867*dbe06d34SDave May 			ierr = _DataBucketLoadFromFileBinary_SEQ(name,db);CHKERRQ(ierr);
86852849c42SDave May 			free(name);
86952849c42SDave May 		}
87052849c42SDave May 	} else {
871a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer requested");
87252849c42SDave May 	}
8737fbf63aeSDave May   PetscFunctionReturn(0);
87452849c42SDave May }
87552849c42SDave May 
8767fbf63aeSDave May #undef __FUNCT__
8777fbf63aeSDave May #define __FUNCT__ "_DataBucketViewBinary"
8787fbf63aeSDave May PetscErrorCode _DataBucketViewBinary(DataBucket db,const char filename[])
87952849c42SDave May {
88052849c42SDave May 	FILE *fp = NULL;
8815c18a9d6SDave May 	PetscInt f;
882*dbe06d34SDave May 	PetscErrorCode ierr;
88352849c42SDave May 
88452849c42SDave May 	fp = fopen(filename,"wb");
885a233d522SDave May 	if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file with name %s", filename);
88652849c42SDave May 
88752849c42SDave May 	/* db header */
888*dbe06d34SDave May 	ierr =_DataBucketViewAscii_HeaderWrite_v00(fp);CHKERRQ(ierr);
88952849c42SDave May 
89052849c42SDave May 	/* meta-data */
89152849c42SDave May 	fprintf(fp,"%d\n%d\n%d\n", db->L,db->buffer,db->nfields);
89252849c42SDave May 
89352849c42SDave May 	for (f=0; f<db->nfields; f++) {
89452849c42SDave May     /* load datafields */
895*dbe06d34SDave May 		ierr = _DataFieldViewBinary(db->field[f],fp);CHKERRQ(ierr);
89652849c42SDave May 	}
89752849c42SDave May 
89852849c42SDave May 	fclose(fp);
8997fbf63aeSDave May   PetscFunctionReturn(0);
90052849c42SDave May }
90152849c42SDave May 
9027fbf63aeSDave May #undef __FUNCT__
9037fbf63aeSDave May #define __FUNCT__ "DataBucketView_SEQ"
9047fbf63aeSDave May PetscErrorCode DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type)
90552849c42SDave May {
906*dbe06d34SDave May   PetscErrorCode ierr;
907*dbe06d34SDave May 
90852849c42SDave May 	switch (type) {
90952849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
91052849c42SDave May 		{
9115c18a9d6SDave May 			PetscInt f;
91252849c42SDave May 			double memory_usage_total = 0.0;
91352849c42SDave May 
914b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"DataBucketView(SEQ): (\"%s\")\n",filename);
915b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  L                  = %D \n", db->L );
916b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  buffer             = %D \n", db->buffer );
917b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  allocated          = %D \n", db->allocated );
91852849c42SDave May 
919b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  nfields registered = %D \n", db->nfields );
92052849c42SDave May 			for (f=0; f<db->nfields; f++) {
92152849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
92252849c42SDave May 
923b3a122caSDave May 				PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) \n", f, db->field[f]->name, memory_usage_f  );
92452849c42SDave May 				memory_usage_total += memory_usage_f;
92552849c42SDave May 			}
926b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) \n", memory_usage_total );
92752849c42SDave May 		}
92852849c42SDave May 			break;
92952849c42SDave May 
93052849c42SDave May 		case DATABUCKET_VIEW_ASCII:
93152849c42SDave May 		{
932a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
93352849c42SDave May 		}
93452849c42SDave May 			break;
93552849c42SDave May 
93652849c42SDave May 		case DATABUCKET_VIEW_BINARY:
93752849c42SDave May 		{
938*dbe06d34SDave May 			ierr = _DataBucketViewBinary(db,filename);CHKERRQ(ierr);
93952849c42SDave May 		}
94052849c42SDave May 			break;
94152849c42SDave May 
94252849c42SDave May 		case DATABUCKET_VIEW_HDF5:
94352849c42SDave May 		{
944a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No HDF5 support");
94552849c42SDave May 		}
94652849c42SDave May 			break;
94752849c42SDave May 
94852849c42SDave May 		default:
949a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
95052849c42SDave May 			break;
95152849c42SDave May 	}
9527fbf63aeSDave May   PetscFunctionReturn(0);
95352849c42SDave May }
95452849c42SDave May 
9557fbf63aeSDave May #undef __FUNCT__
9567fbf63aeSDave May #define __FUNCT__ "DataBucketView_MPI"
9577fbf63aeSDave May PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
95852849c42SDave May {
959997fa542SDave May   PetscErrorCode ierr;
960997fa542SDave May 
96152849c42SDave May 	switch (type) {
96252849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
96352849c42SDave May 		{
9645c18a9d6SDave May 			PetscInt f;
9655c18a9d6SDave May 			PetscInt L,buffer,allocated;
96652849c42SDave May 			double memory_usage_total,memory_usage_total_local = 0.0;
9675c18a9d6SDave May 			PetscMPIInt rank;
96852849c42SDave May 
969997fa542SDave May 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
97052849c42SDave May 
971*dbe06d34SDave May 			ierr = DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated);CHKERRQ(ierr);
97252849c42SDave May 
97352849c42SDave May 			for (f=0; f<db->nfields; f++) {
97452849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
97552849c42SDave May 
97652849c42SDave May 				memory_usage_total_local += memory_usage_f;
97752849c42SDave May 			}
978*dbe06d34SDave May 			ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
97952849c42SDave May 
98052849c42SDave May 			if (rank == 0) {
9815c18a9d6SDave May 				PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename);
9825c18a9d6SDave May 				PetscPrintf(comm,"  L                  = %D \n", L );
9835c18a9d6SDave May 				PetscPrintf(comm,"  buffer (max)       = %D \n", buffer );
9845c18a9d6SDave May 				PetscPrintf(comm,"  allocated          = %D \n", allocated );
98552849c42SDave May 
9865c18a9d6SDave May 				PetscPrintf(comm,"  nfields registered = %D \n", db->nfields );
98752849c42SDave May 				for (f=0; f<db->nfields; f++) {
98852849c42SDave May 					double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
98952849c42SDave May 
990b3a122caSDave May 					PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) : rank0\n", f, db->field[f]->name, memory_usage_f  );
99152849c42SDave May 				}
99252849c42SDave May 
993b3a122caSDave May 				PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) : collective\n", memory_usage_total );
99452849c42SDave May 			}
99552849c42SDave May 
99652849c42SDave May 		}
99752849c42SDave May 			break;
99852849c42SDave May 
99952849c42SDave May 		case DATABUCKET_VIEW_ASCII:
100052849c42SDave May 		{
1001a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying data structure");
100252849c42SDave May 		}
100352849c42SDave May 			break;
100452849c42SDave May 
100552849c42SDave May 		case DATABUCKET_VIEW_BINARY:
100652849c42SDave May 		{
100752849c42SDave May 			char *name;
10085c18a9d6SDave May 			PetscMPIInt rank;
100952849c42SDave May 
101052849c42SDave May 			/* create correct extension */
1011997fa542SDave May 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
101252849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
101352849c42SDave May 
1014*dbe06d34SDave May 			ierr = _DataBucketViewBinary(db,name);CHKERRQ(ierr);
101552849c42SDave May 
101652849c42SDave May 			free(name);
101752849c42SDave May 		}
101852849c42SDave May 			break;
101952849c42SDave May 
102052849c42SDave May 		case DATABUCKET_VIEW_HDF5:
102152849c42SDave May 		{
1022a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5");
102352849c42SDave May 		}
102452849c42SDave May 			break;
102552849c42SDave May 
102652849c42SDave May 		default:
1027a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
102852849c42SDave May 			break;
102952849c42SDave May 	}
10307fbf63aeSDave May   PetscFunctionReturn(0);
103152849c42SDave May }
103252849c42SDave May 
10337fbf63aeSDave May #undef __FUNCT__
10347fbf63aeSDave May #define __FUNCT__ "DataBucketView"
10357fbf63aeSDave May PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
103652849c42SDave May {
10375c18a9d6SDave May 	PetscMPIInt nproc;
1038997fa542SDave May   PetscErrorCode ierr;
103952849c42SDave May 
1040997fa542SDave May 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
104152849c42SDave May 	if (nproc == 1) {
1042*dbe06d34SDave May 		ierr =DataBucketView_SEQ(db,filename,type);CHKERRQ(ierr);
104352849c42SDave May 	} else {
1044*dbe06d34SDave May 		ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
104552849c42SDave May 	}
10467fbf63aeSDave May   PetscFunctionReturn(0);
104752849c42SDave May }
104852849c42SDave May 
10497fbf63aeSDave May #undef __FUNCT__
10507fbf63aeSDave May #define __FUNCT__ "DataBucketDuplicateFields"
10517fbf63aeSDave May PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
105252849c42SDave May {
105352849c42SDave May 	DataBucket db2;
10545c18a9d6SDave May 	PetscInt f;
1055*dbe06d34SDave May 	PetscErrorCode ierr;
105652849c42SDave May 
1057*dbe06d34SDave May 	ierr = DataBucketCreate(&db2);CHKERRQ(ierr);
105852849c42SDave May 
105952849c42SDave May 	/* copy contents from dbA into db2 */
106052849c42SDave May 	for (f=0; f<dbA->nfields; f++) {
106152849c42SDave May 		DataField field;
106252849c42SDave May 		size_t    atomic_size;
106352849c42SDave May 		char      *name;
106452849c42SDave May 
106552849c42SDave May 		field = dbA->field[f];
106652849c42SDave May 
106752849c42SDave May 		atomic_size = field->atomic_size;
106852849c42SDave May 		name        = field->name;
106952849c42SDave May 
1070*dbe06d34SDave May 		ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
107152849c42SDave May 	}
1072*dbe06d34SDave May 	ierr = DataBucketFinalize(db2);CHKERRQ(ierr);
1073*dbe06d34SDave May 	ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
107452849c42SDave May 
107552849c42SDave May 	/* set pointer */
107652849c42SDave May 	*dbB = db2;
10777fbf63aeSDave May   PetscFunctionReturn(0);
107852849c42SDave May }
107952849c42SDave May 
108052849c42SDave May /*
108152849c42SDave May  Insert points from db2 into db1
108252849c42SDave May  db1 <<== db2
108352849c42SDave May  */
10847fbf63aeSDave May #undef __FUNCT__
10857fbf63aeSDave May #define __FUNCT__ "DataBucketInsertValues"
10867fbf63aeSDave May PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2)
108752849c42SDave May {
10885c18a9d6SDave May 	PetscInt n_mp_points1,n_mp_points2;
10895c18a9d6SDave May 	PetscInt n_mp_points1_new,p;
1090*dbe06d34SDave May 	PetscErrorCode ierr;
109152849c42SDave May 
1092*dbe06d34SDave May 	ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr);
1093*dbe06d34SDave May 	ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr);
109452849c42SDave May 
109552849c42SDave May 	n_mp_points1_new = n_mp_points1 + n_mp_points2;
1096*dbe06d34SDave May 	ierr = DataBucketSetSizes(db1,n_mp_points1_new,-1);CHKERRQ(ierr);
109752849c42SDave May 
109852849c42SDave May 	for (p=0; p<n_mp_points2; p++) {
109952849c42SDave May 		// db1 <<== db2 //
1100*dbe06d34SDave May 		ierr =DataBucketCopyPoint( db2,p, db1,(n_mp_points1 + p) );CHKERRQ(ierr);
110152849c42SDave May 	}
11027fbf63aeSDave May   PetscFunctionReturn(0);
110352849c42SDave May }
110452849c42SDave May 
110552849c42SDave May /* helpers for parallel send/recv */
11067fbf63aeSDave May #undef __FUNCT__
11077fbf63aeSDave May #define __FUNCT__ "DataBucketCreatePackedArray"
11087fbf63aeSDave May PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
110952849c42SDave May {
11105c18a9d6SDave May   PetscInt       f;
111152849c42SDave May   size_t    sizeof_marker_contents;
111252849c42SDave May   void      *buffer;
111352849c42SDave May 
111452849c42SDave May   sizeof_marker_contents = 0;
111552849c42SDave May   for (f=0; f<db->nfields; f++) {
111652849c42SDave May     DataField df = db->field[f];
111752849c42SDave May 
111852849c42SDave May     sizeof_marker_contents += df->atomic_size;
111952849c42SDave May   }
112052849c42SDave May 
112152849c42SDave May   buffer = malloc(sizeof_marker_contents);
112252849c42SDave May   memset(buffer,0,sizeof_marker_contents);
112352849c42SDave May 
112452849c42SDave May   if (bytes) { *bytes = sizeof_marker_contents; }
112552849c42SDave May   if (buf)   { *buf   = buffer; }
11267fbf63aeSDave May   PetscFunctionReturn(0);
112752849c42SDave May }
112852849c42SDave May 
11297fbf63aeSDave May #undef __FUNCT__
11307fbf63aeSDave May #define __FUNCT__ "DataBucketDestroyPackedArray"
11317fbf63aeSDave May PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf)
113252849c42SDave May {
113352849c42SDave May   if (buf) {
113452849c42SDave May     free(*buf);
113552849c42SDave May     *buf = NULL;
113652849c42SDave May   }
11377fbf63aeSDave May   PetscFunctionReturn(0);
113852849c42SDave May }
113952849c42SDave May 
11407fbf63aeSDave May #undef __FUNCT__
11417fbf63aeSDave May #define __FUNCT__ "DataBucketFillPackedArray"
11427fbf63aeSDave May PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
114352849c42SDave May {
11445c18a9d6SDave May   PetscInt    f;
114552849c42SDave May   void   *data,*data_p;
114652849c42SDave May   size_t asize,offset;
114752849c42SDave May 
114852849c42SDave May   offset = 0;
114952849c42SDave May   for (f=0; f<db->nfields; f++) {
115052849c42SDave May     DataField df = db->field[f];
115152849c42SDave May 
115252849c42SDave May     asize = df->atomic_size;
115352849c42SDave May 
115452849c42SDave May     data = (void*)( df->data );
115552849c42SDave May     data_p = (void*)( (char*)data + index*asize );
115652849c42SDave May 
115752849c42SDave May     memcpy( (void*)((char*)buf + offset),  data_p,  asize);
115852849c42SDave May     offset = offset + asize;
115952849c42SDave May   }
11607fbf63aeSDave May   PetscFunctionReturn(0);
116152849c42SDave May }
116252849c42SDave May 
11637fbf63aeSDave May #undef __FUNCT__
11647fbf63aeSDave May #define __FUNCT__ "DataBucketInsertPackedArray"
11657fbf63aeSDave May PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
116652849c42SDave May {
11675c18a9d6SDave May   PetscInt f;
116852849c42SDave May   void *data_p;
116952849c42SDave May   size_t offset;
1170*dbe06d34SDave May 	PetscErrorCode ierr;
117152849c42SDave May 
117252849c42SDave May   offset = 0;
117352849c42SDave May   for (f=0; f<db->nfields; f++) {
117452849c42SDave May     DataField df = db->field[f];
117552849c42SDave May 
117652849c42SDave May     data_p = (void*)( (char*)data + offset );
117752849c42SDave May 
1178*dbe06d34SDave May     ierr = DataFieldInsertPoint(df, idx, (void*)data_p );CHKERRQ(ierr);
117952849c42SDave May     offset = offset + df->atomic_size;
118052849c42SDave May   }
11817fbf63aeSDave May   PetscFunctionReturn(0);
118252849c42SDave May }
1183