xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision e78fc626795089c6e4d0f2fa45b895689380b653)
152849c42SDave May 
252849c42SDave May #include "data_bucket.h"
352849c42SDave May 
452849c42SDave May /* string helpers */
52635f519SDave May #undef __FUNCT__
62635f519SDave May #define __FUNCT__ "StringInList"
72eac95f8SDave May PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val)
852849c42SDave May {
95c18a9d6SDave May 	PetscInt i;
1052849c42SDave May 
1152849c42SDave May 	*val = PETSC_FALSE;
1252849c42SDave May 	for (i=0; i<N; i++) {
1352849c42SDave May 		if (strcmp( name, gfield[i]->name) == 0 ) {
1452849c42SDave May 			*val = PETSC_TRUE;
152eac95f8SDave May       PetscFunctionReturn(0);
1652849c42SDave May 		}
1752849c42SDave May 	}
182eac95f8SDave May   PetscFunctionReturn(0);
1952849c42SDave May }
2052849c42SDave May 
212635f519SDave May #undef __FUNCT__
222635f519SDave May #define __FUNCT__ "StringFindInList"
232eac95f8SDave May PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index)
2452849c42SDave May {
255c18a9d6SDave May 	PetscInt i;
2652849c42SDave May 
2752849c42SDave May 	*index = -1;
2852849c42SDave May 	for (i=0; i<N; i++) {
2952849c42SDave May 		if (strcmp( name, gfield[i]->name ) == 0) {
3052849c42SDave May 			*index = i;
312eac95f8SDave May       PetscFunctionReturn(0);
3252849c42SDave May 		}
3352849c42SDave May 	}
342eac95f8SDave May   PetscFunctionReturn(0);
3552849c42SDave May }
3652849c42SDave May 
372eac95f8SDave May #undef __FUNCT__
382eac95f8SDave May #define __FUNCT__ "DataFieldCreate"
392eac95f8SDave May PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF)
4052849c42SDave May {
4152849c42SDave May 	DataField df;
4252849c42SDave May 
4352849c42SDave May 	df = malloc( sizeof(struct _p_DataField) );
4452849c42SDave May 	memset( df, 0, sizeof(struct _p_DataField) );
4552849c42SDave May 
4652849c42SDave May 
4752849c42SDave May 	asprintf( &df->registeration_function, "%s", registeration_function );
4852849c42SDave May 	asprintf( &df->name, "%s", name );
4952849c42SDave May 	df->atomic_size = size;
5052849c42SDave May 	df->L = L;
5152c3ed93SDave May   df->bs = 1;
5252849c42SDave May 
5352849c42SDave May 	df->data = malloc( size * L ); /* allocate something so we don't have to reallocate */
5452849c42SDave May 	memset( df->data, 0, size * L );
5552849c42SDave May 
5652849c42SDave May 	*DF = df;
572eac95f8SDave May   PetscFunctionReturn(0);
5852849c42SDave May }
5952849c42SDave May 
602eac95f8SDave May #undef __FUNCT__
612eac95f8SDave May #define __FUNCT__ "DataFieldDestroy"
622eac95f8SDave May PetscErrorCode DataFieldDestroy(DataField *DF)
6352849c42SDave May {
6452849c42SDave May 	DataField df = *DF;
6552849c42SDave May 
6652849c42SDave May 	free(df->registeration_function);
6752849c42SDave May 	free(df->name);
6852849c42SDave May 	free(df->data);
6952849c42SDave May 	free(df);
7052849c42SDave May 
7152849c42SDave May 	*DF = NULL;
722eac95f8SDave May   PetscFunctionReturn(0);
7352849c42SDave May }
7452849c42SDave May 
7552849c42SDave May /* data bucket */
762eac95f8SDave May #undef __FUNCT__
772eac95f8SDave May #define __FUNCT__ "DataBucketCreate"
782eac95f8SDave May PetscErrorCode DataBucketCreate(DataBucket *DB)
7952849c42SDave May {
8052849c42SDave May 	DataBucket db;
8152849c42SDave May 
8252849c42SDave May 
8352849c42SDave May 	db = malloc( sizeof(struct _p_DataBucket) );
8452849c42SDave May 	memset( db, 0, sizeof(struct _p_DataBucket) );
8552849c42SDave May 
8652849c42SDave May 	db->finalised = PETSC_FALSE;
8752849c42SDave May 
8852849c42SDave May 	/* create empty spaces for fields */
893454631fSDave May 	db->L         = -1;
9052849c42SDave May 	db->buffer    = 1;
9152849c42SDave May 	db->allocated = 1;
9252849c42SDave May 
9352849c42SDave May 	db->nfields   = 0;
9452849c42SDave May 	db->field     = malloc(sizeof(DataField));
9552849c42SDave May 
9652849c42SDave May 	*DB = db;
972eac95f8SDave May   PetscFunctionReturn(0);
9852849c42SDave May }
9952849c42SDave May 
1002eac95f8SDave May #undef __FUNCT__
1012eac95f8SDave May #define __FUNCT__ "DataBucketDestroy"
1022eac95f8SDave May PetscErrorCode DataBucketDestroy(DataBucket *DB)
10352849c42SDave May {
10452849c42SDave May 	DataBucket db = *DB;
1055c18a9d6SDave May 	PetscInt f;
106dbe06d34SDave May 	PetscErrorCode ierr;
10752849c42SDave May 
10852849c42SDave May 	/* release fields */
10952849c42SDave May 	for (f=0; f<db->nfields; f++) {
110dbe06d34SDave May 		ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
11152849c42SDave May 	}
11252849c42SDave May 
11352849c42SDave May 	/* this will catch the initially allocated objects in the event that no fields are registered */
11452849c42SDave May 	if (db->field != NULL) {
11552849c42SDave May 		free(db->field);
11652849c42SDave May 	}
11752849c42SDave May 
11852849c42SDave May 	free(db);
11952849c42SDave May 
12052849c42SDave May 	*DB = NULL;
1212eac95f8SDave May   PetscFunctionReturn(0);
12252849c42SDave May }
12352849c42SDave May 
1242eac95f8SDave May #undef __FUNCT__
1252eac95f8SDave May #define __FUNCT__ "DataBucketRegisterField"
1262eac95f8SDave May PetscErrorCode DataBucketRegisterField(
12752849c42SDave May                               DataBucket db,
12852849c42SDave May                               const char registeration_function[],
12952849c42SDave May                               const char field_name[],
13052849c42SDave May                               size_t atomic_size, DataField *_gfield)
13152849c42SDave May {
13252849c42SDave May 	PetscBool val;
13352849c42SDave May 	DataField *field,fp;
134dbe06d34SDave May 	PetscErrorCode ierr;
13552849c42SDave May 
13652849c42SDave May 	/* check we haven't finalised the registration of fields */
13752849c42SDave May 	/*
13852849c42SDave May    if(db->finalised==PETSC_TRUE) {
13952849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
14052849c42SDave May    ERROR();
14152849c42SDave May    }
14252849c42SDave May    */
14352849c42SDave May 
14452849c42SDave May 	/* check for repeated name */
1452635f519SDave May 	ierr = StringInList( field_name, db->nfields, (const DataField*)db->field, &val );CHKERRQ(ierr);
1462eac95f8SDave May 	if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
14752849c42SDave May 
14852849c42SDave May 	/* create new space for data */
14952849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
15052849c42SDave May 	db->field     = field;
15152849c42SDave May 
15252849c42SDave May 	/* add field */
153dbe06d34SDave May 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, db->allocated, &fp );CHKERRQ(ierr);
15452849c42SDave May 	db->field[ db->nfields ] = fp;
15552849c42SDave May 
15652849c42SDave May 	db->nfields++;
15752849c42SDave May 
15852849c42SDave May 	if (_gfield != NULL) {
15952849c42SDave May 		*_gfield = fp;
16052849c42SDave May 	}
1612eac95f8SDave May   PetscFunctionReturn(0);
16252849c42SDave May }
16352849c42SDave May 
16452849c42SDave May /*
16552849c42SDave May  #define DataBucketRegisterField(db,name,size,k) {\
16652849c42SDave May  char *location;\
16752849c42SDave May  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
16852849c42SDave May  _DataBucketRegisterField( (db), location, (name), (size), (k) );\
16952849c42SDave May  free(location);\
17052849c42SDave May  }
17152849c42SDave May  */
17252849c42SDave May 
1732eac95f8SDave May #undef __FUNCT__
1742eac95f8SDave May #define __FUNCT__ "DataBucketGetDataFieldByName"
1752eac95f8SDave May PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
17652849c42SDave May {
1772635f519SDave May   PetscErrorCode ierr;
1785c18a9d6SDave May 	PetscInt idx;
17952849c42SDave May 	PetscBool found;
18052849c42SDave May 
1812635f519SDave May 	ierr = StringInList(name,db->nfields,(const DataField*)db->field,&found);CHKERRQ(ierr);
1822eac95f8SDave May 	if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name);
1832eac95f8SDave May 
1842635f519SDave May 	ierr = StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);CHKERRQ(ierr);
18552849c42SDave May 
18652849c42SDave May 	*gfield = db->field[idx];
1872eac95f8SDave May   PetscFunctionReturn(0);
18852849c42SDave May }
18952849c42SDave May 
1907fbf63aeSDave May #undef __FUNCT__
1917fbf63aeSDave May #define __FUNCT__ "DataBucketQueryDataFieldByName"
1927fbf63aeSDave May PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
19352849c42SDave May {
1942635f519SDave May   PetscErrorCode ierr;
19552849c42SDave May 	*found = PETSC_FALSE;
1962635f519SDave May 	ierr = StringInList(name,db->nfields,(const DataField*)db->field,found);CHKERRQ(ierr);
1977fbf63aeSDave May   PetscFunctionReturn(0);
19852849c42SDave May }
19952849c42SDave May 
2007fbf63aeSDave May #undef __FUNCT__
2017fbf63aeSDave May #define __FUNCT__ "DataBucketFinalize"
2027fbf63aeSDave May PetscErrorCode DataBucketFinalize(DataBucket db)
20352849c42SDave May {
20452849c42SDave May 	db->finalised = PETSC_TRUE;
2057fbf63aeSDave May   PetscFunctionReturn(0);
20652849c42SDave May }
20752849c42SDave May 
2087fbf63aeSDave May #undef __FUNCT__
2097fbf63aeSDave May #define __FUNCT__ "DataFieldGetNumEntries"
2107fbf63aeSDave May PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum)
21152849c42SDave May {
21252849c42SDave May 	*sum = df->L;
2137fbf63aeSDave May   PetscFunctionReturn(0);
21452849c42SDave May }
21552849c42SDave May 
2167fbf63aeSDave May #undef __FUNCT__
21752c3ed93SDave May #define __FUNCT__ "DataFieldSetBlockSize"
21852c3ed93SDave May PetscErrorCode DataFieldSetBlockSize(DataField df,PetscInt blocksize)
21952c3ed93SDave May {
22052c3ed93SDave May 	df->bs = blocksize;
22152c3ed93SDave May   PetscFunctionReturn(0);
22252c3ed93SDave May }
22352c3ed93SDave May 
22452c3ed93SDave May #undef __FUNCT__
2257fbf63aeSDave May #define __FUNCT__ "DataFieldSetSize"
2267fbf63aeSDave May PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L)
22752849c42SDave May {
22852849c42SDave May 	void *tmp_data;
22952849c42SDave May 
230a233d522SDave May 	if (new_L <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be <= 0");
2317fbf63aeSDave May 
2327fbf63aeSDave May 	if (new_L == df->L) PetscFunctionReturn(0);
23352849c42SDave May 
23452849c42SDave May 	if (new_L > df->L) {
23552849c42SDave May 
23652849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L) );
23752849c42SDave May 		df->data = tmp_data;
23852849c42SDave May 
23952849c42SDave May 		/* init new contents */
24052849c42SDave May 		memset( ( ((char*)df->data)+df->L*df->atomic_size), 0, (new_L-df->L)*df->atomic_size );
24152849c42SDave May 
2427fbf63aeSDave May 	} else {
24352849c42SDave May 		/* reallocate pointer list, add +1 in case new_L = 0 */
24452849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L+1) );
24552849c42SDave May 		df->data = tmp_data;
24652849c42SDave May 	}
24752849c42SDave May 
24852849c42SDave May 	df->L = new_L;
2497fbf63aeSDave May   PetscFunctionReturn(0);
25052849c42SDave May }
25152849c42SDave May 
2527fbf63aeSDave May #undef __FUNCT__
2537fbf63aeSDave May #define __FUNCT__ "DataFieldZeroBlock"
2547fbf63aeSDave May PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end)
25552849c42SDave May {
2567fbf63aeSDave May 	if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
2577fbf63aeSDave May 
2587fbf63aeSDave May 	if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
2597fbf63aeSDave May 
260a233d522SDave 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);
26152849c42SDave May 
26252849c42SDave May 	memset( ( ((char*)df->data)+start*df->atomic_size), 0, (end-start)*df->atomic_size );
2637fbf63aeSDave May   PetscFunctionReturn(0);
26452849c42SDave May }
26552849c42SDave May 
26652849c42SDave May /*
26752849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
26852849c42SDave May  */
2697fbf63aeSDave May #undef __FUNCT__
2707fbf63aeSDave May #define __FUNCT__ "DataBucketSetSizes"
2717fbf63aeSDave May PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
27252849c42SDave May {
2735c18a9d6SDave May 	PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
274dbe06d34SDave May   PetscErrorCode ierr;
27552849c42SDave May 
2767fbf63aeSDave May 	if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()");
27752849c42SDave May 
27852849c42SDave May 	current_allocated = db->allocated;
27952849c42SDave May 
28052849c42SDave May 	new_used   = L;
28152849c42SDave May 	new_unused = current_allocated - new_used;
28252849c42SDave May 	new_buffer = db->buffer;
28352849c42SDave May 	if (buffer >= 0) { /* update the buffer value */
28452849c42SDave May 		new_buffer = buffer;
28552849c42SDave May 	}
28652849c42SDave May 	new_allocated = new_used + new_buffer;
28752849c42SDave May 
28852849c42SDave May 	/* action */
28952849c42SDave May 	if (new_allocated > current_allocated) {
29052849c42SDave May 		/* increase size to new_used + new_buffer */
29152849c42SDave May 		for (f=0; f<db->nfields; f++) {
292dbe06d34SDave May 			ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
29352849c42SDave May 		}
29452849c42SDave May 
29552849c42SDave May 		db->L         = new_used;
29652849c42SDave May 		db->buffer    = new_buffer;
29752849c42SDave May 		db->allocated = new_used + new_buffer;
29852849c42SDave May 	}
29952849c42SDave May 	else {
30052849c42SDave May 		if (new_unused > 2 * new_buffer) {
30152849c42SDave May 
30252849c42SDave May 			/* shrink array to new_used + new_buffer */
30352849c42SDave May 			for (f=0; f<db->nfields; f++) {
304dbe06d34SDave May 				ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
30552849c42SDave May 			}
30652849c42SDave May 
30752849c42SDave May 			db->L         = new_used;
30852849c42SDave May 			db->buffer    = new_buffer;
30952849c42SDave May 			db->allocated = new_used + new_buffer;
31052849c42SDave May 		}
31152849c42SDave May 		else {
31252849c42SDave May 			db->L      = new_used;
31352849c42SDave May 			db->buffer = new_buffer;
31452849c42SDave May 		}
31552849c42SDave May 	}
31652849c42SDave May 
31752849c42SDave May 	/* zero all entries from db->L to db->allocated */
31852849c42SDave May 	for (f=0; f<db->nfields; f++) {
31952849c42SDave May 		DataField field = db->field[f];
320dbe06d34SDave May 		ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
32152849c42SDave May 	}
3227fbf63aeSDave May   PetscFunctionReturn(0);
32352849c42SDave May }
32452849c42SDave May 
3257fbf63aeSDave May #undef __FUNCT__
3267fbf63aeSDave May #define __FUNCT__ "DataBucketSetInitialSizes"
3277fbf63aeSDave May PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
32852849c42SDave May {
3295c18a9d6SDave May 	PetscInt f;
330dbe06d34SDave May   PetscErrorCode ierr;
331dbe06d34SDave May 
332dbe06d34SDave May 	ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
33352849c42SDave May 
33452849c42SDave May 	for (f=0; f<db->nfields; f++) {
33552849c42SDave May 		DataField field = db->field[f];
336dbe06d34SDave May 		ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
33752849c42SDave May 	}
3387fbf63aeSDave May   PetscFunctionReturn(0);
33952849c42SDave May }
34052849c42SDave May 
3417fbf63aeSDave May #undef __FUNCT__
3427fbf63aeSDave May #define __FUNCT__ "DataBucketGetSizes"
3437fbf63aeSDave May PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
34452849c42SDave May {
34552849c42SDave May 	if (L) { *L = db->L; }
34652849c42SDave May 	if (buffer) { *buffer = db->buffer; }
34752849c42SDave May 	if (allocated) { *allocated = db->allocated; }
3487fbf63aeSDave May   PetscFunctionReturn(0);
34952849c42SDave May }
35052849c42SDave May 
3517fbf63aeSDave May #undef __FUNCT__
3527fbf63aeSDave May #define __FUNCT__ "DataBucketGetGlobalSizes"
3537fbf63aeSDave May PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
35452849c42SDave May {
3555c18a9d6SDave May 	PetscInt _L,_buffer,_allocated;
3565c18a9d6SDave May 	PetscInt ierr;
35752849c42SDave May 
3585c18a9d6SDave May 	_L = db->L;
3595c18a9d6SDave May 	_buffer = db->buffer;
3605c18a9d6SDave May 	_allocated = db->allocated;
36152849c42SDave May 
36233564166SDave May 	if (L) {         ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
36333564166SDave May 	if (buffer) {    ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
36433564166SDave May 	if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
3657fbf63aeSDave May   PetscFunctionReturn(0);
36652849c42SDave May }
36752849c42SDave May 
3687fbf63aeSDave May #undef __FUNCT__
36989884300SDave May #define __FUNCT__ "DataBucketGetDataFields"
3707fbf63aeSDave May PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[])
37152849c42SDave May {
37252849c42SDave May 	if (L) {      *L      = db->nfields; }
37352849c42SDave May 	if (fields) { *fields = db->field; }
3747fbf63aeSDave May   PetscFunctionReturn(0);
37552849c42SDave May }
37652849c42SDave May 
3777fbf63aeSDave May #undef __FUNCT__
3787fbf63aeSDave May #define __FUNCT__ "DataFieldGetAccess"
3797fbf63aeSDave May PetscErrorCode DataFieldGetAccess(const DataField gfield)
38052849c42SDave May {
3817fbf63aeSDave May 	if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name);
3827fbf63aeSDave May 
38352849c42SDave May 	gfield->active = PETSC_TRUE;
3847fbf63aeSDave May   PetscFunctionReturn(0);
38552849c42SDave May }
38652849c42SDave May 
3877fbf63aeSDave May #undef __FUNCT__
3887fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPoint"
3897fbf63aeSDave May PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p)
39052849c42SDave May {
39152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
39252849c42SDave May 	/* debug mode */
39384bcda08SDave May 	/* check point is valid */
3947fbf63aeSDave May 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
3957fbf63aeSDave May 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
39652849c42SDave May 
39784bcda08SDave 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);
39852849c42SDave May #endif
39952849c42SDave May 
40052849c42SDave May 	//*ctx_p  = (void*)( ((char*)gfield->data) + pid * gfield->atomic_size);
40152849c42SDave May 	*ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
4027fbf63aeSDave May   PetscFunctionReturn(0);
40352849c42SDave May }
40452849c42SDave May 
4057fbf63aeSDave May #undef __FUNCT__
4067fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPointOffset"
4077fbf63aeSDave May PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
40852849c42SDave May {
40952849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
41052849c42SDave May 	/* debug mode */
41152849c42SDave May 
41284bcda08SDave May 	/* check point is valid */
4137fbf63aeSDave 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 */
4147fbf63aeSDave May 	if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
41552849c42SDave May 
41684bcda08SDave May 	/* check point is valid */
4177fbf63aeSDave May 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
418a233d522SDave May 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
41952849c42SDave May 
420a233d522SDave 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);
42152849c42SDave May #endif
42252849c42SDave May 
42352849c42SDave May 	*ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
4247fbf63aeSDave May   PetscFunctionReturn(0);
42552849c42SDave May }
42652849c42SDave May 
4277fbf63aeSDave May #undef __FUNCT__
42889884300SDave May #define __FUNCT__ "DataFieldRestoreAccess"
4297fbf63aeSDave May PetscErrorCode DataFieldRestoreAccess(DataField gfield)
43052849c42SDave May {
4317fbf63aeSDave May 	if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name );
4327fbf63aeSDave May 
43352849c42SDave May 	gfield->active = PETSC_FALSE;
4347fbf63aeSDave May   PetscFunctionReturn(0);
43552849c42SDave May }
43652849c42SDave May 
4377fbf63aeSDave May #undef __FUNCT__
4387fbf63aeSDave May #define __FUNCT__ "DataFieldVerifyAccess"
4397fbf63aeSDave May PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size)
44052849c42SDave May {
44152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
4427fbf63aeSDave 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.",
44352849c42SDave May            gfield->name, gfield->atomic_size, size );
44452849c42SDave May #endif
4457fbf63aeSDave May   PetscFunctionReturn(0);
44652849c42SDave May }
44752849c42SDave May 
4487fbf63aeSDave May #undef __FUNCT__
4497fbf63aeSDave May #define __FUNCT__ "DataFieldGetAtomicSize"
4507fbf63aeSDave May PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size)
45152849c42SDave May {
45252849c42SDave May   if (size) { *size = gfield->atomic_size; }
4537fbf63aeSDave May   PetscFunctionReturn(0);
45452849c42SDave May }
45552849c42SDave May 
4567fbf63aeSDave May #undef __FUNCT__
4577fbf63aeSDave May #define __FUNCT__ "DataFieldGetEntries"
4587fbf63aeSDave May PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data)
45952849c42SDave May {
46052849c42SDave May   if (data) {
46152849c42SDave May     *data = gfield->data;
46252849c42SDave May   }
4637fbf63aeSDave May   PetscFunctionReturn(0);
46452849c42SDave May }
46552849c42SDave May 
4667fbf63aeSDave May #undef __FUNCT__
4677fbf63aeSDave May #define __FUNCT__ "DataFieldRestoreEntries"
4687fbf63aeSDave May PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data)
46952849c42SDave May {
47052849c42SDave May   if (data) {
47152849c42SDave May     *data = NULL;
47252849c42SDave May   }
4737fbf63aeSDave May   PetscFunctionReturn(0);
47452849c42SDave May }
47552849c42SDave May 
47652849c42SDave May /* y = x */
4777fbf63aeSDave May #undef __FUNCT__
4787fbf63aeSDave May #define __FUNCT__ "DataBucketCopyPoint"
4797fbf63aeSDave May PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,
4805c18a9d6SDave May                          const DataBucket yb,const PetscInt pid_y)
48152849c42SDave May {
4825c18a9d6SDave May 	PetscInt f;
483dbe06d34SDave May 	PetscErrorCode ierr;
484dbe06d34SDave May 
48552849c42SDave May 	for (f=0; f<xb->nfields; f++) {
48652849c42SDave May 		void *dest;
48752849c42SDave May 		void *src;
48852849c42SDave May 
489dbe06d34SDave May 		ierr = DataFieldGetAccess( xb->field[f] );CHKERRQ(ierr);
4906845f8f5SDave May 		if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f] );CHKERRQ(ierr); }
49152849c42SDave May 
492dbe06d34SDave May 		ierr = DataFieldAccessPoint( xb->field[f],pid_x, &src );CHKERRQ(ierr);
493dbe06d34SDave May 		ierr = DataFieldAccessPoint( yb->field[f],pid_y, &dest );CHKERRQ(ierr);
49452849c42SDave May 
49552849c42SDave May 		memcpy( dest, src, xb->field[f]->atomic_size );
49652849c42SDave May 
497dbe06d34SDave May 		ierr = DataFieldRestoreAccess( xb->field[f] );CHKERRQ(ierr);
4986845f8f5SDave May 		if (xb != yb) { ierr = DataFieldRestoreAccess( yb->field[f] );CHKERRQ(ierr); }
49952849c42SDave May 	}
5007fbf63aeSDave May   PetscFunctionReturn(0);
50152849c42SDave May }
50252849c42SDave May 
5037fbf63aeSDave May #undef __FUNCT__
5047fbf63aeSDave May #define __FUNCT__ "DataBucketCreateFromSubset"
5057fbf63aeSDave May PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB)
50652849c42SDave May {
5075c18a9d6SDave May 	PetscInt nfields;
50852849c42SDave May 	DataField *fields;
50952849c42SDave May 	DataBucketCreate(DB);
5105c18a9d6SDave May 	PetscInt f,L,buffer,allocated,p;
511dbe06d34SDave May 	PetscErrorCode ierr;
51252849c42SDave May 
51352849c42SDave May 	/* copy contents of DBIn */
514dbe06d34SDave May 	ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
515dbe06d34SDave May 	ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
51652849c42SDave May 
51752849c42SDave May 	for (f=0; f<nfields; f++) {
518dbe06d34SDave May 		ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
51952849c42SDave May 	}
520dbe06d34SDave May 	ierr = DataBucketFinalize(*DB);CHKERRQ(ierr);
52152849c42SDave May 
522dbe06d34SDave May 	ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
52352849c42SDave May 
52452849c42SDave May 	/* now copy the desired guys from DBIn => DB */
52552849c42SDave May 	for (p=0; p<N; p++) {
526dbe06d34SDave May 		ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
52752849c42SDave May 	}
5287fbf63aeSDave May   PetscFunctionReturn(0);
52952849c42SDave May }
53052849c42SDave May 
53152849c42SDave May // insert into an exisitng location
5327fbf63aeSDave May #undef __FUNCT__
5337fbf63aeSDave May #define __FUNCT__ "DataFieldInsertPoint"
5347fbf63aeSDave May PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx)
53552849c42SDave May {
53652849c42SDave May 
53752849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
53884bcda08SDave May 	/* check point is valid */
539a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
540a233d522SDave May 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
54152849c42SDave May #endif
54252849c42SDave May 
54352849c42SDave May   //	memcpy( (void*)((char*)field->data + index*field->atomic_size), ctx, field->atomic_size );
54452849c42SDave May 	memcpy( __DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size );
5457fbf63aeSDave May   PetscFunctionReturn(0);
54652849c42SDave May }
54752849c42SDave May 
54852849c42SDave May // remove data at index - replace with last point
549a233d522SDave May #undef __FUNCT__
550a233d522SDave May #define __FUNCT__ "DataBucketRemovePointAtIndex"
5517fbf63aeSDave May PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index)
55252849c42SDave May {
5535c18a9d6SDave May 	PetscInt f;
554dbe06d34SDave May 	PetscErrorCode ierr;
55552849c42SDave May 
55652849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
55784bcda08SDave May 	/* check point is valid */
558a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
559a233d522SDave May 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
56052849c42SDave May #endif
56152849c42SDave May 
562a233d522SDave May 	if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
563a233d522SDave 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 );
56452849c42SDave May 	}
56552849c42SDave May 
566a233d522SDave May 	if (index != db->L-1) { /* not last point in list */
56752849c42SDave May 		for (f=0; f<db->nfields; f++) {
56852849c42SDave May 			DataField field = db->field[f];
56952849c42SDave May 
57052849c42SDave May 			/* copy then remove */
571dbe06d34SDave May 			ierr = DataFieldCopyPoint( db->L-1,field, index,field );CHKERRQ(ierr);
57252849c42SDave May 
57352849c42SDave May 			//DataFieldZeroPoint(field,index);
57452849c42SDave May 		}
57552849c42SDave May 	}
57652849c42SDave May 
57752849c42SDave May 	/* decrement size */
57852849c42SDave May 	/* this will zero out an crap at the end of the list */
579dbe06d34SDave May 	ierr = DataBucketRemovePoint(db);CHKERRQ(ierr);
5807fbf63aeSDave May   PetscFunctionReturn(0);
58152849c42SDave May }
58252849c42SDave May 
58352849c42SDave May /* copy x into y */
5847fbf63aeSDave May #undef __FUNCT__
5857fbf63aeSDave May #define __FUNCT__ "DataFieldCopyPoint"
5867fbf63aeSDave May PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,
5875c18a9d6SDave May                         const PetscInt pid_y,const DataField field_y )
58852849c42SDave May {
58952849c42SDave May 
59052849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
59184bcda08SDave May 	/* check point is valid */
592a233d522SDave May 	if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
593a233d522SDave May 	if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
59452849c42SDave May 
595a233d522SDave May 	if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
596a233d522SDave May 	if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
59752849c42SDave May 
59852849c42SDave May 	if( field_y->atomic_size != field_x->atomic_size ) {
599a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
60052849c42SDave May 	}
60152849c42SDave May #endif
60252849c42SDave May 	/*
60352849c42SDave May    memcpy( (void*)((char*)field_y->data + pid_y*field_y->atomic_size),
60452849c42SDave May    (void*)((char*)field_x->data + pid_x*field_x->atomic_size),
60552849c42SDave May    field_x->atomic_size );
60652849c42SDave May    */
60752849c42SDave May 	memcpy(		__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
60852849c42SDave May          __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
60952849c42SDave May          field_y->atomic_size );
6107fbf63aeSDave May   PetscFunctionReturn(0);
61152849c42SDave May }
61252849c42SDave May 
61352849c42SDave May 
61452849c42SDave May // zero only the datafield at this point
6157fbf63aeSDave May #undef __FUNCT__
6167fbf63aeSDave May #define __FUNCT__ "DataFieldZeroPoint"
6177fbf63aeSDave May PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index)
61852849c42SDave May {
61952849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
62084bcda08SDave May 	/* check point is valid */
621a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
622a233d522SDave May 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
62352849c42SDave May #endif
62452849c42SDave May 
62552849c42SDave May   //	memset( (void*)((char*)field->data + index*field->atomic_size), 0, field->atomic_size );
62652849c42SDave May 	memset( __DATATFIELD_point_access(field->data,index,field->atomic_size), 0, field->atomic_size );
6277fbf63aeSDave May   PetscFunctionReturn(0);
62852849c42SDave May }
62952849c42SDave May 
63052849c42SDave May // zero ALL data for this point
6317fbf63aeSDave May #undef __FUNCT__
6327fbf63aeSDave May #define __FUNCT__ "DataBucketZeroPoint"
6337fbf63aeSDave May PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index)
63452849c42SDave May {
6355c18a9d6SDave May 	PetscInt f;
636dbe06d34SDave May 	PetscErrorCode ierr;
63752849c42SDave May 
63884bcda08SDave May 	/* check point is valid */
639a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
640a233d522SDave May 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
64152849c42SDave May 
64252849c42SDave May 	for (f=0; f<db->nfields; f++) {
64352849c42SDave May 		DataField field = db->field[f];
64452849c42SDave May 
645dbe06d34SDave May 		ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr);
64652849c42SDave May 	}
6477fbf63aeSDave May   PetscFunctionReturn(0);
64852849c42SDave May }
64952849c42SDave May 
65052849c42SDave May /* increment */
6517fbf63aeSDave May #undef __FUNCT__
6527fbf63aeSDave May #define __FUNCT__ "DataBucketAddPoint"
6537fbf63aeSDave May PetscErrorCode DataBucketAddPoint(DataBucket db)
65452849c42SDave May {
655dbe06d34SDave May 	PetscErrorCode ierr;
656dbe06d34SDave May 
657dbe06d34SDave May 	ierr = DataBucketSetSizes(db,db->L+1,-1);CHKERRQ(ierr);
6587fbf63aeSDave May   PetscFunctionReturn(0);
65952849c42SDave May }
66052849c42SDave May 
6617fbf63aeSDave May /* decrement */
6627fbf63aeSDave May #undef __FUNCT__
6637fbf63aeSDave May #define __FUNCT__ "DataBucketRemovePoint"
6647fbf63aeSDave May PetscErrorCode DataBucketRemovePoint(DataBucket db)
6657fbf63aeSDave May {
666dbe06d34SDave May 	PetscErrorCode ierr;
667dbe06d34SDave May 
668dbe06d34SDave May 	ierr = DataBucketSetSizes(db,db->L-1,-1);CHKERRQ(ierr);
6697fbf63aeSDave May   PetscFunctionReturn(0);
6707fbf63aeSDave May }
6717fbf63aeSDave May 
6727fbf63aeSDave May #undef __FUNCT__
6737fbf63aeSDave May #define __FUNCT__ "_DataFieldViewBinary"
6747fbf63aeSDave May PetscErrorCode _DataFieldViewBinary(DataField field,FILE *fp)
67552849c42SDave May {
67652849c42SDave May 	fprintf(fp,"<DataField>\n");
67752849c42SDave May 	fprintf(fp,"%d\n", field->L);
67852849c42SDave May 	fprintf(fp,"%zu\n",field->atomic_size);
67952849c42SDave May 	fprintf(fp,"%s\n", field->registeration_function);
68052849c42SDave May 	fprintf(fp,"%s\n", field->name);
68152849c42SDave May 
68252849c42SDave May 	fwrite(field->data, field->atomic_size, field->L, fp);
68352849c42SDave May   /*
68452849c42SDave May    printf("  ** wrote %zu bytes for DataField \"%s\" \n", field->atomic_size * field->L, field->name );
68552849c42SDave May    */
68652849c42SDave May 	fprintf(fp,"\n</DataField>\n");
6877fbf63aeSDave May   PetscFunctionReturn(0);
68852849c42SDave May }
68952849c42SDave May 
6907fbf63aeSDave May #undef __FUNCT__
6917fbf63aeSDave May #define __FUNCT__ "_DataBucketRegisterFieldFromFile"
6927fbf63aeSDave May PetscErrorCode _DataBucketRegisterFieldFromFile(FILE *fp,DataBucket db)
69352849c42SDave May {
69452849c42SDave May 	PetscBool val;
69552849c42SDave May 	DataField *field;
69652849c42SDave May 
69752849c42SDave May 	DataField gfield;
69852849c42SDave May 	char dummy[100];
69952849c42SDave May 	char registeration_function[5000];
70052849c42SDave May 	char field_name[5000];
7015c18a9d6SDave May 	PetscInt L;
70252849c42SDave May 	size_t atomic_size,strL;
703dbe06d34SDave May 	PetscErrorCode ierr;
70452849c42SDave May 
70552849c42SDave May 	/* check we haven't finalised the registration of fields */
70652849c42SDave May 	/*
70752849c42SDave May    if(db->finalised==PETSC_TRUE) {
70852849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
70952849c42SDave May    ERROR();
71052849c42SDave May    }
71152849c42SDave May    */
71252849c42SDave May 
71352849c42SDave May 
71452849c42SDave May 	/* read file contents */
71552849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
71652849c42SDave May 
71752849c42SDave May 	fscanf( fp, "%d\n",&L); //printf("read(L): %d\n", L);
71852849c42SDave May 
71952849c42SDave May 	fscanf( fp, "%zu\n",&atomic_size); //printf("read(size): %zu\n",atomic_size);
72052849c42SDave May 
72152849c42SDave May 	fgets(registeration_function,4999,fp); //printf("read(reg func): %s", registeration_function );
72252849c42SDave May 	strL = strlen(registeration_function);
72352849c42SDave May 	if (strL > 1) {
72452849c42SDave May 		registeration_function[strL-1] = 0;
72552849c42SDave May 	}
72652849c42SDave May 
72752849c42SDave May 	fgets(field_name,4999,fp); //printf("read(name): %s", field_name );
72852849c42SDave May 	strL = strlen(field_name);
72952849c42SDave May 	if (strL > 1) {
73052849c42SDave May 		field_name[strL-1] = 0;
73152849c42SDave May 	}
73252849c42SDave May 
7334b46c5e1SDave May #ifdef DATA_BUCKET_LOG
734b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"  ** read L=%D; atomic_size=%zu; reg_func=\"%s\"; name=\"%s\" \n", L,atomic_size,registeration_function,field_name);
73552849c42SDave May #endif
73652849c42SDave May 
73752849c42SDave May 
73852849c42SDave May 	/* check for repeated name */
7392635f519SDave May 	ierr = StringInList( field_name, db->nfields, (const DataField*)db->field, &val );CHKERRQ(ierr);
740a233d522SDave May 	if (val == PETSC_TRUE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot add same field twice");
74152849c42SDave May 
74252849c42SDave May 	/* create new space for data */
74352849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
74452849c42SDave May 	db->field     = field;
74552849c42SDave May 
74652849c42SDave May 	/* add field */
747dbe06d34SDave May 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield );CHKERRQ(ierr);
74852849c42SDave May 
74952849c42SDave May 	/* copy contents of file */
75052849c42SDave May 	fread(gfield->data, gfield->atomic_size, gfield->L, fp);
7514b46c5e1SDave May #ifdef DATA_BUCKET_LOG
752b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"  ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name );
75352849c42SDave May #endif
75452849c42SDave May 	/* finish reading meta data */
75552849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
75652849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
75752849c42SDave May 
75852849c42SDave May 	db->field[ db->nfields ] = gfield;
75952849c42SDave May 
76052849c42SDave May 	db->nfields++;
7617fbf63aeSDave May   PetscFunctionReturn(0);
76252849c42SDave May }
76352849c42SDave May 
7647fbf63aeSDave May #undef __FUNCT__
7657fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderWrite_v00"
7667fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderWrite_v00(FILE *fp)
76752849c42SDave May {
76852849c42SDave May 	fprintf(fp,"<DataBucketHeader>\n");
76952849c42SDave May 	fprintf(fp,"type=DataBucket\n");
77052849c42SDave May 	fprintf(fp,"format=ascii\n");
77152849c42SDave May 	fprintf(fp,"version=0.0\n");
77252849c42SDave May 	fprintf(fp,"options=\n");
77352849c42SDave May 	fprintf(fp,"</DataBucketHeader>\n");
7747fbf63aeSDave May   PetscFunctionReturn(0);
77552849c42SDave May }
7767fbf63aeSDave May 
7777fbf63aeSDave May #undef __FUNCT__
7787fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderRead_v00"
7797fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderRead_v00(FILE *fp)
78052849c42SDave May {
78152849c42SDave May 	char dummy[100];
78252849c42SDave May 	size_t strL;
78352849c42SDave May 
78452849c42SDave May 	// header open
78552849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
78652849c42SDave May 
78752849c42SDave May 	// type
78852849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
78952849c42SDave May 	strL = strlen(dummy);
79052849c42SDave May 	if (strL > 1) { dummy[strL-1] = 0; }
79152849c42SDave May 	if (strcmp(dummy,"type=DataBucket") != 0) {
792a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Data file doesn't contain a DataBucket type");
79352849c42SDave May 	}
79452849c42SDave May 
79552849c42SDave May 	// format
79652849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
79752849c42SDave May 
79852849c42SDave May 	// version
79952849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
80052849c42SDave May 	strL = strlen(dummy);
80152849c42SDave May 	if (strL > 1) { dummy[strL-1] = 0; }
80252849c42SDave May 	if (strcmp(dummy,"version=0.0") != 0) {
803a233d522SDave May 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"DataBucket file must be parsed with version=0.0 : You tried %s", dummy);
80452849c42SDave May 	}
80552849c42SDave May 
80652849c42SDave May 	// options
80752849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
80852849c42SDave May 	// header close
80952849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
8107fbf63aeSDave May   PetscFunctionReturn(0);
81152849c42SDave May }
81252849c42SDave May 
8137fbf63aeSDave May #undef __FUNCT__
8147fbf63aeSDave May #define __FUNCT__ "_DataBucketLoadFromFileBinary_SEQ"
8157fbf63aeSDave May PetscErrorCode _DataBucketLoadFromFileBinary_SEQ(const char filename[],DataBucket *_db)
81652849c42SDave May {
81752849c42SDave May 	DataBucket db;
81852849c42SDave May 	FILE *fp;
8195c18a9d6SDave May 	PetscInt L,buffer,f,nfields;
820dbe06d34SDave May 	PetscErrorCode ierr;
82152849c42SDave May 
82252849c42SDave May 
8234b46c5e1SDave May #ifdef DATA_BUCKET_LOG
824b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
82552849c42SDave May #endif
82652849c42SDave May 
82752849c42SDave May 	/* open file */
82852849c42SDave May 	fp = fopen(filename,"rb");
82952849c42SDave May 	if (fp == NULL) {
830a233d522SDave May 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file with name %s", filename);
83152849c42SDave May 	}
83252849c42SDave May 
83352849c42SDave May 	/* read header */
834dbe06d34SDave May 	ierr = _DataBucketViewAscii_HeaderRead_v00(fp);CHKERRQ(ierr);
83552849c42SDave May 
83652849c42SDave May 	fscanf(fp,"%d\n%d\n%d\n",&L,&buffer,&nfields);
83752849c42SDave May 
838dbe06d34SDave May 	ierr = DataBucketCreate(&db);CHKERRQ(ierr);
83952849c42SDave May 
84052849c42SDave May 	for (f=0; f<nfields; f++) {
841dbe06d34SDave May 		ierr = _DataBucketRegisterFieldFromFile(fp,db);CHKERRQ(ierr);
84252849c42SDave May 	}
84352849c42SDave May 	fclose(fp);
84452849c42SDave May 
845dbe06d34SDave May 	ierr = DataBucketFinalize(db);CHKERRQ(ierr);
84652849c42SDave May 
84752849c42SDave May   /*
84852849c42SDave May    DataBucketSetSizes(db,L,buffer);
84952849c42SDave May    */
85052849c42SDave May 	db->L = L;
85152849c42SDave May 	db->buffer = buffer;
85252849c42SDave May 	db->allocated = L + buffer;
85352849c42SDave May 
85452849c42SDave May 	*_db = db;
8557fbf63aeSDave May   PetscFunctionReturn(0);
85652849c42SDave May }
85752849c42SDave May 
8587fbf63aeSDave May #undef __FUNCT__
8597fbf63aeSDave May #define __FUNCT__ "DataBucketLoadFromFile"
8607fbf63aeSDave May PetscErrorCode DataBucketLoadFromFile(MPI_Comm comm,const char filename[],DataBucketViewType type,DataBucket *db)
86152849c42SDave May {
8625c18a9d6SDave May 	PetscMPIInt nproc,rank;
863997fa542SDave May 	PetscErrorCode ierr;
86452849c42SDave May 
865997fa542SDave May 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
866997fa542SDave May 	ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
86752849c42SDave May 
8684b46c5e1SDave May #ifdef DATA_BUCKET_LOG
869b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
87052849c42SDave May #endif
87152849c42SDave May 	if (type == DATABUCKET_VIEW_STDOUT) {
87252849c42SDave May 
87352849c42SDave May 	} else if (type == DATABUCKET_VIEW_ASCII) {
874a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
87552849c42SDave May 	} else if (type == DATABUCKET_VIEW_BINARY) {
87652849c42SDave May 		if (nproc == 1) {
877dbe06d34SDave May 			ierr = _DataBucketLoadFromFileBinary_SEQ(filename,db);CHKERRQ(ierr);
87852849c42SDave May 		} else {
87952849c42SDave May 			char *name;
88052849c42SDave May 
88152849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
882dbe06d34SDave May 			ierr = _DataBucketLoadFromFileBinary_SEQ(name,db);CHKERRQ(ierr);
88352849c42SDave May 			free(name);
88452849c42SDave May 		}
88552849c42SDave May 	} else {
886a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer requested");
88752849c42SDave May 	}
8887fbf63aeSDave May   PetscFunctionReturn(0);
88952849c42SDave May }
89052849c42SDave May 
8917fbf63aeSDave May #undef __FUNCT__
8927fbf63aeSDave May #define __FUNCT__ "_DataBucketViewBinary"
8937fbf63aeSDave May PetscErrorCode _DataBucketViewBinary(DataBucket db,const char filename[])
89452849c42SDave May {
89552849c42SDave May 	FILE *fp = NULL;
8965c18a9d6SDave May 	PetscInt f;
897dbe06d34SDave May 	PetscErrorCode ierr;
89852849c42SDave May 
89952849c42SDave May 	fp = fopen(filename,"wb");
900a233d522SDave May 	if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file with name %s", filename);
90152849c42SDave May 
90252849c42SDave May 	/* db header */
903dbe06d34SDave May 	ierr =_DataBucketViewAscii_HeaderWrite_v00(fp);CHKERRQ(ierr);
90452849c42SDave May 
90552849c42SDave May 	/* meta-data */
90652849c42SDave May 	fprintf(fp,"%d\n%d\n%d\n", db->L,db->buffer,db->nfields);
90752849c42SDave May 
90852849c42SDave May 	for (f=0; f<db->nfields; f++) {
90952849c42SDave May     /* load datafields */
910dbe06d34SDave May 		ierr = _DataFieldViewBinary(db->field[f],fp);CHKERRQ(ierr);
91152849c42SDave May 	}
91252849c42SDave May 
91352849c42SDave May 	fclose(fp);
9147fbf63aeSDave May   PetscFunctionReturn(0);
91552849c42SDave May }
91652849c42SDave May 
9177fbf63aeSDave May #undef __FUNCT__
9187fbf63aeSDave May #define __FUNCT__ "DataBucketView_SEQ"
9197fbf63aeSDave May PetscErrorCode DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type)
92052849c42SDave May {
921dbe06d34SDave May   PetscErrorCode ierr;
922dbe06d34SDave May 
92352849c42SDave May 	switch (type) {
92452849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
92552849c42SDave May 		{
9265c18a9d6SDave May 			PetscInt f;
92752849c42SDave May 			double memory_usage_total = 0.0;
92852849c42SDave May 
929b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"DataBucketView(SEQ): (\"%s\")\n",filename);
930b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  L                  = %D \n", db->L );
931b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  buffer             = %D \n", db->buffer );
932b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  allocated          = %D \n", db->allocated );
93352849c42SDave May 
934b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  nfields registered = %D \n", db->nfields );
93552849c42SDave May 			for (f=0; f<db->nfields; f++) {
93652849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
93752849c42SDave May 
938b3a122caSDave May 				PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) \n", f, db->field[f]->name, memory_usage_f  );
93952c3ed93SDave May         PetscPrintf(PETSC_COMM_SELF,"           blocksize          = %D \n", db->field[f]->bs );
94052c3ed93SDave May         if (db->field[f]->bs != 1) {
941*e78fc626SDave May           PetscPrintf(PETSC_COMM_SELF,"           atomic size        = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs );
94252c3ed93SDave May           PetscPrintf(PETSC_COMM_SELF,"           atomic size/item   = %zu \n", db->field[f]->atomic_size/db->field[f]->bs );
94352c3ed93SDave May         } else {
94452c3ed93SDave May           PetscPrintf(PETSC_COMM_SELF,"           atomic size        = %zu \n", db->field[f]->atomic_size );
94552c3ed93SDave May         }
94652849c42SDave May 				memory_usage_total += memory_usage_f;
94752849c42SDave May 			}
948b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) \n", memory_usage_total );
94952849c42SDave May 		}
95052849c42SDave May 			break;
95152849c42SDave May 
95252849c42SDave May 		case DATABUCKET_VIEW_ASCII:
95352849c42SDave May 		{
954a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
95552849c42SDave May 		}
95652849c42SDave May 			break;
95752849c42SDave May 
95852849c42SDave May 		case DATABUCKET_VIEW_BINARY:
95952849c42SDave May 		{
960dbe06d34SDave May 			ierr = _DataBucketViewBinary(db,filename);CHKERRQ(ierr);
96152849c42SDave May 		}
96252849c42SDave May 			break;
96352849c42SDave May 
96452849c42SDave May 		case DATABUCKET_VIEW_HDF5:
96552849c42SDave May 		{
966a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No HDF5 support");
96752849c42SDave May 		}
96852849c42SDave May 			break;
96952849c42SDave May 
97052849c42SDave May 		default:
971a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
97252849c42SDave May 			break;
97352849c42SDave May 	}
9747fbf63aeSDave May   PetscFunctionReturn(0);
97552849c42SDave May }
97652849c42SDave May 
9777fbf63aeSDave May #undef __FUNCT__
9787fbf63aeSDave May #define __FUNCT__ "DataBucketView_MPI"
9797fbf63aeSDave May PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
98052849c42SDave May {
981997fa542SDave May   PetscErrorCode ierr;
982997fa542SDave May 
98352849c42SDave May 	switch (type) {
98452849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
98552849c42SDave May 		{
9865c18a9d6SDave May 			PetscInt f;
9875c18a9d6SDave May 			PetscInt L,buffer,allocated;
98852849c42SDave May 			double memory_usage_total,memory_usage_total_local = 0.0;
9895c18a9d6SDave May 			PetscMPIInt rank;
99052849c42SDave May 
991997fa542SDave May 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
99252849c42SDave May 
993dbe06d34SDave May 			ierr = DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated);CHKERRQ(ierr);
99452849c42SDave May 
99552849c42SDave May 			for (f=0; f<db->nfields; f++) {
99652849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
99752849c42SDave May 
99852849c42SDave May 				memory_usage_total_local += memory_usage_f;
99952849c42SDave May 			}
1000dbe06d34SDave May 			ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
100152849c42SDave May 
100252849c42SDave May 			if (rank == 0) {
10035c18a9d6SDave May 				PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename);
10045c18a9d6SDave May 				PetscPrintf(comm,"  L                  = %D \n", L );
10055c18a9d6SDave May 				PetscPrintf(comm,"  buffer (max)       = %D \n", buffer );
10065c18a9d6SDave May 				PetscPrintf(comm,"  allocated          = %D \n", allocated );
100752849c42SDave May 
10085c18a9d6SDave May 				PetscPrintf(comm,"  nfields registered = %D \n", db->nfields );
100952849c42SDave May 				for (f=0; f<db->nfields; f++) {
101052849c42SDave May 					double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
101152849c42SDave May 
1012b3a122caSDave May 					PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) : rank0\n", f, db->field[f]->name, memory_usage_f  );
101352849c42SDave May 				}
101452849c42SDave May 
1015b3a122caSDave May 				PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) : collective\n", memory_usage_total );
101652849c42SDave May 			}
101752849c42SDave May 
101852849c42SDave May 		}
101952849c42SDave May 			break;
102052849c42SDave May 
102152849c42SDave May 		case DATABUCKET_VIEW_ASCII:
102252849c42SDave May 		{
1023a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying data structure");
102452849c42SDave May 		}
102552849c42SDave May 			break;
102652849c42SDave May 
102752849c42SDave May 		case DATABUCKET_VIEW_BINARY:
102852849c42SDave May 		{
102952849c42SDave May 			char *name;
10305c18a9d6SDave May 			PetscMPIInt rank;
103152849c42SDave May 
103252849c42SDave May 			/* create correct extension */
1033997fa542SDave May 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
103452849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
103552849c42SDave May 
1036dbe06d34SDave May 			ierr = _DataBucketViewBinary(db,name);CHKERRQ(ierr);
103752849c42SDave May 
103852849c42SDave May 			free(name);
103952849c42SDave May 		}
104052849c42SDave May 			break;
104152849c42SDave May 
104252849c42SDave May 		case DATABUCKET_VIEW_HDF5:
104352849c42SDave May 		{
1044a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5");
104552849c42SDave May 		}
104652849c42SDave May 			break;
104752849c42SDave May 
104852849c42SDave May 		default:
1049a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
105052849c42SDave May 			break;
105152849c42SDave May 	}
10527fbf63aeSDave May   PetscFunctionReturn(0);
105352849c42SDave May }
105452849c42SDave May 
10557fbf63aeSDave May #undef __FUNCT__
10567fbf63aeSDave May #define __FUNCT__ "DataBucketView"
10577fbf63aeSDave May PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
105852849c42SDave May {
10595c18a9d6SDave May 	PetscMPIInt nproc;
1060997fa542SDave May   PetscErrorCode ierr;
106152849c42SDave May 
1062997fa542SDave May 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
106352849c42SDave May 	if (nproc == 1) {
1064dbe06d34SDave May 		ierr =DataBucketView_SEQ(db,filename,type);CHKERRQ(ierr);
106552849c42SDave May 	} else {
1066dbe06d34SDave May 		ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
106752849c42SDave May 	}
10687fbf63aeSDave May   PetscFunctionReturn(0);
106952849c42SDave May }
107052849c42SDave May 
10717fbf63aeSDave May #undef __FUNCT__
10727fbf63aeSDave May #define __FUNCT__ "DataBucketDuplicateFields"
10737fbf63aeSDave May PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
107452849c42SDave May {
107552849c42SDave May 	DataBucket db2;
10765c18a9d6SDave May 	PetscInt f;
1077dbe06d34SDave May 	PetscErrorCode ierr;
107852849c42SDave May 
1079dbe06d34SDave May 	ierr = DataBucketCreate(&db2);CHKERRQ(ierr);
108052849c42SDave May 
108152849c42SDave May 	/* copy contents from dbA into db2 */
108252849c42SDave May 	for (f=0; f<dbA->nfields; f++) {
108352849c42SDave May 		DataField field;
108452849c42SDave May 		size_t    atomic_size;
108552849c42SDave May 		char      *name;
108652849c42SDave May 
108752849c42SDave May 		field = dbA->field[f];
108852849c42SDave May 
108952849c42SDave May 		atomic_size = field->atomic_size;
109052849c42SDave May 		name        = field->name;
109152849c42SDave May 
1092dbe06d34SDave May 		ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
109352849c42SDave May 	}
1094dbe06d34SDave May 	ierr = DataBucketFinalize(db2);CHKERRQ(ierr);
1095dbe06d34SDave May 	ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
109652849c42SDave May 
109752849c42SDave May 	/* set pointer */
109852849c42SDave May 	*dbB = db2;
10997fbf63aeSDave May   PetscFunctionReturn(0);
110052849c42SDave May }
110152849c42SDave May 
110252849c42SDave May /*
110352849c42SDave May  Insert points from db2 into db1
110452849c42SDave May  db1 <<== db2
110552849c42SDave May  */
11067fbf63aeSDave May #undef __FUNCT__
11077fbf63aeSDave May #define __FUNCT__ "DataBucketInsertValues"
11087fbf63aeSDave May PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2)
110952849c42SDave May {
11105c18a9d6SDave May 	PetscInt n_mp_points1,n_mp_points2;
11115c18a9d6SDave May 	PetscInt n_mp_points1_new,p;
1112dbe06d34SDave May 	PetscErrorCode ierr;
111352849c42SDave May 
1114dbe06d34SDave May 	ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr);
1115dbe06d34SDave May 	ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr);
111652849c42SDave May 
111752849c42SDave May 	n_mp_points1_new = n_mp_points1 + n_mp_points2;
1118dbe06d34SDave May 	ierr = DataBucketSetSizes(db1,n_mp_points1_new,-1);CHKERRQ(ierr);
111952849c42SDave May 
112052849c42SDave May 	for (p=0; p<n_mp_points2; p++) {
112152849c42SDave May 		// db1 <<== db2 //
1122dbe06d34SDave May 		ierr =DataBucketCopyPoint( db2,p, db1,(n_mp_points1 + p) );CHKERRQ(ierr);
112352849c42SDave May 	}
11247fbf63aeSDave May   PetscFunctionReturn(0);
112552849c42SDave May }
112652849c42SDave May 
112752849c42SDave May /* helpers for parallel send/recv */
11287fbf63aeSDave May #undef __FUNCT__
11297fbf63aeSDave May #define __FUNCT__ "DataBucketCreatePackedArray"
11307fbf63aeSDave May PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
113152849c42SDave May {
11325c18a9d6SDave May   PetscInt       f;
113352849c42SDave May   size_t    sizeof_marker_contents;
113452849c42SDave May   void      *buffer;
113552849c42SDave May 
113652849c42SDave May   sizeof_marker_contents = 0;
113752849c42SDave May   for (f=0; f<db->nfields; f++) {
113852849c42SDave May     DataField df = db->field[f];
113952849c42SDave May 
114052849c42SDave May     sizeof_marker_contents += df->atomic_size;
114152849c42SDave May   }
114252849c42SDave May 
114352849c42SDave May   buffer = malloc(sizeof_marker_contents);
114452849c42SDave May   memset(buffer,0,sizeof_marker_contents);
114552849c42SDave May 
114652849c42SDave May   if (bytes) { *bytes = sizeof_marker_contents; }
114752849c42SDave May   if (buf)   { *buf   = buffer; }
11487fbf63aeSDave May   PetscFunctionReturn(0);
114952849c42SDave May }
115052849c42SDave May 
11517fbf63aeSDave May #undef __FUNCT__
11527fbf63aeSDave May #define __FUNCT__ "DataBucketDestroyPackedArray"
11537fbf63aeSDave May PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf)
115452849c42SDave May {
115552849c42SDave May   if (buf) {
115652849c42SDave May     free(*buf);
115752849c42SDave May     *buf = NULL;
115852849c42SDave May   }
11597fbf63aeSDave May   PetscFunctionReturn(0);
116052849c42SDave May }
116152849c42SDave May 
11627fbf63aeSDave May #undef __FUNCT__
11637fbf63aeSDave May #define __FUNCT__ "DataBucketFillPackedArray"
11647fbf63aeSDave May PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
116552849c42SDave May {
11665c18a9d6SDave May   PetscInt    f;
116752849c42SDave May   void   *data,*data_p;
116852849c42SDave May   size_t asize,offset;
116952849c42SDave May 
117052849c42SDave May   offset = 0;
117152849c42SDave May   for (f=0; f<db->nfields; f++) {
117252849c42SDave May     DataField df = db->field[f];
117352849c42SDave May 
117452849c42SDave May     asize = df->atomic_size;
117552849c42SDave May 
117652849c42SDave May     data = (void*)( df->data );
117752849c42SDave May     data_p = (void*)( (char*)data + index*asize );
117852849c42SDave May 
117952849c42SDave May     memcpy( (void*)((char*)buf + offset),  data_p,  asize);
118052849c42SDave May     offset = offset + asize;
118152849c42SDave May   }
11827fbf63aeSDave May   PetscFunctionReturn(0);
118352849c42SDave May }
118452849c42SDave May 
11857fbf63aeSDave May #undef __FUNCT__
11867fbf63aeSDave May #define __FUNCT__ "DataBucketInsertPackedArray"
11877fbf63aeSDave May PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
118852849c42SDave May {
11895c18a9d6SDave May   PetscInt f;
119052849c42SDave May   void *data_p;
119152849c42SDave May   size_t offset;
1192dbe06d34SDave May 	PetscErrorCode ierr;
119352849c42SDave May 
119452849c42SDave May   offset = 0;
119552849c42SDave May   for (f=0; f<db->nfields; f++) {
119652849c42SDave May     DataField df = db->field[f];
119752849c42SDave May 
119852849c42SDave May     data_p = (void*)( (char*)data + offset );
119952849c42SDave May 
1200dbe06d34SDave May     ierr = DataFieldInsertPoint(df, idx, (void*)data_p );CHKERRQ(ierr);
120152849c42SDave May     offset = offset + df->atomic_size;
120252849c42SDave May   }
12037fbf63aeSDave May   PetscFunctionReturn(0);
120452849c42SDave May }
1205