xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision 4b46c5e15be7c29a43df8a7ac989acd1e2bd59f6)
152849c42SDave May 
252849c42SDave May #include "data_bucket.h"
352849c42SDave May 
452849c42SDave May /* string helpers */
55c18a9d6SDave May void StringInList( const char name[], const PetscInt N, const DataField gfield[], PetscBool *val )
652849c42SDave May {
75c18a9d6SDave May 	PetscInt i;
852849c42SDave May 
952849c42SDave May 	*val = PETSC_FALSE;
1052849c42SDave May 	for( i=0; i<N; i++ ) {
1152849c42SDave May 		if( strcmp( name, gfield[i]->name ) == 0 ) {
1252849c42SDave May 			*val = PETSC_TRUE;
1352849c42SDave May 			return;
1452849c42SDave May 		}
1552849c42SDave May 	}
1652849c42SDave May }
1752849c42SDave May 
185c18a9d6SDave May void StringFindInList( const char name[], const PetscInt N, const DataField gfield[], PetscInt *index )
1952849c42SDave May {
205c18a9d6SDave May 	PetscInt i;
2152849c42SDave May 
2252849c42SDave May 	*index = -1;
2352849c42SDave May 	for( i=0; i<N; i++ ) {
2452849c42SDave May 		if( strcmp( name, gfield[i]->name ) == 0 ) {
2552849c42SDave May 			*index = i;
2652849c42SDave May 			return;
2752849c42SDave May 		}
2852849c42SDave May 	}
2952849c42SDave May }
3052849c42SDave May 
315c18a9d6SDave May void DataFieldCreate( const char registeration_function[], const char name[], const size_t size, const PetscInt L, DataField *DF )
3252849c42SDave May {
3352849c42SDave May 	DataField df;
3452849c42SDave May 
3552849c42SDave May 	df = malloc( sizeof(struct _p_DataField) );
3652849c42SDave May 	memset( df, 0, sizeof(struct _p_DataField) );
3752849c42SDave May 
3852849c42SDave May 
3952849c42SDave May 	asprintf( &df->registeration_function, "%s", registeration_function );
4052849c42SDave May 	asprintf( &df->name, "%s", name );
4152849c42SDave May 	df->atomic_size = size;
4252849c42SDave May 	df->L = L;
4352849c42SDave May 
4452849c42SDave May 	df->data = malloc( size * L ); /* allocate something so we don't have to reallocate */
4552849c42SDave May 	memset( df->data, 0, size * L );
4652849c42SDave May 
4752849c42SDave May 	*DF = df;
4852849c42SDave May }
4952849c42SDave May 
5052849c42SDave May void DataFieldDestroy( DataField *DF )
5152849c42SDave May {
5252849c42SDave May 	DataField df = *DF;
5352849c42SDave May 
5452849c42SDave May 	free( df->registeration_function );
5552849c42SDave May 	free( df->name );
5652849c42SDave May 	free( df->data );
5752849c42SDave May 	free(df);
5852849c42SDave May 
5952849c42SDave May 	*DF = NULL;
6052849c42SDave May }
6152849c42SDave May 
6252849c42SDave May /* data bucket */
6352849c42SDave May void DataBucketCreate( DataBucket *DB )
6452849c42SDave May {
6552849c42SDave May 	DataBucket db;
6652849c42SDave May 
6752849c42SDave May 
6852849c42SDave May 	db = malloc( sizeof(struct _p_DataBucket) );
6952849c42SDave May 	memset( db, 0, sizeof(struct _p_DataBucket) );
7052849c42SDave May 
7152849c42SDave May 	db->finalised = PETSC_FALSE;
7252849c42SDave May 
7352849c42SDave May 	/* create empty spaces for fields */
7452849c42SDave May 	db->L         = 0;
7552849c42SDave May 	db->buffer    = 1;
7652849c42SDave May 	db->allocated = 1;
7752849c42SDave May 
7852849c42SDave May 	db->nfields   = 0;
7952849c42SDave May 	db->field     = malloc(sizeof(DataField));
8052849c42SDave May 
8152849c42SDave May 	*DB = db;
8252849c42SDave May }
8352849c42SDave May 
8452849c42SDave May void DataBucketDestroy( DataBucket *DB )
8552849c42SDave May {
8652849c42SDave May 	DataBucket db = *DB;
875c18a9d6SDave May 	PetscInt f;
8852849c42SDave May 
8952849c42SDave May 	/* release fields */
9052849c42SDave May 	for( f=0; f<db->nfields; f++ ) {
9152849c42SDave May 		DataFieldDestroy(&db->field[f]);
9252849c42SDave May 	}
9352849c42SDave May 
9452849c42SDave May 	/* this will catch the initially allocated objects in the event that no fields are registered */
9552849c42SDave May 	if(db->field!=NULL) {
9652849c42SDave May 		free(db->field);
9752849c42SDave May 	}
9852849c42SDave May 
9952849c42SDave May 	free(db);
10052849c42SDave May 
10152849c42SDave May 	*DB = NULL;
10252849c42SDave May }
10352849c42SDave May 
10452849c42SDave May void _DataBucketRegisterField(
10552849c42SDave May 						DataBucket db,
10652849c42SDave May 						const char registeration_function[],
10752849c42SDave May 						const char field_name[],
10852849c42SDave May 						size_t atomic_size, DataField *_gfield )
10952849c42SDave May {
11052849c42SDave May 	PetscBool val;
11152849c42SDave May 	DataField *field,fp;
11252849c42SDave May 
11352849c42SDave May 	/* check we haven't finalised the registration of fields */
11452849c42SDave May 	/*
11552849c42SDave May 	if(db->finalised==PETSC_TRUE) {
11652849c42SDave May 		printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
11752849c42SDave May 		ERROR();
11852849c42SDave May 	}
11952849c42SDave May 	*/
12052849c42SDave May 
12152849c42SDave May 	/* check for repeated name */
12252849c42SDave May 	StringInList( field_name, db->nfields, (const DataField*)db->field, &val );
12352849c42SDave May 	if(val == PETSC_TRUE ) {
12452849c42SDave May 		printf("ERROR: Cannot add same field twice\n");
12552849c42SDave May 		ERROR();
12652849c42SDave May 	}
12752849c42SDave May 
12852849c42SDave May 	/* create new space for data */
12952849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
13052849c42SDave May 	db->field     = field;
13152849c42SDave May 
13252849c42SDave May 	/* add field */
13352849c42SDave May 	DataFieldCreate( registeration_function, field_name, atomic_size, db->allocated, &fp );
13452849c42SDave May 	db->field[ db->nfields ] = fp;
13552849c42SDave May 
13652849c42SDave May 	db->nfields++;
13752849c42SDave May 
13852849c42SDave May 	if(_gfield!=NULL){
13952849c42SDave May 		*_gfield = fp;
14052849c42SDave May 	}
14152849c42SDave May }
14252849c42SDave May 
14352849c42SDave May /*
14452849c42SDave May #define DataBucketRegisterField(db,name,size,k) {\
14552849c42SDave May   char *location;\
14652849c42SDave May   asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
14752849c42SDave May   _DataBucketRegisterField( (db), location, (name), (size), (k) );\
14852849c42SDave May   free(location);\
14952849c42SDave May }
15052849c42SDave May */
15152849c42SDave May 
15252849c42SDave May void DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
15352849c42SDave May {
1545c18a9d6SDave May 	PetscInt idx;
15552849c42SDave May 	PetscBool found;
15652849c42SDave May 
15752849c42SDave May 	StringInList(name,db->nfields,(const DataField*)db->field,&found);
15852849c42SDave May 	if(found==PETSC_FALSE) {
15952849c42SDave May 		printf("ERROR: Cannot find DataField with name %s \n", name );
16052849c42SDave May 		ERROR();
16152849c42SDave May 	}
16252849c42SDave May 	StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);
16352849c42SDave May 
16452849c42SDave May 	*gfield = db->field[idx];
16552849c42SDave May }
16652849c42SDave May 
16752849c42SDave May void DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
16852849c42SDave May {
16952849c42SDave May 	*found = PETSC_FALSE;
17052849c42SDave May 	StringInList(name,db->nfields,(const DataField*)db->field,found);
17152849c42SDave May }
17252849c42SDave May 
17352849c42SDave May void DataBucketFinalize(DataBucket db)
17452849c42SDave May {
17552849c42SDave May 	db->finalised = PETSC_TRUE;
17652849c42SDave May }
17752849c42SDave May 
1785c18a9d6SDave May void DataFieldGetNumEntries(DataField df, PetscInt *sum)
17952849c42SDave May {
18052849c42SDave May 	*sum = df->L;
18152849c42SDave May }
18252849c42SDave May 
1835c18a9d6SDave May void DataFieldSetSize( DataField df, const PetscInt new_L )
18452849c42SDave May {
18552849c42SDave May 	void *tmp_data;
18652849c42SDave May 
18752849c42SDave May 	if( new_L <= 0 ) {
18852849c42SDave May 		printf("ERROR: Cannot set size of DataField to be <= 0 \n");
18952849c42SDave May 		ERROR();
19052849c42SDave May 	}
19152849c42SDave May 	if( new_L == df->L ) return;
19252849c42SDave May 
19352849c42SDave May 	if( new_L > df->L ) {
19452849c42SDave May 
19552849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L) );
19652849c42SDave May 		df->data = tmp_data;
19752849c42SDave May 
19852849c42SDave May 		/* init new contents */
19952849c42SDave May 		memset( ( ((char*)df->data)+df->L*df->atomic_size), 0, (new_L-df->L)*df->atomic_size );
20052849c42SDave May 
20152849c42SDave May 	}
20252849c42SDave May 	else {
20352849c42SDave May 		/* reallocate pointer list, add +1 in case new_L = 0 */
20452849c42SDave May 		tmp_data = realloc( df->data, df->atomic_size * (new_L+1) );
20552849c42SDave May 		df->data = tmp_data;
20652849c42SDave May 	}
20752849c42SDave May 
20852849c42SDave May 	df->L = new_L;
20952849c42SDave May }
21052849c42SDave May 
2115c18a9d6SDave May void DataFieldZeroBlock( DataField df, const PetscInt start, const PetscInt end )
21252849c42SDave May {
21352849c42SDave May 	if( start > end ) {
21452849c42SDave May 		printf("ERROR: Cannot zero a block of entries if start(%d) > end(%d) \n",start,end);
21552849c42SDave May 		ERROR();
21652849c42SDave May 	}
21752849c42SDave May 	if( start < 0 ) {
21852849c42SDave May 		printf("ERROR: Cannot zero a block of entries if start(%d) < 0 \n",start);
21952849c42SDave May 		ERROR();
22052849c42SDave May 	}
22152849c42SDave May 	if( end > df->L ) {
22252849c42SDave May 		printf("ERROR: Cannot zero a block of entries if end(%d) >= array size(%d) \n",end,df->L);
22352849c42SDave May 		ERROR();
22452849c42SDave May 	}
22552849c42SDave May 
22652849c42SDave May 	memset( ( ((char*)df->data)+start*df->atomic_size), 0, (end-start)*df->atomic_size );
22752849c42SDave May }
22852849c42SDave May 
22952849c42SDave May /*
23052849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
23152849c42SDave May  */
2325c18a9d6SDave May void DataBucketSetSizes( DataBucket db, const PetscInt L, const PetscInt buffer )
23352849c42SDave May {
2345c18a9d6SDave May 	PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
23552849c42SDave May 
23652849c42SDave May 
23752849c42SDave May 	if( db->finalised == PETSC_FALSE ) {
23852849c42SDave May 		printf("ERROR: You must call DataBucketFinalize() before DataBucketSetSizes() \n");
23952849c42SDave May 		ERROR();
24052849c42SDave May 	}
24152849c42SDave May 
24252849c42SDave May 	current_allocated = db->allocated;
24352849c42SDave May 
24452849c42SDave May 	new_used   = L;
24552849c42SDave May 	new_unused = current_allocated - new_used;
24652849c42SDave May 	new_buffer = db->buffer;
24752849c42SDave May 	if( buffer >= 0 ) { /* update the buffer value */
24852849c42SDave May 		new_buffer = buffer;
24952849c42SDave May 	}
25052849c42SDave May 	new_allocated = new_used + new_buffer;
25152849c42SDave May 
25252849c42SDave May 	/* action */
25352849c42SDave May 	if ( new_allocated > current_allocated ) {
25452849c42SDave May 		/* increase size to new_used + new_buffer */
25552849c42SDave May 		for( f=0; f<db->nfields; f++ ) {
25652849c42SDave May 			DataFieldSetSize( db->field[f], new_allocated );
25752849c42SDave May 		}
25852849c42SDave May 
25952849c42SDave May 		db->L         = new_used;
26052849c42SDave May 		db->buffer    = new_buffer;
26152849c42SDave May 		db->allocated = new_used + new_buffer;
26252849c42SDave May 	}
26352849c42SDave May 	else {
26452849c42SDave May 		if( new_unused > 2 * new_buffer ) {
26552849c42SDave May 
26652849c42SDave May 			/* shrink array to new_used + new_buffer */
26752849c42SDave May 			for( f=0; f<db->nfields; f++ ) {
26852849c42SDave May 				DataFieldSetSize( db->field[f], new_allocated );
26952849c42SDave May 			}
27052849c42SDave May 
27152849c42SDave May 			db->L         = new_used;
27252849c42SDave May 			db->buffer    = new_buffer;
27352849c42SDave May 			db->allocated = new_used + new_buffer;
27452849c42SDave May 		}
27552849c42SDave May 		else {
27652849c42SDave May 			db->L      = new_used;
27752849c42SDave May 			db->buffer = new_buffer;
27852849c42SDave May 		}
27952849c42SDave May 	}
28052849c42SDave May 
28152849c42SDave May 	/* zero all entries from db->L to db->allocated */
28252849c42SDave May 	for( f=0; f<db->nfields; f++ ) {
28352849c42SDave May 		DataField field = db->field[f];
28452849c42SDave May 		DataFieldZeroBlock(field, db->L,db->allocated);
28552849c42SDave May 	}
28652849c42SDave May }
28752849c42SDave May 
2885c18a9d6SDave May void DataBucketSetInitialSizes( DataBucket db, const PetscInt L, const PetscInt buffer )
28952849c42SDave May {
2905c18a9d6SDave May 	PetscInt f;
29152849c42SDave May 	DataBucketSetSizes(db,L,buffer);
29252849c42SDave May 
29352849c42SDave May 	for( f=0; f<db->nfields; f++ ) {
29452849c42SDave May 		DataField field = db->field[f];
29552849c42SDave May 		DataFieldZeroBlock(field,0,db->allocated);
29652849c42SDave May 	}
29752849c42SDave May }
29852849c42SDave May 
2995c18a9d6SDave May void DataBucketGetSizes( DataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated )
30052849c42SDave May {
30152849c42SDave May 	if (L) { *L = db->L; }
30252849c42SDave May 	if (buffer) { *buffer = db->buffer; }
30352849c42SDave May 	if (allocated) { *allocated = db->allocated; }
30452849c42SDave May }
30552849c42SDave May 
3065c18a9d6SDave May void DataBucketGetGlobalSizes(MPI_Comm comm, DataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated )
30752849c42SDave May {
3085c18a9d6SDave May 	PetscInt _L,_buffer,_allocated;
3095c18a9d6SDave May 	PetscInt ierr;
31052849c42SDave May 
3115c18a9d6SDave May 	_L = db->L;
3125c18a9d6SDave May 	_buffer = db->buffer;
3135c18a9d6SDave May 	_allocated = db->allocated;
31452849c42SDave May 
3155c18a9d6SDave May 	if (L) {         ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm); }
3165c18a9d6SDave May 	if (buffer) {    ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm); }
3175c18a9d6SDave May 	if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm); }
31852849c42SDave May }
31952849c42SDave May 
3205c18a9d6SDave May void DataBucketGetDataFields( DataBucket db, PetscInt *L, DataField *fields[] )
32152849c42SDave May {
32252849c42SDave May 	if(L){      *L      = db->nfields; }
32352849c42SDave May 	if(fields){ *fields = db->field; }
32452849c42SDave May }
32552849c42SDave May 
32652849c42SDave May void DataFieldGetAccess( const DataField gfield )
32752849c42SDave May {
32852849c42SDave May 	if(gfield->active==PETSC_TRUE) {
32952849c42SDave May 		printf("ERROR: Field \"%s\" is already active. You must call DataFieldRestoreAccess()\n", gfield->name );
33052849c42SDave May 		ERROR();
33152849c42SDave May 	}
33252849c42SDave May 	gfield->active = PETSC_TRUE;
33352849c42SDave May }
33452849c42SDave May 
3355c18a9d6SDave May void DataFieldAccessPoint( const DataField gfield, const PetscInt pid, void **ctx_p )
33652849c42SDave May {
33752849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
33852849c42SDave May 	/* debug mode */
3395c18a9d6SDave May 	/* check poPetscInt is valid */
34052849c42SDave May 	if( pid < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR();  }
34152849c42SDave May 	if( pid >= gfield->L ){ printf("ERROR: index must be < %d\n",gfield->L); ERROR(); }
34252849c42SDave May 
34352849c42SDave May 	if(gfield->active==PETSC_FALSE) {
3445c18a9d6SDave May 		printf("ERROR: Field \"%s\" is not active. You must call DataFieldGetAccess() before poPetscInt data can be retrivied\n",gfield->name);
34552849c42SDave May 		ERROR();
34652849c42SDave May 	}
34752849c42SDave May #endif
34852849c42SDave May 
34952849c42SDave May 	//*ctx_p  = (void*)( ((char*)gfield->data) + pid * gfield->atomic_size);
35052849c42SDave May 	*ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
35152849c42SDave May }
35252849c42SDave May 
3535c18a9d6SDave May void DataFieldAccessPointOffset( const DataField gfield, const size_t offset, const PetscInt pid, void **ctx_p )
35452849c42SDave May {
35552849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
35652849c42SDave May 	/* debug mode */
35752849c42SDave May 
3585c18a9d6SDave May 	/* check poPetscInt is valid */
3595c18a9d6SDave May 	/* if( offset < 0 ){ printf("ERROR: offset must be >= 0\n"); ERROR();  } *//* Note compiler realizes this can never happen with an unsigned PetscInt */
36052849c42SDave May 	if( offset >= gfield->atomic_size ){ printf("ERROR: offset must be < %zu\n",gfield->atomic_size); ERROR(); }
36152849c42SDave May 
3625c18a9d6SDave May 	/* check poPetscInt is valid */
36352849c42SDave May 	if( pid < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR();  }
36452849c42SDave May 	if( pid >= gfield->L ){ printf("ERROR: index must be < %d\n",gfield->L); ERROR(); }
36552849c42SDave May 
36652849c42SDave May 	if(gfield->active==PETSC_FALSE) {
3675c18a9d6SDave May 		printf("ERROR: Field \"%s\" is not active. You must call DataFieldGetAccess() before poPetscInt data can be retrivied\n",gfield->name);
36852849c42SDave May 		ERROR();
36952849c42SDave May 	}
37052849c42SDave May #endif
37152849c42SDave May 
37252849c42SDave May 	*ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
37352849c42SDave May }
37452849c42SDave May 
37552849c42SDave May void DataFieldRestoreAccess( DataField gfield )
37652849c42SDave May {
37752849c42SDave May 	if(gfield->active==PETSC_FALSE) {
37852849c42SDave May 		printf("ERROR: Field \"%s\" is not active. You must call DataFieldGetAccess()\n", gfield->name );
37952849c42SDave May 		ERROR();
38052849c42SDave May 	}
38152849c42SDave May 	gfield->active = PETSC_FALSE;
38252849c42SDave May }
38352849c42SDave May 
38452849c42SDave May void DataFieldVerifyAccess( const DataField gfield, const size_t size)
38552849c42SDave May {
38652849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
38752849c42SDave May 	if(gfield->atomic_size != size ) {
38852849c42SDave May         printf("ERROR: Field \"%s\" must be mapped to %zu bytes, your intended structure is %zu bytes in length.\n",
38952849c42SDave May                gfield->name, gfield->atomic_size, size );
39052849c42SDave May 		ERROR();
39152849c42SDave May 	}
39252849c42SDave May #endif
39352849c42SDave May }
39452849c42SDave May 
39552849c42SDave May void DataFieldGetAtomicSize(const DataField gfield,size_t *size)
39652849c42SDave May {
39752849c42SDave May     if (size) { *size = gfield->atomic_size; }
39852849c42SDave May }
39952849c42SDave May 
40052849c42SDave May void DataFieldGetEntries(const DataField gfield,void **data)
40152849c42SDave May {
40252849c42SDave May     if (data) {
40352849c42SDave May         *data = gfield->data;
40452849c42SDave May     }
40552849c42SDave May }
40652849c42SDave May 
40752849c42SDave May void DataFieldRestoreEntries(const DataField gfield,void **data)
40852849c42SDave May {
40952849c42SDave May     if (data) {
41052849c42SDave May         *data = NULL;
41152849c42SDave May     }
41252849c42SDave May }
41352849c42SDave May 
41452849c42SDave May /* y = x */
4155c18a9d6SDave May void DataBucketCopyPoint( const DataBucket xb, const PetscInt pid_x,
4165c18a9d6SDave May 												  const DataBucket yb, const PetscInt pid_y )
41752849c42SDave May {
4185c18a9d6SDave May 	PetscInt f;
41952849c42SDave May 	for( f=0; f<xb->nfields; f++ ) {
42052849c42SDave May 		void *dest;
42152849c42SDave May 		void *src;
42252849c42SDave May 
42352849c42SDave May 		DataFieldGetAccess( xb->field[f] );
42452849c42SDave May 		if (xb!=yb) { DataFieldGetAccess( yb->field[f] ); }
42552849c42SDave May 
42652849c42SDave May 		DataFieldAccessPoint( xb->field[f],pid_x, &src );
42752849c42SDave May 		DataFieldAccessPoint( yb->field[f],pid_y, &dest );
42852849c42SDave May 
42952849c42SDave May 		memcpy( dest, src, xb->field[f]->atomic_size );
43052849c42SDave May 
43152849c42SDave May 		DataFieldRestoreAccess( xb->field[f] );
43252849c42SDave May 		if (xb!=yb) { DataFieldRestoreAccess( yb->field[f] ); }
43352849c42SDave May 	}
43452849c42SDave May 
43552849c42SDave May }
43652849c42SDave May 
4375c18a9d6SDave May void DataBucketCreateFromSubset( DataBucket DBIn, const PetscInt N, const PetscInt list[], DataBucket *DB )
43852849c42SDave May {
4395c18a9d6SDave May 	PetscInt nfields;
44052849c42SDave May 	DataField *fields;
44152849c42SDave May 	DataBucketCreate(DB);
4425c18a9d6SDave May 	PetscInt f,L,buffer,allocated,p;
44352849c42SDave May 
44452849c42SDave May 	/* copy contents of DBIn */
44552849c42SDave May 	DataBucketGetDataFields(DBIn,&nfields,&fields);
44652849c42SDave May 	DataBucketGetSizes(DBIn,&L,&buffer,&allocated);
44752849c42SDave May 
44852849c42SDave May 	for(f=0;f<nfields;f++) {
44952849c42SDave May 		DataBucketRegisterField(*DB,fields[f]->name,fields[f]->atomic_size,NULL);
45052849c42SDave May 	}
45152849c42SDave May 	DataBucketFinalize(*DB);
45252849c42SDave May 
45352849c42SDave May 	DataBucketSetSizes(*DB,L,buffer);
45452849c42SDave May 
45552849c42SDave May 	/* now copy the desired guys from DBIn => DB */
45652849c42SDave May 	for( p=0; p<N; p++ ) {
45752849c42SDave May 		DataBucketCopyPoint(DBIn,list[p], *DB,p);
45852849c42SDave May 	}
45952849c42SDave May 
46052849c42SDave May }
46152849c42SDave May 
46252849c42SDave May // insert into an exisitng location
4635c18a9d6SDave May void DataFieldInsertPoint( const DataField field, const PetscInt index, const void *ctx )
46452849c42SDave May {
46552849c42SDave May 
46652849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
4675c18a9d6SDave May 	/* check poPetscInt is valid */
46852849c42SDave May 	if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR();  }
46952849c42SDave May 	if( index >= field->L ){ printf("ERROR: index must be < %d\n",field->L); ERROR(); }
47052849c42SDave May #endif
47152849c42SDave May 
47252849c42SDave May //	memcpy( (void*)((char*)field->data + index*field->atomic_size), ctx, field->atomic_size );
47352849c42SDave May 	memcpy( __DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size );
47452849c42SDave May }
47552849c42SDave May 
47652849c42SDave May // remove data at index - replace with last point
4775c18a9d6SDave May void DataBucketRemovePointAtIndex( const DataBucket db, const PetscInt index )
47852849c42SDave May {
4795c18a9d6SDave May 	PetscInt f;
48052849c42SDave May 
48152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
4825c18a9d6SDave May 	/* check poPetscInt is valid */
48352849c42SDave May 	if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); }
48452849c42SDave May 	if( index >= db->allocated ){ printf("ERROR: index must be < %d\n",db->L+db->buffer); ERROR(); }
48552849c42SDave May #endif
48652849c42SDave May 
4875c18a9d6SDave May 	if (index >= db->L) { /* this poPetscInt is not in the list - no need to error, but I will anyway */
4885c18a9d6SDave May 		printf("ERROR: You should not be trying to remove poPetscInt at index=%d since it's < db->L = %d \n", index, db->L );
48952849c42SDave May 		ERROR();
49052849c42SDave May 	}
49152849c42SDave May 
49252849c42SDave May #if 0
4935c18a9d6SDave May 	if (index == db->L-1) { /* last poPetscInt in list */
49452849c42SDave May 		for( f=0; f<db->nfields; f++ ) {
49552849c42SDave May 			DataField field = db->field[f];
49652849c42SDave May 
49752849c42SDave May 			DataFieldZeroPoint(field,index);
49852849c42SDave May 		}
49952849c42SDave May 	}
50052849c42SDave May 	else {
50152849c42SDave May 		for( f=0; f<db->nfields; f++ ) {
50252849c42SDave May 			DataField field = db->field[f];
50352849c42SDave May 
50452849c42SDave May 			/* copy then remove */
50552849c42SDave May 			DataFieldCopyPoint( db->L-1,field, index,field );
50652849c42SDave May 
50752849c42SDave May 			DataFieldZeroPoint(field,index);
50852849c42SDave May 		}
50952849c42SDave May 	}
51052849c42SDave May #endif
51152849c42SDave May 
5125c18a9d6SDave May 	if (index != db->L-1) { /* not last poPetscInt in list */
51352849c42SDave May 		for( f=0; f<db->nfields; f++ ) {
51452849c42SDave May 			DataField field = db->field[f];
51552849c42SDave May 
51652849c42SDave May 			/* copy then remove */
51752849c42SDave May 			DataFieldCopyPoint( db->L-1,field, index,field );
51852849c42SDave May 
51952849c42SDave May 			//DataFieldZeroPoint(field,index);
52052849c42SDave May 		}
52152849c42SDave May 	}
52252849c42SDave May 
52352849c42SDave May 	/* decrement size */
52452849c42SDave May 	/* this will zero out an crap at the end of the list */
52552849c42SDave May 	DataBucketRemovePoint(db);
52652849c42SDave May 
52752849c42SDave May }
52852849c42SDave May 
52952849c42SDave May /* copy x into y */
5305c18a9d6SDave May void DataFieldCopyPoint( const PetscInt pid_x, const DataField field_x,
5315c18a9d6SDave May 												 const PetscInt pid_y, const DataField field_y )
53252849c42SDave May {
53352849c42SDave May 
53452849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
5355c18a9d6SDave May 	/* check poPetscInt is valid */
53652849c42SDave May 	if( pid_x < 0 ){ printf("ERROR: (IN) index must be >= 0\n"); ERROR(); }
53752849c42SDave May 	if( pid_x >= field_x->L ){ printf("ERROR: (IN) index must be < %d\n",field_x->L); ERROR(); }
53852849c42SDave May 
53952849c42SDave May 	if( pid_y < 0 ){ printf("ERROR: (OUT) index must be >= 0\n"); ERROR(); }
54052849c42SDave May 	if( pid_y >= field_y->L ){ printf("ERROR: (OUT) index must be < %d\n",field_y->L); ERROR(); }
54152849c42SDave May 
54252849c42SDave May 	if( field_y->atomic_size != field_x->atomic_size ) {
54352849c42SDave May 		printf("ERROR: atomic size must match \n"); ERROR();
54452849c42SDave May 	}
54552849c42SDave May #endif
54652849c42SDave May 	/*
54752849c42SDave May 	memcpy( (void*)((char*)field_y->data + pid_y*field_y->atomic_size),
54852849c42SDave May 					(void*)((char*)field_x->data + pid_x*field_x->atomic_size),
54952849c42SDave May 				  field_x->atomic_size );
55052849c42SDave May 	*/
55152849c42SDave May 	memcpy(		__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
55252849c42SDave May 						__DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
55352849c42SDave May 						field_y->atomic_size );
55452849c42SDave May 
55552849c42SDave May }
55652849c42SDave May 
55752849c42SDave May 
55852849c42SDave May // zero only the datafield at this point
5595c18a9d6SDave May void DataFieldZeroPoint( const DataField field, const PetscInt index )
56052849c42SDave May {
56152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
5625c18a9d6SDave May 	/* check poPetscInt is valid */
56352849c42SDave May 	if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); }
56452849c42SDave May 	if( index >= field->L ){ printf("ERROR: index must be < %d\n",field->L); ERROR(); }
56552849c42SDave May #endif
56652849c42SDave May 
56752849c42SDave May //	memset( (void*)((char*)field->data + index*field->atomic_size), 0, field->atomic_size );
56852849c42SDave May 	memset( __DATATFIELD_point_access(field->data,index,field->atomic_size), 0, field->atomic_size );
56952849c42SDave May }
57052849c42SDave May 
57152849c42SDave May // zero ALL data for this point
5725c18a9d6SDave May void DataBucketZeroPoint( const DataBucket db, const PetscInt index )
57352849c42SDave May {
5745c18a9d6SDave May 	PetscInt f;
57552849c42SDave May 
5765c18a9d6SDave May 	/* check poPetscInt is valid */
57752849c42SDave May 	if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); }
57852849c42SDave May 	if( index >= db->allocated ){ printf("ERROR: index must be < %d\n",db->allocated); ERROR(); }
57952849c42SDave May 
58052849c42SDave May 	for(f=0;f<db->nfields;f++){
58152849c42SDave May 		DataField field = db->field[f];
58252849c42SDave May 
58352849c42SDave May 		DataFieldZeroPoint(field,index);
58452849c42SDave May 	}
58552849c42SDave May }
58652849c42SDave May 
58752849c42SDave May /* increment */
58852849c42SDave May void DataBucketAddPoint( DataBucket db )
58952849c42SDave May {
59052849c42SDave May 	DataBucketSetSizes( db, db->L+1, -1 );
59152849c42SDave May }
59252849c42SDave May /* decrement */
59352849c42SDave May void DataBucketRemovePoint( DataBucket db )
59452849c42SDave May {
59552849c42SDave May 	DataBucketSetSizes( db, db->L-1, -1 );
59652849c42SDave May }
59752849c42SDave May 
59852849c42SDave May void _DataFieldViewBinary(DataField field, FILE *fp )
59952849c42SDave May {
60052849c42SDave May 	fprintf(fp,"<DataField>\n");
60152849c42SDave May 	fprintf(fp,"%d\n", field->L);
60252849c42SDave May 	fprintf(fp,"%zu\n",field->atomic_size);
60352849c42SDave May 	fprintf(fp,"%s\n", field->registeration_function);
60452849c42SDave May 	fprintf(fp,"%s\n", field->name);
60552849c42SDave May 
60652849c42SDave May 	fwrite(field->data, field->atomic_size, field->L, fp);
60752849c42SDave May /*
60852849c42SDave May 	printf("  ** wrote %zu bytes for DataField \"%s\" \n", field->atomic_size * field->L, field->name );
60952849c42SDave May */
61052849c42SDave May 	fprintf(fp,"\n</DataField>\n");
61152849c42SDave May }
61252849c42SDave May 
61352849c42SDave May void _DataBucketRegisterFieldFromFile( FILE *fp, DataBucket db )
61452849c42SDave May {
61552849c42SDave May 	PetscBool val;
61652849c42SDave May 	DataField *field;
61752849c42SDave May 
61852849c42SDave May 	DataField gfield;
61952849c42SDave May 	char dummy[100];
62052849c42SDave May 	char registeration_function[5000];
62152849c42SDave May 	char field_name[5000];
6225c18a9d6SDave May 	PetscInt L;
62352849c42SDave May 	size_t atomic_size,strL;
62452849c42SDave May 
62552849c42SDave May 
62652849c42SDave May 	/* check we haven't finalised the registration of fields */
62752849c42SDave May 	/*
62852849c42SDave May 	if(db->finalised==PETSC_TRUE) {
62952849c42SDave May 		printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
63052849c42SDave May 		ERROR();
63152849c42SDave May 	}
63252849c42SDave May 	*/
63352849c42SDave May 
63452849c42SDave May 
63552849c42SDave May 	/* read file contents */
63652849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
63752849c42SDave May 
63852849c42SDave May 	fscanf( fp, "%d\n",&L); //printf("read(L): %d\n", L);
63952849c42SDave May 
64052849c42SDave May 	fscanf( fp, "%zu\n",&atomic_size); //printf("read(size): %zu\n",atomic_size);
64152849c42SDave May 
64252849c42SDave May 	fgets(registeration_function,4999,fp); //printf("read(reg func): %s", registeration_function );
64352849c42SDave May 	strL = strlen(registeration_function);
64452849c42SDave May 	if(strL>1){
64552849c42SDave May 		registeration_function[strL-1] = 0;
64652849c42SDave May 	}
64752849c42SDave May 
64852849c42SDave May 	fgets(field_name,4999,fp); //printf("read(name): %s", field_name );
64952849c42SDave May 	strL = strlen(field_name);
65052849c42SDave May 	if(strL>1){
65152849c42SDave May 		field_name[strL-1] = 0;
65252849c42SDave May 	}
65352849c42SDave May 
654*4b46c5e1SDave May #ifdef DATA_BUCKET_LOG
65552849c42SDave May 	printf("  ** read L=%d; atomic_size=%zu; reg_func=\"%s\"; name=\"%s\" \n", L,atomic_size,registeration_function,field_name);
65652849c42SDave May #endif
65752849c42SDave May 
65852849c42SDave May 
65952849c42SDave May 	/* check for repeated name */
66052849c42SDave May 	StringInList( field_name, db->nfields, (const DataField*)db->field, &val );
66152849c42SDave May 	if(val == PETSC_TRUE ) {
66252849c42SDave May 		printf("ERROR: Cannot add same field twice\n");
66352849c42SDave May 		ERROR();
66452849c42SDave May 	}
66552849c42SDave May 
66652849c42SDave May 	/* create new space for data */
66752849c42SDave May 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
66852849c42SDave May 	db->field     = field;
66952849c42SDave May 
67052849c42SDave May 	/* add field */
67152849c42SDave May 	DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield );
67252849c42SDave May 
67352849c42SDave May 	/* copy contents of file */
67452849c42SDave May 	fread(gfield->data, gfield->atomic_size, gfield->L, fp);
675*4b46c5e1SDave May #ifdef DATA_BUCKET_LOG
67652849c42SDave May 	printf("  ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name );
67752849c42SDave May #endif
67852849c42SDave May 	/* finish reading meta data */
67952849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
68052849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
68152849c42SDave May 
68252849c42SDave May 	db->field[ db->nfields ] = gfield;
68352849c42SDave May 
68452849c42SDave May 	db->nfields++;
68552849c42SDave May 
68652849c42SDave May }
68752849c42SDave May 
68852849c42SDave May void _DataBucketViewAscii_HeaderWrite_v00(FILE *fp)
68952849c42SDave May {
69052849c42SDave May 	fprintf(fp,"<DataBucketHeader>\n");
69152849c42SDave May 	fprintf(fp,"type=DataBucket\n");
69252849c42SDave May 	fprintf(fp,"format=ascii\n");
69352849c42SDave May 	fprintf(fp,"version=0.0\n");
69452849c42SDave May 	fprintf(fp,"options=\n");
69552849c42SDave May 	fprintf(fp,"</DataBucketHeader>\n");
69652849c42SDave May }
69752849c42SDave May void _DataBucketViewAscii_HeaderRead_v00(FILE *fp)
69852849c42SDave May {
69952849c42SDave May 	char dummy[100];
70052849c42SDave May 	size_t strL;
70152849c42SDave May 
70252849c42SDave May 	// header open
70352849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
70452849c42SDave May 
70552849c42SDave May 	// type
70652849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
70752849c42SDave May 	strL = strlen(dummy);
70852849c42SDave May 	if(strL>1) { dummy[strL-1] = 0; }
70952849c42SDave May 	if(strcmp(dummy,"type=DataBucket")!=0) {
71052849c42SDave May 		printf("ERROR: Data file doesn't contain a DataBucket type\n");
71152849c42SDave May 		ERROR();
71252849c42SDave May 	}
71352849c42SDave May 
71452849c42SDave May 	// format
71552849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
71652849c42SDave May 
71752849c42SDave May 	// version
71852849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
71952849c42SDave May 	strL = strlen(dummy);
72052849c42SDave May 	if(strL>1) { dummy[strL-1] = 0; }
72152849c42SDave May 	if(strcmp(dummy,"version=0.0")!=0) {
72252849c42SDave May 		printf("ERROR: DataBucket file must be parsed with version=0.0 : You tried %s \n", dummy);
72352849c42SDave May 		ERROR();
72452849c42SDave May 	}
72552849c42SDave May 
72652849c42SDave May 	// options
72752849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
72852849c42SDave May 	// header close
72952849c42SDave May 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
73052849c42SDave May }
73152849c42SDave May 
73252849c42SDave May 
73352849c42SDave May void _DataBucketLoadFromFileBinary_SEQ(const char filename[], DataBucket *_db)
73452849c42SDave May {
73552849c42SDave May 	DataBucket db;
73652849c42SDave May 	FILE *fp;
7375c18a9d6SDave May 	PetscInt L,buffer,f,nfields;
73852849c42SDave May 
73952849c42SDave May 
740*4b46c5e1SDave May #ifdef DATA_BUCKET_LOG
74152849c42SDave May 	printf("** DataBucketLoadFromFile **\n");
74252849c42SDave May #endif
74352849c42SDave May 
74452849c42SDave May 	/* open file */
74552849c42SDave May 	fp = fopen(filename,"rb");
74652849c42SDave May 	if(fp==NULL){
74752849c42SDave May 		printf("ERROR: Cannot open file with name %s \n", filename);
74852849c42SDave May 		ERROR();
74952849c42SDave May 	}
75052849c42SDave May 
75152849c42SDave May 	/* read header */
75252849c42SDave May 	_DataBucketViewAscii_HeaderRead_v00(fp);
75352849c42SDave May 
75452849c42SDave May 	fscanf(fp,"%d\n%d\n%d\n",&L,&buffer,&nfields);
75552849c42SDave May 
75652849c42SDave May 	DataBucketCreate(&db);
75752849c42SDave May 
75852849c42SDave May 	for( f=0; f<nfields; f++ ) {
75952849c42SDave May 		_DataBucketRegisterFieldFromFile(fp,db);
76052849c42SDave May 	}
76152849c42SDave May 	fclose(fp);
76252849c42SDave May 
76352849c42SDave May 	DataBucketFinalize(db);
76452849c42SDave May 
76552849c42SDave May 
76652849c42SDave May /*
76752849c42SDave May   DataBucketSetSizes(db,L,buffer);
76852849c42SDave May */
76952849c42SDave May 	db->L = L;
77052849c42SDave May 	db->buffer = buffer;
77152849c42SDave May 	db->allocated = L + buffer;
77252849c42SDave May 
77352849c42SDave May 	*_db = db;
77452849c42SDave May }
77552849c42SDave May 
77652849c42SDave May void DataBucketLoadFromFile(MPI_Comm comm,const char filename[], DataBucketViewType type, DataBucket *db)
77752849c42SDave May {
7785c18a9d6SDave May 	PetscMPIInt nproc,rank;
77952849c42SDave May 
78052849c42SDave May 	MPI_Comm_size(comm,&nproc);
78152849c42SDave May 	MPI_Comm_rank(comm,&rank);
78252849c42SDave May 
783*4b46c5e1SDave May #ifdef DATA_BUCKET_LOG
78452849c42SDave May 	printf("** DataBucketLoadFromFile **\n");
78552849c42SDave May #endif
78652849c42SDave May 	if(type==DATABUCKET_VIEW_STDOUT) {
78752849c42SDave May 
78852849c42SDave May 	} else if(type==DATABUCKET_VIEW_ASCII) {
78952849c42SDave May 		printf("ERROR: Cannot be implemented as we don't know the underlying particle data structure\n");
79052849c42SDave May 		ERROR();
79152849c42SDave May 	} else if(type==DATABUCKET_VIEW_BINARY) {
79252849c42SDave May 		if (nproc==1) {
79352849c42SDave May 			_DataBucketLoadFromFileBinary_SEQ(filename,db);
79452849c42SDave May 		} else {
79552849c42SDave May 			char *name;
79652849c42SDave May 
79752849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
79852849c42SDave May 			_DataBucketLoadFromFileBinary_SEQ(name,db);
79952849c42SDave May 			free(name);
80052849c42SDave May 		}
80152849c42SDave May 	} else {
80252849c42SDave May 		printf("ERROR: Not implemented\n");
80352849c42SDave May 		ERROR();
80452849c42SDave May 	}
80552849c42SDave May }
80652849c42SDave May 
80752849c42SDave May 
80852849c42SDave May void _DataBucketViewBinary(DataBucket db,const char filename[])
80952849c42SDave May {
81052849c42SDave May 	FILE *fp = NULL;
8115c18a9d6SDave May 	PetscInt f;
81252849c42SDave May 
81352849c42SDave May 	fp = fopen(filename,"wb");
81452849c42SDave May 	if(fp==NULL){
81552849c42SDave May 		printf("ERROR: Cannot open file with name %s \n", filename);
81652849c42SDave May 		ERROR();
81752849c42SDave May 	}
81852849c42SDave May 
81952849c42SDave May 	/* db header */
82052849c42SDave May 	_DataBucketViewAscii_HeaderWrite_v00(fp);
82152849c42SDave May 
82252849c42SDave May 	/* meta-data */
82352849c42SDave May 	fprintf(fp,"%d\n%d\n%d\n", db->L,db->buffer,db->nfields);
82452849c42SDave May 
82552849c42SDave May 	for( f=0; f<db->nfields; f++ ) {
82652849c42SDave May 			/* load datafields */
82752849c42SDave May 		_DataFieldViewBinary(db->field[f],fp);
82852849c42SDave May 	}
82952849c42SDave May 
83052849c42SDave May 	fclose(fp);
83152849c42SDave May }
83252849c42SDave May 
83352849c42SDave May void DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type)
83452849c42SDave May {
83552849c42SDave May 	switch (type) {
83652849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
83752849c42SDave May 		{
8385c18a9d6SDave May 			PetscInt f;
83952849c42SDave May 			double memory_usage_total = 0.0;
84052849c42SDave May 
84152849c42SDave May 			printf("DataBucketView(SEQ): (\"%s\")\n",filename);
84252849c42SDave May 			printf("  L                  = %d \n", db->L );
84352849c42SDave May 			printf("  buffer             = %d \n", db->buffer );
84452849c42SDave May 			printf("  allocated          = %d \n", db->allocated );
84552849c42SDave May 
84652849c42SDave May 			printf("  nfields registered = %d \n", db->nfields );
84752849c42SDave May 			for( f=0; f<db->nfields; f++ ) {
84852849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
84952849c42SDave May 
85052849c42SDave May 				printf("    [%3d]: field name  ==>> %30s : Mem. usage = %1.2e (MB) \n", f, db->field[f]->name, memory_usage_f  );
85152849c42SDave May 				memory_usage_total += memory_usage_f;
85252849c42SDave May 			}
85352849c42SDave May 			printf("  Total mem. usage                                                      = %1.2e (MB) \n", memory_usage_total );
85452849c42SDave May 		}
85552849c42SDave May 			break;
85652849c42SDave May 
85752849c42SDave May 		case DATABUCKET_VIEW_ASCII:
85852849c42SDave May 		{
85952849c42SDave May 			printf("ERROR: Cannot be implemented as we don't know the underlying particle data structure\n");
86052849c42SDave May 			ERROR();
86152849c42SDave May 		}
86252849c42SDave May 			break;
86352849c42SDave May 
86452849c42SDave May 		case DATABUCKET_VIEW_BINARY:
86552849c42SDave May 		{
86652849c42SDave May 			_DataBucketViewBinary(db,filename);
86752849c42SDave May 		}
86852849c42SDave May 			break;
86952849c42SDave May 
87052849c42SDave May 		case DATABUCKET_VIEW_HDF5:
87152849c42SDave May 		{
87252849c42SDave May 			printf("ERROR: Has not been implemented \n");
87352849c42SDave May 			ERROR();
87452849c42SDave May 		}
87552849c42SDave May 			break;
87652849c42SDave May 
87752849c42SDave May 		default:
87852849c42SDave May 			printf("ERROR: Unknown method requested \n");
87952849c42SDave May 			ERROR();
88052849c42SDave May 			break;
88152849c42SDave May 	}
88252849c42SDave May }
88352849c42SDave May 
88452849c42SDave May void DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
88552849c42SDave May {
88652849c42SDave May 	switch (type) {
88752849c42SDave May 		case DATABUCKET_VIEW_STDOUT:
88852849c42SDave May 		{
8895c18a9d6SDave May 			PetscInt f;
8905c18a9d6SDave May 			PetscInt L,buffer,allocated;
89152849c42SDave May 			double memory_usage_total,memory_usage_total_local = 0.0;
8925c18a9d6SDave May 			PetscMPIInt rank;
8935c18a9d6SDave May 			PetscInt ierr;
89452849c42SDave May 
89552849c42SDave May 			ierr = MPI_Comm_rank(comm,&rank);
89652849c42SDave May 
89752849c42SDave May 			DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated);
89852849c42SDave May 
89952849c42SDave May 			for( f=0; f<db->nfields; f++ ) {
90052849c42SDave May 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
90152849c42SDave May 
90252849c42SDave May 				memory_usage_total_local += memory_usage_f;
90352849c42SDave May 			}
90452849c42SDave May 			MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);
90552849c42SDave May 
90652849c42SDave May 			if (rank==0) {
9075c18a9d6SDave May 				PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename);
9085c18a9d6SDave May 				PetscPrintf(comm,"  L                  = %D \n", L );
9095c18a9d6SDave May 				PetscPrintf(comm,"  buffer (max)       = %D \n", buffer );
9105c18a9d6SDave May 				PetscPrintf(comm,"  allocated          = %D \n", allocated );
91152849c42SDave May 
9125c18a9d6SDave May 				PetscPrintf(comm,"  nfields registered = %D \n", db->nfields );
91352849c42SDave May 				for( f=0; f<db->nfields; f++ ) {
91452849c42SDave May 					double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
91552849c42SDave May 
91652849c42SDave May 					printf("    [%3d]: field name  ==>> %30s : Mem. usage = %1.2e (MB) : rank0\n", f, db->field[f]->name, memory_usage_f  );
91752849c42SDave May 				}
91852849c42SDave May 
91952849c42SDave May 				printf("  Total mem. usage                                                      = %1.2e (MB) : collective\n", memory_usage_total );
92052849c42SDave May 			}
92152849c42SDave May 
92252849c42SDave May 		}
92352849c42SDave May 			break;
92452849c42SDave May 
92552849c42SDave May 		case DATABUCKET_VIEW_ASCII:
92652849c42SDave May 		{
92752849c42SDave May 			printf("ERROR: Cannot be implemented as we don't know the underlying particle data structure\n");
92852849c42SDave May 			ERROR();
92952849c42SDave May 		}
93052849c42SDave May 			break;
93152849c42SDave May 
93252849c42SDave May 		case DATABUCKET_VIEW_BINARY:
93352849c42SDave May 		{
93452849c42SDave May 			char *name;
9355c18a9d6SDave May 			PetscMPIInt rank;
93652849c42SDave May 
93752849c42SDave May 			/* create correct extension */
93852849c42SDave May 			MPI_Comm_rank(comm,&rank);
93952849c42SDave May 			asprintf(&name,"%s_p%1.5d",filename, rank );
94052849c42SDave May 
94152849c42SDave May 			_DataBucketViewBinary(db,name);
94252849c42SDave May 
94352849c42SDave May 			free(name);
94452849c42SDave May 		}
94552849c42SDave May 			break;
94652849c42SDave May 
94752849c42SDave May 		case DATABUCKET_VIEW_HDF5:
94852849c42SDave May 		{
94952849c42SDave May 			printf("ERROR: Has not been implemented \n");
95052849c42SDave May 			ERROR();
95152849c42SDave May 		}
95252849c42SDave May 			break;
95352849c42SDave May 
95452849c42SDave May 		default:
95552849c42SDave May 			printf("ERROR: Unknown method requested \n");
95652849c42SDave May 			ERROR();
95752849c42SDave May 			break;
95852849c42SDave May 	}
95952849c42SDave May }
96052849c42SDave May 
96152849c42SDave May 
96252849c42SDave May void DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
96352849c42SDave May {
9645c18a9d6SDave May 	PetscMPIInt nproc;
96552849c42SDave May 
96652849c42SDave May 	MPI_Comm_size(comm,&nproc);
96752849c42SDave May 	if (nproc==1) {
96852849c42SDave May 		DataBucketView_SEQ(db,filename,type);
96952849c42SDave May 	} else {
97052849c42SDave May 		DataBucketView_MPI(comm,db,filename,type);
97152849c42SDave May 	}
97252849c42SDave May }
97352849c42SDave May 
97452849c42SDave May void DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
97552849c42SDave May {
97652849c42SDave May 	DataBucket db2;
9775c18a9d6SDave May 	PetscInt f;
97852849c42SDave May 
97952849c42SDave May 	DataBucketCreate(&db2);
98052849c42SDave May 
98152849c42SDave May 	/* copy contents from dbA into db2 */
98252849c42SDave May 	for (f=0; f<dbA->nfields; f++) {
98352849c42SDave May 		DataField field;
98452849c42SDave May 		size_t    atomic_size;
98552849c42SDave May 		char      *name;
98652849c42SDave May 
98752849c42SDave May 		field = dbA->field[f];
98852849c42SDave May 
98952849c42SDave May 		atomic_size = field->atomic_size;
99052849c42SDave May 		name        = field->name;
99152849c42SDave May 
99252849c42SDave May 		DataBucketRegisterField(db2,name,atomic_size,NULL);
99352849c42SDave May 	}
99452849c42SDave May 	DataBucketFinalize(db2);
99552849c42SDave May 	DataBucketSetInitialSizes(db2,0,1000);
99652849c42SDave May 
99752849c42SDave May 	/* set pointer */
99852849c42SDave May 	*dbB = db2;
99952849c42SDave May }
100052849c42SDave May 
100152849c42SDave May /*
100252849c42SDave May  Insert points from db2 into db1
100352849c42SDave May  db1 <<== db2
100452849c42SDave May  */
100552849c42SDave May void DataBucketInsertValues(DataBucket db1,DataBucket db2)
100652849c42SDave May {
10075c18a9d6SDave May 	PetscInt n_mp_points1,n_mp_points2;
10085c18a9d6SDave May 	PetscInt n_mp_points1_new,p;
100952849c42SDave May 
101052849c42SDave May 	DataBucketGetSizes(db1,&n_mp_points1,0,0);
101152849c42SDave May 	DataBucketGetSizes(db2,&n_mp_points2,0,0);
101252849c42SDave May 
101352849c42SDave May 	n_mp_points1_new = n_mp_points1 + n_mp_points2;
101452849c42SDave May 	DataBucketSetSizes(db1,n_mp_points1_new,-1);
101552849c42SDave May 
101652849c42SDave May 	for (p=0; p<n_mp_points2; p++) {
101752849c42SDave May 		// db1 <<== db2 //
101852849c42SDave May 		DataBucketCopyPoint( db2,p, db1,(n_mp_points1 + p) );
101952849c42SDave May 	}
102052849c42SDave May }
102152849c42SDave May 
102252849c42SDave May /* helpers for parallel send/recv */
102352849c42SDave May void DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
102452849c42SDave May {
10255c18a9d6SDave May     PetscInt       f;
102652849c42SDave May     size_t    sizeof_marker_contents;
102752849c42SDave May     void      *buffer;
102852849c42SDave May 
102952849c42SDave May     sizeof_marker_contents = 0;
103052849c42SDave May     for (f=0; f<db->nfields; f++) {
103152849c42SDave May         DataField df = db->field[f];
103252849c42SDave May 
103352849c42SDave May         sizeof_marker_contents += df->atomic_size;
103452849c42SDave May     }
103552849c42SDave May 
103652849c42SDave May     buffer = malloc(sizeof_marker_contents);
103752849c42SDave May     memset(buffer,0,sizeof_marker_contents);
103852849c42SDave May 
103952849c42SDave May     if (bytes) { *bytes = sizeof_marker_contents; }
104052849c42SDave May     if (buf)   { *buf   = buffer; }
104152849c42SDave May }
104252849c42SDave May 
104352849c42SDave May void DataBucketDestroyPackedArray(DataBucket db,void **buf)
104452849c42SDave May {
104552849c42SDave May     if (buf) {
104652849c42SDave May         free(*buf);
104752849c42SDave May         *buf = NULL;
104852849c42SDave May     }
104952849c42SDave May }
105052849c42SDave May 
10515c18a9d6SDave May void DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
105252849c42SDave May {
10535c18a9d6SDave May     PetscInt    f;
105452849c42SDave May     void   *data,*data_p;
105552849c42SDave May     size_t asize,offset;
105652849c42SDave May 
105752849c42SDave May     offset = 0;
105852849c42SDave May     for (f=0; f<db->nfields; f++) {
105952849c42SDave May         DataField df = db->field[f];
106052849c42SDave May 
106152849c42SDave May         asize = df->atomic_size;
106252849c42SDave May 
106352849c42SDave May         data = (void*)( df->data );
106452849c42SDave May         data_p = (void*)( (char*)data + index*asize );
106552849c42SDave May 
106652849c42SDave May         memcpy( (void*)((char*)buf + offset),  data_p,  asize);
106752849c42SDave May         offset = offset + asize;
106852849c42SDave May     }
106952849c42SDave May }
107052849c42SDave May 
10715c18a9d6SDave May void DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
107252849c42SDave May {
10735c18a9d6SDave May     PetscInt f;
107452849c42SDave May     void *data_p;
107552849c42SDave May     size_t offset;
107652849c42SDave May 
107752849c42SDave May     offset = 0;
107852849c42SDave May     for (f=0; f<db->nfields; f++) {
107952849c42SDave May         DataField df = db->field[f];
108052849c42SDave May 
108152849c42SDave May         data_p = (void*)( (char*)data + offset );
108252849c42SDave May 
108352849c42SDave May         DataFieldInsertPoint(df, idx, (void*)data_p );
108452849c42SDave May         offset = offset + df->atomic_size;
108552849c42SDave May     }
108652849c42SDave May }
1087