xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision 2635f51937deeb34105c00bbd4aa264ff3f215f6)
152849c42SDave May 
252849c42SDave May #include "data_bucket.h"
352849c42SDave May 
452849c42SDave May /* string helpers */
5*2635f519SDave May #undef __FUNCT__
6*2635f519SDave 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 
21*2635f519SDave May #undef __FUNCT__
22*2635f519SDave 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;
5152849c42SDave May 
5252849c42SDave May 	df->data = malloc( size * L ); /* allocate something so we don't have to reallocate */
5352849c42SDave May 	memset( df->data, 0, size * L );
5452849c42SDave May 
5552849c42SDave May 	*DF = df;
562eac95f8SDave May   PetscFunctionReturn(0);
5752849c42SDave May }
5852849c42SDave May 
592eac95f8SDave May #undef __FUNCT__
602eac95f8SDave May #define __FUNCT__ "DataFieldDestroy"
612eac95f8SDave May PetscErrorCode DataFieldDestroy(DataField *DF)
6252849c42SDave May {
6352849c42SDave May 	DataField df = *DF;
6452849c42SDave May 
6552849c42SDave May 	free(df->registeration_function);
6652849c42SDave May 	free(df->name);
6752849c42SDave May 	free(df->data);
6852849c42SDave May 	free(df);
6952849c42SDave May 
7052849c42SDave May 	*DF = NULL;
712eac95f8SDave May   PetscFunctionReturn(0);
7252849c42SDave May }
7352849c42SDave May 
7452849c42SDave May /* data bucket */
752eac95f8SDave May #undef __FUNCT__
762eac95f8SDave May #define __FUNCT__ "DataBucketCreate"
772eac95f8SDave May PetscErrorCode DataBucketCreate(DataBucket *DB)
7852849c42SDave May {
7952849c42SDave May 	DataBucket db;
8052849c42SDave May 
8152849c42SDave May 
8252849c42SDave May 	db = malloc( sizeof(struct _p_DataBucket) );
8352849c42SDave May 	memset( db, 0, sizeof(struct _p_DataBucket) );
8452849c42SDave May 
8552849c42SDave May 	db->finalised = PETSC_FALSE;
8652849c42SDave May 
8752849c42SDave May 	/* create empty spaces for fields */
883454631fSDave May 	db->L         = -1;
8952849c42SDave May 	db->buffer    = 1;
9052849c42SDave May 	db->allocated = 1;
9152849c42SDave May 
9252849c42SDave May 	db->nfields   = 0;
9352849c42SDave May 	db->field     = malloc(sizeof(DataField));
9452849c42SDave May 
9552849c42SDave May 	*DB = db;
962eac95f8SDave May   PetscFunctionReturn(0);
9752849c42SDave May }
9852849c42SDave May 
992eac95f8SDave May #undef __FUNCT__
1002eac95f8SDave May #define __FUNCT__ "DataBucketDestroy"
1012eac95f8SDave May PetscErrorCode DataBucketDestroy(DataBucket *DB)
10252849c42SDave May {
10352849c42SDave May 	DataBucket db = *DB;
1045c18a9d6SDave May 	PetscInt f;
105dbe06d34SDave May 	PetscErrorCode ierr;
10652849c42SDave May 
10752849c42SDave May 	/* release fields */
10852849c42SDave May 	for (f=0; f<db->nfields; f++) {
109dbe06d34SDave May 		ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
11052849c42SDave May 	}
11152849c42SDave May 
11252849c42SDave May 	/* this will catch the initially allocated objects in the event that no fields are registered */
11352849c42SDave May 	if (db->field != NULL) {
11452849c42SDave May 		free(db->field);
11552849c42SDave May 	}
11652849c42SDave May 
11752849c42SDave May 	free(db);
11852849c42SDave May 
11952849c42SDave May 	*DB = NULL;
1202eac95f8SDave May   PetscFunctionReturn(0);
12152849c42SDave May }
12252849c42SDave May 
1232eac95f8SDave May #undef __FUNCT__
1242eac95f8SDave May #define __FUNCT__ "DataBucketRegisterField"
1252eac95f8SDave May PetscErrorCode DataBucketRegisterField(
12652849c42SDave May                               DataBucket db,
12752849c42SDave May                               const char registeration_function[],
12852849c42SDave May                               const char field_name[],
12952849c42SDave May                               size_t atomic_size, DataField *_gfield)
13052849c42SDave May {
13152849c42SDave May 	PetscBool val;
13252849c42SDave May 	DataField *field,fp;
133dbe06d34SDave May 	PetscErrorCode ierr;
13452849c42SDave May 
13552849c42SDave May 	/* check we haven't finalised the registration of fields */
13652849c42SDave May 	/*
13752849c42SDave May    if(db->finalised==PETSC_TRUE) {
13852849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
13952849c42SDave May    ERROR();
14052849c42SDave May    }
14152849c42SDave May    */
14252849c42SDave May 
14352849c42SDave May 	/* check for repeated name */
144*2635f519SDave May 	ierr = StringInList( field_name, db->nfields, (const DataField*)db->field, &val );CHKERRQ(ierr);
1452eac95f8SDave May 	if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
14652849c42SDave May 
14752849c42SDave May 	/* create new space for data */
14852849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
14952849c42SDave May 	db->field     = field;
15052849c42SDave May 
15152849c42SDave May 	/* add field */
152dbe06d34SDave May 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, db->allocated, &fp );CHKERRQ(ierr);
15352849c42SDave May 	db->field[ db->nfields ] = fp;
15452849c42SDave May 
15552849c42SDave May 	db->nfields++;
15652849c42SDave May 
15752849c42SDave May 	if (_gfield != NULL) {
15852849c42SDave May 		*_gfield = fp;
15952849c42SDave May 	}
1602eac95f8SDave May   PetscFunctionReturn(0);
16152849c42SDave May }
16252849c42SDave May 
16352849c42SDave May /*
16452849c42SDave May  #define DataBucketRegisterField(db,name,size,k) {\
16552849c42SDave May  char *location;\
16652849c42SDave May  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
16752849c42SDave May  _DataBucketRegisterField( (db), location, (name), (size), (k) );\
16852849c42SDave May  free(location);\
16952849c42SDave May  }
17052849c42SDave May  */
17152849c42SDave May 
1722eac95f8SDave May #undef __FUNCT__
1732eac95f8SDave May #define __FUNCT__ "DataBucketGetDataFieldByName"
1742eac95f8SDave May PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
17552849c42SDave May {
176*2635f519SDave May   PetscErrorCode ierr;
1775c18a9d6SDave May 	PetscInt idx;
17852849c42SDave May 	PetscBool found;
17952849c42SDave May 
180*2635f519SDave May 	ierr = StringInList(name,db->nfields,(const DataField*)db->field,&found);CHKERRQ(ierr);
1812eac95f8SDave May 	if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name);
1822eac95f8SDave May 
183*2635f519SDave May 	ierr = StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);CHKERRQ(ierr);
18452849c42SDave May 
18552849c42SDave May 	*gfield = db->field[idx];
1862eac95f8SDave May   PetscFunctionReturn(0);
18752849c42SDave May }
18852849c42SDave May 
1897fbf63aeSDave May #undef __FUNCT__
1907fbf63aeSDave May #define __FUNCT__ "DataBucketQueryDataFieldByName"
1917fbf63aeSDave May PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
19252849c42SDave May {
193*2635f519SDave May   PetscErrorCode ierr;
19452849c42SDave May 	*found = PETSC_FALSE;
195*2635f519SDave May 	ierr = StringInList(name,db->nfields,(const DataField*)db->field,found);CHKERRQ(ierr);
1967fbf63aeSDave May   PetscFunctionReturn(0);
19752849c42SDave May }
19852849c42SDave May 
1997fbf63aeSDave May #undef __FUNCT__
2007fbf63aeSDave May #define __FUNCT__ "DataBucketFinalize"
2017fbf63aeSDave May PetscErrorCode DataBucketFinalize(DataBucket db)
20252849c42SDave May {
20352849c42SDave May 	db->finalised = PETSC_TRUE;
2047fbf63aeSDave May   PetscFunctionReturn(0);
20552849c42SDave May }
20652849c42SDave May 
2077fbf63aeSDave May #undef __FUNCT__
2087fbf63aeSDave May #define __FUNCT__ "DataFieldGetNumEntries"
2097fbf63aeSDave May PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum)
21052849c42SDave May {
21152849c42SDave May 	*sum = df->L;
2127fbf63aeSDave May   PetscFunctionReturn(0);
21352849c42SDave May }
21452849c42SDave May 
2157fbf63aeSDave May #undef __FUNCT__
2167fbf63aeSDave May #define __FUNCT__ "DataFieldSetSize"
2177fbf63aeSDave May PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L)
21852849c42SDave May {
21952849c42SDave May 	void *tmp_data;
22052849c42SDave May 
221a233d522SDave May 	if (new_L <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be <= 0");
2227fbf63aeSDave May 
2237fbf63aeSDave May 	if (new_L == df->L) PetscFunctionReturn(0);
22452849c42SDave May 
22552849c42SDave May 	if (new_L > df->L) {
22652849c42SDave May 
22752849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L) );
22852849c42SDave May 		df->data = tmp_data;
22952849c42SDave May 
23052849c42SDave May 		/* init new contents */
23152849c42SDave May 		memset( ( ((char*)df->data)+df->L*df->atomic_size), 0, (new_L-df->L)*df->atomic_size );
23252849c42SDave May 
2337fbf63aeSDave May 	} else {
23452849c42SDave May 		/* reallocate pointer list, add +1 in case new_L = 0 */
23552849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L+1) );
23652849c42SDave May 		df->data = tmp_data;
23752849c42SDave May 	}
23852849c42SDave May 
23952849c42SDave May 	df->L = new_L;
2407fbf63aeSDave May   PetscFunctionReturn(0);
24152849c42SDave May }
24252849c42SDave May 
2437fbf63aeSDave May #undef __FUNCT__
2447fbf63aeSDave May #define __FUNCT__ "DataFieldZeroBlock"
2457fbf63aeSDave May PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end)
24652849c42SDave May {
2477fbf63aeSDave May 	if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
2487fbf63aeSDave May 
2497fbf63aeSDave May 	if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
2507fbf63aeSDave May 
251a233d522SDave 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);
25252849c42SDave May 
25352849c42SDave May 	memset( ( ((char*)df->data)+start*df->atomic_size), 0, (end-start)*df->atomic_size );
2547fbf63aeSDave May   PetscFunctionReturn(0);
25552849c42SDave May }
25652849c42SDave May 
25752849c42SDave May /*
25852849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
25952849c42SDave May  */
2607fbf63aeSDave May #undef __FUNCT__
2617fbf63aeSDave May #define __FUNCT__ "DataBucketSetSizes"
2627fbf63aeSDave May PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
26352849c42SDave May {
2645c18a9d6SDave May 	PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
265dbe06d34SDave May   PetscErrorCode ierr;
26652849c42SDave May 
2677fbf63aeSDave May 	if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()");
26852849c42SDave May 
26952849c42SDave May 	current_allocated = db->allocated;
27052849c42SDave May 
27152849c42SDave May 	new_used   = L;
27252849c42SDave May 	new_unused = current_allocated - new_used;
27352849c42SDave May 	new_buffer = db->buffer;
27452849c42SDave May 	if (buffer >= 0) { /* update the buffer value */
27552849c42SDave May 		new_buffer = buffer;
27652849c42SDave May 	}
27752849c42SDave May 	new_allocated = new_used + new_buffer;
27852849c42SDave May 
27952849c42SDave May 	/* action */
28052849c42SDave May 	if (new_allocated > current_allocated) {
28152849c42SDave May 		/* increase size to new_used + new_buffer */
28252849c42SDave May 		for (f=0; f<db->nfields; f++) {
283dbe06d34SDave May 			ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
28452849c42SDave May 		}
28552849c42SDave May 
28652849c42SDave May 		db->L         = new_used;
28752849c42SDave May 		db->buffer    = new_buffer;
28852849c42SDave May 		db->allocated = new_used + new_buffer;
28952849c42SDave May 	}
29052849c42SDave May 	else {
29152849c42SDave May 		if (new_unused > 2 * new_buffer) {
29252849c42SDave May 
29352849c42SDave May 			/* shrink array to new_used + new_buffer */
29452849c42SDave May 			for (f=0; f<db->nfields; f++) {
295dbe06d34SDave May 				ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
29652849c42SDave May 			}
29752849c42SDave May 
29852849c42SDave May 			db->L         = new_used;
29952849c42SDave May 			db->buffer    = new_buffer;
30052849c42SDave May 			db->allocated = new_used + new_buffer;
30152849c42SDave May 		}
30252849c42SDave May 		else {
30352849c42SDave May 			db->L      = new_used;
30452849c42SDave May 			db->buffer = new_buffer;
30552849c42SDave May 		}
30652849c42SDave May 	}
30752849c42SDave May 
30852849c42SDave May 	/* zero all entries from db->L to db->allocated */
30952849c42SDave May 	for (f=0; f<db->nfields; f++) {
31052849c42SDave May 		DataField field = db->field[f];
311dbe06d34SDave May 		ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
31252849c42SDave May 	}
3137fbf63aeSDave May   PetscFunctionReturn(0);
31452849c42SDave May }
31552849c42SDave May 
3167fbf63aeSDave May #undef __FUNCT__
3177fbf63aeSDave May #define __FUNCT__ "DataBucketSetInitialSizes"
3187fbf63aeSDave May PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
31952849c42SDave May {
3205c18a9d6SDave May 	PetscInt f;
321dbe06d34SDave May   PetscErrorCode ierr;
322dbe06d34SDave May 
323dbe06d34SDave May 	ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
32452849c42SDave May 
32552849c42SDave May 	for (f=0; f<db->nfields; f++) {
32652849c42SDave May 		DataField field = db->field[f];
327dbe06d34SDave May 		ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
32852849c42SDave May 	}
3297fbf63aeSDave May   PetscFunctionReturn(0);
33052849c42SDave May }
33152849c42SDave May 
3327fbf63aeSDave May #undef __FUNCT__
3337fbf63aeSDave May #define __FUNCT__ "DataBucketGetSizes"
3347fbf63aeSDave May PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
33552849c42SDave May {
33652849c42SDave May 	if (L) { *L = db->L; }
33752849c42SDave May 	if (buffer) { *buffer = db->buffer; }
33852849c42SDave May 	if (allocated) { *allocated = db->allocated; }
3397fbf63aeSDave May   PetscFunctionReturn(0);
34052849c42SDave May }
34152849c42SDave May 
3427fbf63aeSDave May #undef __FUNCT__
3437fbf63aeSDave May #define __FUNCT__ "DataBucketGetGlobalSizes"
3447fbf63aeSDave May PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
34552849c42SDave May {
3465c18a9d6SDave May 	PetscInt _L,_buffer,_allocated;
3475c18a9d6SDave May 	PetscInt ierr;
34852849c42SDave May 
3495c18a9d6SDave May 	_L = db->L;
3505c18a9d6SDave May 	_buffer = db->buffer;
3515c18a9d6SDave May 	_allocated = db->allocated;
35252849c42SDave May 
35333564166SDave May 	if (L) {         ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
35433564166SDave May 	if (buffer) {    ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
35533564166SDave May 	if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
3567fbf63aeSDave May   PetscFunctionReturn(0);
35752849c42SDave May }
35852849c42SDave May 
3597fbf63aeSDave May #undef __FUNCT__
36089884300SDave May #define __FUNCT__ "DataBucketGetDataFields"
3617fbf63aeSDave May PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[])
36252849c42SDave May {
36352849c42SDave May 	if (L) {      *L      = db->nfields; }
36452849c42SDave May 	if (fields) { *fields = db->field; }
3657fbf63aeSDave May   PetscFunctionReturn(0);
36652849c42SDave May }
36752849c42SDave May 
3687fbf63aeSDave May #undef __FUNCT__
3697fbf63aeSDave May #define __FUNCT__ "DataFieldGetAccess"
3707fbf63aeSDave May PetscErrorCode DataFieldGetAccess(const DataField gfield)
37152849c42SDave May {
3727fbf63aeSDave May 	if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name);
3737fbf63aeSDave May 
37452849c42SDave May 	gfield->active = PETSC_TRUE;
3757fbf63aeSDave May   PetscFunctionReturn(0);
37652849c42SDave May }
37752849c42SDave May 
3787fbf63aeSDave May #undef __FUNCT__
3797fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPoint"
3807fbf63aeSDave May PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p)
38152849c42SDave May {
38252849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
38352849c42SDave May 	/* debug mode */
38484bcda08SDave May 	/* check point is valid */
3857fbf63aeSDave May 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
3867fbf63aeSDave May 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
38752849c42SDave May 
38884bcda08SDave 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);
38952849c42SDave May #endif
39052849c42SDave May 
39152849c42SDave May 	//*ctx_p  = (void*)( ((char*)gfield->data) + pid * gfield->atomic_size);
39252849c42SDave May 	*ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
3937fbf63aeSDave May   PetscFunctionReturn(0);
39452849c42SDave May }
39552849c42SDave May 
3967fbf63aeSDave May #undef __FUNCT__
3977fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPointOffset"
3987fbf63aeSDave May PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
39952849c42SDave May {
40052849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
40152849c42SDave May 	/* debug mode */
40252849c42SDave May 
40384bcda08SDave May 	/* check point is valid */
4047fbf63aeSDave 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 */
4057fbf63aeSDave May 	if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
40652849c42SDave May 
40784bcda08SDave May 	/* check point is valid */
4087fbf63aeSDave May 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
409a233d522SDave May 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
41052849c42SDave May 
411a233d522SDave 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);
41252849c42SDave May #endif
41352849c42SDave May 
41452849c42SDave May 	*ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
4157fbf63aeSDave May   PetscFunctionReturn(0);
41652849c42SDave May }
41752849c42SDave May 
4187fbf63aeSDave May #undef __FUNCT__
41989884300SDave May #define __FUNCT__ "DataFieldRestoreAccess"
4207fbf63aeSDave May PetscErrorCode DataFieldRestoreAccess(DataField gfield)
42152849c42SDave May {
4227fbf63aeSDave May 	if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name );
4237fbf63aeSDave May 
42452849c42SDave May 	gfield->active = PETSC_FALSE;
4257fbf63aeSDave May   PetscFunctionReturn(0);
42652849c42SDave May }
42752849c42SDave May 
4287fbf63aeSDave May #undef __FUNCT__
4297fbf63aeSDave May #define __FUNCT__ "DataFieldVerifyAccess"
4307fbf63aeSDave May PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size)
43152849c42SDave May {
43252849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
4337fbf63aeSDave 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.",
43452849c42SDave May            gfield->name, gfield->atomic_size, size );
43552849c42SDave May #endif
4367fbf63aeSDave May   PetscFunctionReturn(0);
43752849c42SDave May }
43852849c42SDave May 
4397fbf63aeSDave May #undef __FUNCT__
4407fbf63aeSDave May #define __FUNCT__ "DataFieldGetAtomicSize"
4417fbf63aeSDave May PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size)
44252849c42SDave May {
44352849c42SDave May   if (size) { *size = gfield->atomic_size; }
4447fbf63aeSDave May   PetscFunctionReturn(0);
44552849c42SDave May }
44652849c42SDave May 
4477fbf63aeSDave May #undef __FUNCT__
4487fbf63aeSDave May #define __FUNCT__ "DataFieldGetEntries"
4497fbf63aeSDave May PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data)
45052849c42SDave May {
45152849c42SDave May   if (data) {
45252849c42SDave May     *data = gfield->data;
45352849c42SDave May   }
4547fbf63aeSDave May   PetscFunctionReturn(0);
45552849c42SDave May }
45652849c42SDave May 
4577fbf63aeSDave May #undef __FUNCT__
4587fbf63aeSDave May #define __FUNCT__ "DataFieldRestoreEntries"
4597fbf63aeSDave May PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data)
46052849c42SDave May {
46152849c42SDave May   if (data) {
46252849c42SDave May     *data = NULL;
46352849c42SDave May   }
4647fbf63aeSDave May   PetscFunctionReturn(0);
46552849c42SDave May }
46652849c42SDave May 
46752849c42SDave May /* y = x */
4687fbf63aeSDave May #undef __FUNCT__
4697fbf63aeSDave May #define __FUNCT__ "DataBucketCopyPoint"
4707fbf63aeSDave May PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,
4715c18a9d6SDave May                          const DataBucket yb,const PetscInt pid_y)
47252849c42SDave May {
4735c18a9d6SDave May 	PetscInt f;
474dbe06d34SDave May 	PetscErrorCode ierr;
475dbe06d34SDave May 
47652849c42SDave May 	for (f=0; f<xb->nfields; f++) {
47752849c42SDave May 		void *dest;
47852849c42SDave May 		void *src;
47952849c42SDave May 
480dbe06d34SDave May 		ierr = DataFieldGetAccess( xb->field[f] );CHKERRQ(ierr);
4816845f8f5SDave May 		if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f] );CHKERRQ(ierr); }
48252849c42SDave May 
483dbe06d34SDave May 		ierr = DataFieldAccessPoint( xb->field[f],pid_x, &src );CHKERRQ(ierr);
484dbe06d34SDave May 		ierr = DataFieldAccessPoint( yb->field[f],pid_y, &dest );CHKERRQ(ierr);
48552849c42SDave May 
48652849c42SDave May 		memcpy( dest, src, xb->field[f]->atomic_size );
48752849c42SDave May 
488dbe06d34SDave May 		ierr = DataFieldRestoreAccess( xb->field[f] );CHKERRQ(ierr);
4896845f8f5SDave May 		if (xb != yb) { ierr = DataFieldRestoreAccess( yb->field[f] );CHKERRQ(ierr); }
49052849c42SDave May 	}
4917fbf63aeSDave May   PetscFunctionReturn(0);
49252849c42SDave May }
49352849c42SDave May 
4947fbf63aeSDave May #undef __FUNCT__
4957fbf63aeSDave May #define __FUNCT__ "DataBucketCreateFromSubset"
4967fbf63aeSDave May PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB)
49752849c42SDave May {
4985c18a9d6SDave May 	PetscInt nfields;
49952849c42SDave May 	DataField *fields;
50052849c42SDave May 	DataBucketCreate(DB);
5015c18a9d6SDave May 	PetscInt f,L,buffer,allocated,p;
502dbe06d34SDave May 	PetscErrorCode ierr;
50352849c42SDave May 
50452849c42SDave May 	/* copy contents of DBIn */
505dbe06d34SDave May 	ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
506dbe06d34SDave May 	ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
50752849c42SDave May 
50852849c42SDave May 	for (f=0; f<nfields; f++) {
509dbe06d34SDave May 		ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
51052849c42SDave May 	}
511dbe06d34SDave May 	ierr = DataBucketFinalize(*DB);CHKERRQ(ierr);
51252849c42SDave May 
513dbe06d34SDave May 	ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
51452849c42SDave May 
51552849c42SDave May 	/* now copy the desired guys from DBIn => DB */
51652849c42SDave May 	for (p=0; p<N; p++) {
517dbe06d34SDave May 		ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
51852849c42SDave May 	}
5197fbf63aeSDave May   PetscFunctionReturn(0);
52052849c42SDave May }
52152849c42SDave May 
52252849c42SDave May // insert into an exisitng location
5237fbf63aeSDave May #undef __FUNCT__
5247fbf63aeSDave May #define __FUNCT__ "DataFieldInsertPoint"
5257fbf63aeSDave May PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx)
52652849c42SDave May {
52752849c42SDave May 
52852849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
52984bcda08SDave May 	/* check point is valid */
530a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
531a233d522SDave May 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
53252849c42SDave May #endif
53352849c42SDave May 
53452849c42SDave May   //	memcpy( (void*)((char*)field->data + index*field->atomic_size), ctx, field->atomic_size );
53552849c42SDave May 	memcpy( __DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size );
5367fbf63aeSDave May   PetscFunctionReturn(0);
53752849c42SDave May }
53852849c42SDave May 
53952849c42SDave May // remove data at index - replace with last point
540a233d522SDave May #undef __FUNCT__
541a233d522SDave May #define __FUNCT__ "DataBucketRemovePointAtIndex"
5427fbf63aeSDave May PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index)
54352849c42SDave May {
5445c18a9d6SDave May 	PetscInt f;
545dbe06d34SDave May 	PetscErrorCode ierr;
54652849c42SDave May 
54752849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
54884bcda08SDave May 	/* check point is valid */
549a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
550a233d522SDave May 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
55152849c42SDave May #endif
55252849c42SDave May 
553a233d522SDave May 	if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
554a233d522SDave 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 );
55552849c42SDave May 	}
55652849c42SDave May 
557a233d522SDave May 	if (index != db->L-1) { /* not last point in list */
55852849c42SDave May 		for (f=0; f<db->nfields; f++) {
55952849c42SDave May 			DataField field = db->field[f];
56052849c42SDave May 
56152849c42SDave May 			/* copy then remove */
562dbe06d34SDave May 			ierr = DataFieldCopyPoint( db->L-1,field, index,field );CHKERRQ(ierr);
56352849c42SDave May 
56452849c42SDave May 			//DataFieldZeroPoint(field,index);
56552849c42SDave May 		}
56652849c42SDave May 	}
56752849c42SDave May 
56852849c42SDave May 	/* decrement size */
56952849c42SDave May 	/* this will zero out an crap at the end of the list */
570dbe06d34SDave May 	ierr = DataBucketRemovePoint(db);CHKERRQ(ierr);
5717fbf63aeSDave May   PetscFunctionReturn(0);
57252849c42SDave May }
57352849c42SDave May 
57452849c42SDave May /* copy x into y */
5757fbf63aeSDave May #undef __FUNCT__
5767fbf63aeSDave May #define __FUNCT__ "DataFieldCopyPoint"
5777fbf63aeSDave May PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,
5785c18a9d6SDave May                         const PetscInt pid_y,const DataField field_y )
57952849c42SDave May {
58052849c42SDave May 
58152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
58284bcda08SDave May 	/* check point is valid */
583a233d522SDave May 	if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
584a233d522SDave May 	if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
58552849c42SDave May 
586a233d522SDave May 	if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
587a233d522SDave May 	if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
58852849c42SDave May 
58952849c42SDave May 	if( field_y->atomic_size != field_x->atomic_size ) {
590a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
59152849c42SDave May 	}
59252849c42SDave May #endif
59352849c42SDave May 	/*
59452849c42SDave May    memcpy( (void*)((char*)field_y->data + pid_y*field_y->atomic_size),
59552849c42SDave May    (void*)((char*)field_x->data + pid_x*field_x->atomic_size),
59652849c42SDave May    field_x->atomic_size );
59752849c42SDave May    */
59852849c42SDave May 	memcpy(		__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
59952849c42SDave May          __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
60052849c42SDave May          field_y->atomic_size );
6017fbf63aeSDave May   PetscFunctionReturn(0);
60252849c42SDave May }
60352849c42SDave May 
60452849c42SDave May 
60552849c42SDave May // zero only the datafield at this point
6067fbf63aeSDave May #undef __FUNCT__
6077fbf63aeSDave May #define __FUNCT__ "DataFieldZeroPoint"
6087fbf63aeSDave May PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index)
60952849c42SDave May {
61052849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
61184bcda08SDave May 	/* check point is valid */
612a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
613a233d522SDave May 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
61452849c42SDave May #endif
61552849c42SDave May 
61652849c42SDave May   //	memset( (void*)((char*)field->data + index*field->atomic_size), 0, field->atomic_size );
61752849c42SDave May 	memset( __DATATFIELD_point_access(field->data,index,field->atomic_size), 0, field->atomic_size );
6187fbf63aeSDave May   PetscFunctionReturn(0);
61952849c42SDave May }
62052849c42SDave May 
62152849c42SDave May // zero ALL data for this point
6227fbf63aeSDave May #undef __FUNCT__
6237fbf63aeSDave May #define __FUNCT__ "DataBucketZeroPoint"
6247fbf63aeSDave May PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index)
62552849c42SDave May {
6265c18a9d6SDave May 	PetscInt f;
627dbe06d34SDave May 	PetscErrorCode ierr;
62852849c42SDave May 
62984bcda08SDave May 	/* check point is valid */
630a233d522SDave May 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
631a233d522SDave May 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
63252849c42SDave May 
63352849c42SDave May 	for (f=0; f<db->nfields; f++) {
63452849c42SDave May 		DataField field = db->field[f];
63552849c42SDave May 
636dbe06d34SDave May 		ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr);
63752849c42SDave May 	}
6387fbf63aeSDave May   PetscFunctionReturn(0);
63952849c42SDave May }
64052849c42SDave May 
64152849c42SDave May /* increment */
6427fbf63aeSDave May #undef __FUNCT__
6437fbf63aeSDave May #define __FUNCT__ "DataBucketAddPoint"
6447fbf63aeSDave May PetscErrorCode DataBucketAddPoint(DataBucket db)
64552849c42SDave May {
646dbe06d34SDave May 	PetscErrorCode ierr;
647dbe06d34SDave May 
648dbe06d34SDave May 	ierr = DataBucketSetSizes(db,db->L+1,-1);CHKERRQ(ierr);
6497fbf63aeSDave May   PetscFunctionReturn(0);
65052849c42SDave May }
65152849c42SDave May 
6527fbf63aeSDave May /* decrement */
6537fbf63aeSDave May #undef __FUNCT__
6547fbf63aeSDave May #define __FUNCT__ "DataBucketRemovePoint"
6557fbf63aeSDave May PetscErrorCode DataBucketRemovePoint(DataBucket db)
6567fbf63aeSDave May {
657dbe06d34SDave May 	PetscErrorCode ierr;
658dbe06d34SDave May 
659dbe06d34SDave May 	ierr = DataBucketSetSizes(db,db->L-1,-1);CHKERRQ(ierr);
6607fbf63aeSDave May   PetscFunctionReturn(0);
6617fbf63aeSDave May }
6627fbf63aeSDave May 
6637fbf63aeSDave May #undef __FUNCT__
6647fbf63aeSDave May #define __FUNCT__ "_DataFieldViewBinary"
6657fbf63aeSDave May PetscErrorCode _DataFieldViewBinary(DataField field,FILE *fp)
66652849c42SDave May {
66752849c42SDave May 	fprintf(fp,"<DataField>\n");
66852849c42SDave May 	fprintf(fp,"%d\n", field->L);
66952849c42SDave May 	fprintf(fp,"%zu\n",field->atomic_size);
67052849c42SDave May 	fprintf(fp,"%s\n", field->registeration_function);
67152849c42SDave May 	fprintf(fp,"%s\n", field->name);
67252849c42SDave May 
67352849c42SDave May 	fwrite(field->data, field->atomic_size, field->L, fp);
67452849c42SDave May   /*
67552849c42SDave May    printf("  ** wrote %zu bytes for DataField \"%s\" \n", field->atomic_size * field->L, field->name );
67652849c42SDave May    */
67752849c42SDave May 	fprintf(fp,"\n</DataField>\n");
6787fbf63aeSDave May   PetscFunctionReturn(0);
67952849c42SDave May }
68052849c42SDave May 
6817fbf63aeSDave May #undef __FUNCT__
6827fbf63aeSDave May #define __FUNCT__ "_DataBucketRegisterFieldFromFile"
6837fbf63aeSDave May PetscErrorCode _DataBucketRegisterFieldFromFile(FILE *fp,DataBucket db)
68452849c42SDave May {
68552849c42SDave May 	PetscBool val;
68652849c42SDave May 	DataField *field;
68752849c42SDave May 
68852849c42SDave May 	DataField gfield;
68952849c42SDave May 	char dummy[100];
69052849c42SDave May 	char registeration_function[5000];
69152849c42SDave May 	char field_name[5000];
6925c18a9d6SDave May 	PetscInt L;
69352849c42SDave May 	size_t atomic_size,strL;
694dbe06d34SDave May 	PetscErrorCode ierr;
69552849c42SDave May 
69652849c42SDave May 	/* check we haven't finalised the registration of fields */
69752849c42SDave May 	/*
69852849c42SDave May    if(db->finalised==PETSC_TRUE) {
69952849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
70052849c42SDave May    ERROR();
70152849c42SDave May    }
70252849c42SDave May    */
70352849c42SDave May 
70452849c42SDave May 
70552849c42SDave May 	/* read file contents */
70652849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
70752849c42SDave May 
70852849c42SDave May 	fscanf( fp, "%d\n",&L); //printf("read(L): %d\n", L);
70952849c42SDave May 
71052849c42SDave May 	fscanf( fp, "%zu\n",&atomic_size); //printf("read(size): %zu\n",atomic_size);
71152849c42SDave May 
71252849c42SDave May 	fgets(registeration_function,4999,fp); //printf("read(reg func): %s", registeration_function );
71352849c42SDave May 	strL = strlen(registeration_function);
71452849c42SDave May 	if (strL > 1) {
71552849c42SDave May 		registeration_function[strL-1] = 0;
71652849c42SDave May 	}
71752849c42SDave May 
71852849c42SDave May 	fgets(field_name,4999,fp); //printf("read(name): %s", field_name );
71952849c42SDave May 	strL = strlen(field_name);
72052849c42SDave May 	if (strL > 1) {
72152849c42SDave May 		field_name[strL-1] = 0;
72252849c42SDave May 	}
72352849c42SDave May 
7244b46c5e1SDave May #ifdef DATA_BUCKET_LOG
725b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"  ** read L=%D; atomic_size=%zu; reg_func=\"%s\"; name=\"%s\" \n", L,atomic_size,registeration_function,field_name);
72652849c42SDave May #endif
72752849c42SDave May 
72852849c42SDave May 
72952849c42SDave May 	/* check for repeated name */
730*2635f519SDave May 	ierr = StringInList( field_name, db->nfields, (const DataField*)db->field, &val );CHKERRQ(ierr);
731a233d522SDave May 	if (val == PETSC_TRUE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot add same field twice");
73252849c42SDave May 
73352849c42SDave May 	/* create new space for data */
73452849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
73552849c42SDave May 	db->field     = field;
73652849c42SDave May 
73752849c42SDave May 	/* add field */
738dbe06d34SDave May 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield );CHKERRQ(ierr);
73952849c42SDave May 
74052849c42SDave May 	/* copy contents of file */
74152849c42SDave May 	fread(gfield->data, gfield->atomic_size, gfield->L, fp);
7424b46c5e1SDave May #ifdef DATA_BUCKET_LOG
743b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"  ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name );
74452849c42SDave May #endif
74552849c42SDave May 	/* finish reading meta data */
74652849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
74752849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
74852849c42SDave May 
74952849c42SDave May 	db->field[ db->nfields ] = gfield;
75052849c42SDave May 
75152849c42SDave May 	db->nfields++;
7527fbf63aeSDave May   PetscFunctionReturn(0);
75352849c42SDave May }
75452849c42SDave May 
7557fbf63aeSDave May #undef __FUNCT__
7567fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderWrite_v00"
7577fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderWrite_v00(FILE *fp)
75852849c42SDave May {
75952849c42SDave May 	fprintf(fp,"<DataBucketHeader>\n");
76052849c42SDave May 	fprintf(fp,"type=DataBucket\n");
76152849c42SDave May 	fprintf(fp,"format=ascii\n");
76252849c42SDave May 	fprintf(fp,"version=0.0\n");
76352849c42SDave May 	fprintf(fp,"options=\n");
76452849c42SDave May 	fprintf(fp,"</DataBucketHeader>\n");
7657fbf63aeSDave May   PetscFunctionReturn(0);
76652849c42SDave May }
7677fbf63aeSDave May 
7687fbf63aeSDave May #undef __FUNCT__
7697fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderRead_v00"
7707fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderRead_v00(FILE *fp)
77152849c42SDave May {
77252849c42SDave May 	char dummy[100];
77352849c42SDave May 	size_t strL;
77452849c42SDave May 
77552849c42SDave May 	// header open
77652849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
77752849c42SDave May 
77852849c42SDave May 	// type
77952849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
78052849c42SDave May 	strL = strlen(dummy);
78152849c42SDave May 	if (strL > 1) { dummy[strL-1] = 0; }
78252849c42SDave May 	if (strcmp(dummy,"type=DataBucket") != 0) {
783a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Data file doesn't contain a DataBucket type");
78452849c42SDave May 	}
78552849c42SDave May 
78652849c42SDave May 	// format
78752849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
78852849c42SDave May 
78952849c42SDave May 	// version
79052849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
79152849c42SDave May 	strL = strlen(dummy);
79252849c42SDave May 	if (strL > 1) { dummy[strL-1] = 0; }
79352849c42SDave May 	if (strcmp(dummy,"version=0.0") != 0) {
794a233d522SDave May 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"DataBucket file must be parsed with version=0.0 : You tried %s", dummy);
79552849c42SDave May 	}
79652849c42SDave May 
79752849c42SDave May 	// options
79852849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
79952849c42SDave May 	// header close
80052849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
8017fbf63aeSDave May   PetscFunctionReturn(0);
80252849c42SDave May }
80352849c42SDave May 
8047fbf63aeSDave May #undef __FUNCT__
8057fbf63aeSDave May #define __FUNCT__ "_DataBucketLoadFromFileBinary_SEQ"
8067fbf63aeSDave May PetscErrorCode _DataBucketLoadFromFileBinary_SEQ(const char filename[],DataBucket *_db)
80752849c42SDave May {
80852849c42SDave May 	DataBucket db;
80952849c42SDave May 	FILE *fp;
8105c18a9d6SDave May 	PetscInt L,buffer,f,nfields;
811dbe06d34SDave May 	PetscErrorCode ierr;
81252849c42SDave May 
81352849c42SDave May 
8144b46c5e1SDave May #ifdef DATA_BUCKET_LOG
815b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
81652849c42SDave May #endif
81752849c42SDave May 
81852849c42SDave May 	/* open file */
81952849c42SDave May 	fp = fopen(filename,"rb");
82052849c42SDave May 	if (fp == NULL) {
821a233d522SDave May 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file with name %s", filename);
82252849c42SDave May 	}
82352849c42SDave May 
82452849c42SDave May 	/* read header */
825dbe06d34SDave May 	ierr = _DataBucketViewAscii_HeaderRead_v00(fp);CHKERRQ(ierr);
82652849c42SDave May 
82752849c42SDave May 	fscanf(fp,"%d\n%d\n%d\n",&L,&buffer,&nfields);
82852849c42SDave May 
829dbe06d34SDave May 	ierr = DataBucketCreate(&db);CHKERRQ(ierr);
83052849c42SDave May 
83152849c42SDave May 	for (f=0; f<nfields; f++) {
832dbe06d34SDave May 		ierr = _DataBucketRegisterFieldFromFile(fp,db);CHKERRQ(ierr);
83352849c42SDave May 	}
83452849c42SDave May 	fclose(fp);
83552849c42SDave May 
836dbe06d34SDave May 	ierr = DataBucketFinalize(db);CHKERRQ(ierr);
83752849c42SDave May 
83852849c42SDave May   /*
83952849c42SDave May    DataBucketSetSizes(db,L,buffer);
84052849c42SDave May    */
84152849c42SDave May 	db->L = L;
84252849c42SDave May 	db->buffer = buffer;
84352849c42SDave May 	db->allocated = L + buffer;
84452849c42SDave May 
84552849c42SDave May 	*_db = db;
8467fbf63aeSDave May   PetscFunctionReturn(0);
84752849c42SDave May }
84852849c42SDave May 
8497fbf63aeSDave May #undef __FUNCT__
8507fbf63aeSDave May #define __FUNCT__ "DataBucketLoadFromFile"
8517fbf63aeSDave May PetscErrorCode DataBucketLoadFromFile(MPI_Comm comm,const char filename[],DataBucketViewType type,DataBucket *db)
85252849c42SDave May {
8535c18a9d6SDave May 	PetscMPIInt nproc,rank;
854997fa542SDave May 	PetscErrorCode ierr;
85552849c42SDave May 
856997fa542SDave May 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
857997fa542SDave May 	ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
85852849c42SDave May 
8594b46c5e1SDave May #ifdef DATA_BUCKET_LOG
860b3a122caSDave May 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
86152849c42SDave May #endif
86252849c42SDave May 	if (type == DATABUCKET_VIEW_STDOUT) {
86352849c42SDave May 
86452849c42SDave May 	} else if (type == DATABUCKET_VIEW_ASCII) {
865a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
86652849c42SDave May 	} else if (type == DATABUCKET_VIEW_BINARY) {
86752849c42SDave May 		if (nproc == 1) {
868dbe06d34SDave May 			ierr = _DataBucketLoadFromFileBinary_SEQ(filename,db);CHKERRQ(ierr);
86952849c42SDave May 		} else {
87052849c42SDave May 			char *name;
87152849c42SDave May 
87252849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
873dbe06d34SDave May 			ierr = _DataBucketLoadFromFileBinary_SEQ(name,db);CHKERRQ(ierr);
87452849c42SDave May 			free(name);
87552849c42SDave May 		}
87652849c42SDave May 	} else {
877a233d522SDave May 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer requested");
87852849c42SDave May 	}
8797fbf63aeSDave May   PetscFunctionReturn(0);
88052849c42SDave May }
88152849c42SDave May 
8827fbf63aeSDave May #undef __FUNCT__
8837fbf63aeSDave May #define __FUNCT__ "_DataBucketViewBinary"
8847fbf63aeSDave May PetscErrorCode _DataBucketViewBinary(DataBucket db,const char filename[])
88552849c42SDave May {
88652849c42SDave May 	FILE *fp = NULL;
8875c18a9d6SDave May 	PetscInt f;
888dbe06d34SDave May 	PetscErrorCode ierr;
88952849c42SDave May 
89052849c42SDave May 	fp = fopen(filename,"wb");
891a233d522SDave May 	if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file with name %s", filename);
89252849c42SDave May 
89352849c42SDave May 	/* db header */
894dbe06d34SDave May 	ierr =_DataBucketViewAscii_HeaderWrite_v00(fp);CHKERRQ(ierr);
89552849c42SDave May 
89652849c42SDave May 	/* meta-data */
89752849c42SDave May 	fprintf(fp,"%d\n%d\n%d\n", db->L,db->buffer,db->nfields);
89852849c42SDave May 
89952849c42SDave May 	for (f=0; f<db->nfields; f++) {
90052849c42SDave May     /* load datafields */
901dbe06d34SDave May 		ierr = _DataFieldViewBinary(db->field[f],fp);CHKERRQ(ierr);
90252849c42SDave May 	}
90352849c42SDave May 
90452849c42SDave May 	fclose(fp);
9057fbf63aeSDave May   PetscFunctionReturn(0);
90652849c42SDave May }
90752849c42SDave May 
9087fbf63aeSDave May #undef __FUNCT__
9097fbf63aeSDave May #define __FUNCT__ "DataBucketView_SEQ"
9107fbf63aeSDave May PetscErrorCode DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type)
91152849c42SDave May {
912dbe06d34SDave May   PetscErrorCode ierr;
913dbe06d34SDave May 
91452849c42SDave May 	switch (type) {
91552849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
91652849c42SDave May 		{
9175c18a9d6SDave May 			PetscInt f;
91852849c42SDave May 			double memory_usage_total = 0.0;
91952849c42SDave May 
920b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"DataBucketView(SEQ): (\"%s\")\n",filename);
921b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  L                  = %D \n", db->L );
922b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  buffer             = %D \n", db->buffer );
923b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  allocated          = %D \n", db->allocated );
92452849c42SDave May 
925b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  nfields registered = %D \n", db->nfields );
92652849c42SDave May 			for (f=0; f<db->nfields; f++) {
92752849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
92852849c42SDave May 
929b3a122caSDave May 				PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) \n", f, db->field[f]->name, memory_usage_f  );
93052849c42SDave May 				memory_usage_total += memory_usage_f;
93152849c42SDave May 			}
932b3a122caSDave May 			PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) \n", memory_usage_total );
93352849c42SDave May 		}
93452849c42SDave May 			break;
93552849c42SDave May 
93652849c42SDave May 		case DATABUCKET_VIEW_ASCII:
93752849c42SDave May 		{
938a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
93952849c42SDave May 		}
94052849c42SDave May 			break;
94152849c42SDave May 
94252849c42SDave May 		case DATABUCKET_VIEW_BINARY:
94352849c42SDave May 		{
944dbe06d34SDave May 			ierr = _DataBucketViewBinary(db,filename);CHKERRQ(ierr);
94552849c42SDave May 		}
94652849c42SDave May 			break;
94752849c42SDave May 
94852849c42SDave May 		case DATABUCKET_VIEW_HDF5:
94952849c42SDave May 		{
950a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No HDF5 support");
95152849c42SDave May 		}
95252849c42SDave May 			break;
95352849c42SDave May 
95452849c42SDave May 		default:
955a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
95652849c42SDave May 			break;
95752849c42SDave May 	}
9587fbf63aeSDave May   PetscFunctionReturn(0);
95952849c42SDave May }
96052849c42SDave May 
9617fbf63aeSDave May #undef __FUNCT__
9627fbf63aeSDave May #define __FUNCT__ "DataBucketView_MPI"
9637fbf63aeSDave May PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
96452849c42SDave May {
965997fa542SDave May   PetscErrorCode ierr;
966997fa542SDave May 
96752849c42SDave May 	switch (type) {
96852849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
96952849c42SDave May 		{
9705c18a9d6SDave May 			PetscInt f;
9715c18a9d6SDave May 			PetscInt L,buffer,allocated;
97252849c42SDave May 			double memory_usage_total,memory_usage_total_local = 0.0;
9735c18a9d6SDave May 			PetscMPIInt rank;
97452849c42SDave May 
975997fa542SDave May 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
97652849c42SDave May 
977dbe06d34SDave May 			ierr = DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated);CHKERRQ(ierr);
97852849c42SDave May 
97952849c42SDave May 			for (f=0; f<db->nfields; f++) {
98052849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
98152849c42SDave May 
98252849c42SDave May 				memory_usage_total_local += memory_usage_f;
98352849c42SDave May 			}
984dbe06d34SDave May 			ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
98552849c42SDave May 
98652849c42SDave May 			if (rank == 0) {
9875c18a9d6SDave May 				PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename);
9885c18a9d6SDave May 				PetscPrintf(comm,"  L                  = %D \n", L );
9895c18a9d6SDave May 				PetscPrintf(comm,"  buffer (max)       = %D \n", buffer );
9905c18a9d6SDave May 				PetscPrintf(comm,"  allocated          = %D \n", allocated );
99152849c42SDave May 
9925c18a9d6SDave May 				PetscPrintf(comm,"  nfields registered = %D \n", db->nfields );
99352849c42SDave May 				for (f=0; f<db->nfields; f++) {
99452849c42SDave May 					double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
99552849c42SDave May 
996b3a122caSDave May 					PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) : rank0\n", f, db->field[f]->name, memory_usage_f  );
99752849c42SDave May 				}
99852849c42SDave May 
999b3a122caSDave May 				PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) : collective\n", memory_usage_total );
100052849c42SDave May 			}
100152849c42SDave May 
100252849c42SDave May 		}
100352849c42SDave May 			break;
100452849c42SDave May 
100552849c42SDave May 		case DATABUCKET_VIEW_ASCII:
100652849c42SDave May 		{
1007a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying data structure");
100852849c42SDave May 		}
100952849c42SDave May 			break;
101052849c42SDave May 
101152849c42SDave May 		case DATABUCKET_VIEW_BINARY:
101252849c42SDave May 		{
101352849c42SDave May 			char *name;
10145c18a9d6SDave May 			PetscMPIInt rank;
101552849c42SDave May 
101652849c42SDave May 			/* create correct extension */
1017997fa542SDave May 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
101852849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
101952849c42SDave May 
1020dbe06d34SDave May 			ierr = _DataBucketViewBinary(db,name);CHKERRQ(ierr);
102152849c42SDave May 
102252849c42SDave May 			free(name);
102352849c42SDave May 		}
102452849c42SDave May 			break;
102552849c42SDave May 
102652849c42SDave May 		case DATABUCKET_VIEW_HDF5:
102752849c42SDave May 		{
1028a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5");
102952849c42SDave May 		}
103052849c42SDave May 			break;
103152849c42SDave May 
103252849c42SDave May 		default:
1033a233d522SDave May 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
103452849c42SDave May 			break;
103552849c42SDave May 	}
10367fbf63aeSDave May   PetscFunctionReturn(0);
103752849c42SDave May }
103852849c42SDave May 
10397fbf63aeSDave May #undef __FUNCT__
10407fbf63aeSDave May #define __FUNCT__ "DataBucketView"
10417fbf63aeSDave May PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
104252849c42SDave May {
10435c18a9d6SDave May 	PetscMPIInt nproc;
1044997fa542SDave May   PetscErrorCode ierr;
104552849c42SDave May 
1046997fa542SDave May 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
104752849c42SDave May 	if (nproc == 1) {
1048dbe06d34SDave May 		ierr =DataBucketView_SEQ(db,filename,type);CHKERRQ(ierr);
104952849c42SDave May 	} else {
1050dbe06d34SDave May 		ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
105152849c42SDave May 	}
10527fbf63aeSDave May   PetscFunctionReturn(0);
105352849c42SDave May }
105452849c42SDave May 
10557fbf63aeSDave May #undef __FUNCT__
10567fbf63aeSDave May #define __FUNCT__ "DataBucketDuplicateFields"
10577fbf63aeSDave May PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
105852849c42SDave May {
105952849c42SDave May 	DataBucket db2;
10605c18a9d6SDave May 	PetscInt f;
1061dbe06d34SDave May 	PetscErrorCode ierr;
106252849c42SDave May 
1063dbe06d34SDave May 	ierr = DataBucketCreate(&db2);CHKERRQ(ierr);
106452849c42SDave May 
106552849c42SDave May 	/* copy contents from dbA into db2 */
106652849c42SDave May 	for (f=0; f<dbA->nfields; f++) {
106752849c42SDave May 		DataField field;
106852849c42SDave May 		size_t    atomic_size;
106952849c42SDave May 		char      *name;
107052849c42SDave May 
107152849c42SDave May 		field = dbA->field[f];
107252849c42SDave May 
107352849c42SDave May 		atomic_size = field->atomic_size;
107452849c42SDave May 		name        = field->name;
107552849c42SDave May 
1076dbe06d34SDave May 		ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
107752849c42SDave May 	}
1078dbe06d34SDave May 	ierr = DataBucketFinalize(db2);CHKERRQ(ierr);
1079dbe06d34SDave May 	ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
108052849c42SDave May 
108152849c42SDave May 	/* set pointer */
108252849c42SDave May 	*dbB = db2;
10837fbf63aeSDave May   PetscFunctionReturn(0);
108452849c42SDave May }
108552849c42SDave May 
108652849c42SDave May /*
108752849c42SDave May  Insert points from db2 into db1
108852849c42SDave May  db1 <<== db2
108952849c42SDave May  */
10907fbf63aeSDave May #undef __FUNCT__
10917fbf63aeSDave May #define __FUNCT__ "DataBucketInsertValues"
10927fbf63aeSDave May PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2)
109352849c42SDave May {
10945c18a9d6SDave May 	PetscInt n_mp_points1,n_mp_points2;
10955c18a9d6SDave May 	PetscInt n_mp_points1_new,p;
1096dbe06d34SDave May 	PetscErrorCode ierr;
109752849c42SDave May 
1098dbe06d34SDave May 	ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr);
1099dbe06d34SDave May 	ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr);
110052849c42SDave May 
110152849c42SDave May 	n_mp_points1_new = n_mp_points1 + n_mp_points2;
1102dbe06d34SDave May 	ierr = DataBucketSetSizes(db1,n_mp_points1_new,-1);CHKERRQ(ierr);
110352849c42SDave May 
110452849c42SDave May 	for (p=0; p<n_mp_points2; p++) {
110552849c42SDave May 		// db1 <<== db2 //
1106dbe06d34SDave May 		ierr =DataBucketCopyPoint( db2,p, db1,(n_mp_points1 + p) );CHKERRQ(ierr);
110752849c42SDave May 	}
11087fbf63aeSDave May   PetscFunctionReturn(0);
110952849c42SDave May }
111052849c42SDave May 
111152849c42SDave May /* helpers for parallel send/recv */
11127fbf63aeSDave May #undef __FUNCT__
11137fbf63aeSDave May #define __FUNCT__ "DataBucketCreatePackedArray"
11147fbf63aeSDave May PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
111552849c42SDave May {
11165c18a9d6SDave May   PetscInt       f;
111752849c42SDave May   size_t    sizeof_marker_contents;
111852849c42SDave May   void      *buffer;
111952849c42SDave May 
112052849c42SDave May   sizeof_marker_contents = 0;
112152849c42SDave May   for (f=0; f<db->nfields; f++) {
112252849c42SDave May     DataField df = db->field[f];
112352849c42SDave May 
112452849c42SDave May     sizeof_marker_contents += df->atomic_size;
112552849c42SDave May   }
112652849c42SDave May 
112752849c42SDave May   buffer = malloc(sizeof_marker_contents);
112852849c42SDave May   memset(buffer,0,sizeof_marker_contents);
112952849c42SDave May 
113052849c42SDave May   if (bytes) { *bytes = sizeof_marker_contents; }
113152849c42SDave May   if (buf)   { *buf   = buffer; }
11327fbf63aeSDave May   PetscFunctionReturn(0);
113352849c42SDave May }
113452849c42SDave May 
11357fbf63aeSDave May #undef __FUNCT__
11367fbf63aeSDave May #define __FUNCT__ "DataBucketDestroyPackedArray"
11377fbf63aeSDave May PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf)
113852849c42SDave May {
113952849c42SDave May   if (buf) {
114052849c42SDave May     free(*buf);
114152849c42SDave May     *buf = NULL;
114252849c42SDave May   }
11437fbf63aeSDave May   PetscFunctionReturn(0);
114452849c42SDave May }
114552849c42SDave May 
11467fbf63aeSDave May #undef __FUNCT__
11477fbf63aeSDave May #define __FUNCT__ "DataBucketFillPackedArray"
11487fbf63aeSDave May PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
114952849c42SDave May {
11505c18a9d6SDave May   PetscInt    f;
115152849c42SDave May   void   *data,*data_p;
115252849c42SDave May   size_t asize,offset;
115352849c42SDave May 
115452849c42SDave May   offset = 0;
115552849c42SDave May   for (f=0; f<db->nfields; f++) {
115652849c42SDave May     DataField df = db->field[f];
115752849c42SDave May 
115852849c42SDave May     asize = df->atomic_size;
115952849c42SDave May 
116052849c42SDave May     data = (void*)( df->data );
116152849c42SDave May     data_p = (void*)( (char*)data + index*asize );
116252849c42SDave May 
116352849c42SDave May     memcpy( (void*)((char*)buf + offset),  data_p,  asize);
116452849c42SDave May     offset = offset + asize;
116552849c42SDave May   }
11667fbf63aeSDave May   PetscFunctionReturn(0);
116752849c42SDave May }
116852849c42SDave May 
11697fbf63aeSDave May #undef __FUNCT__
11707fbf63aeSDave May #define __FUNCT__ "DataBucketInsertPackedArray"
11717fbf63aeSDave May PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
117252849c42SDave May {
11735c18a9d6SDave May   PetscInt f;
117452849c42SDave May   void *data_p;
117552849c42SDave May   size_t offset;
1176dbe06d34SDave May 	PetscErrorCode ierr;
117752849c42SDave May 
117852849c42SDave May   offset = 0;
117952849c42SDave May   for (f=0; f<db->nfields; f++) {
118052849c42SDave May     DataField df = db->field[f];
118152849c42SDave May 
118252849c42SDave May     data_p = (void*)( (char*)data + offset );
118352849c42SDave May 
1184dbe06d34SDave May     ierr = DataFieldInsertPoint(df, idx, (void*)data_p );CHKERRQ(ierr);
118552849c42SDave May     offset = offset + df->atomic_size;
118652849c42SDave May   }
11877fbf63aeSDave May   PetscFunctionReturn(0);
118852849c42SDave May }
1189