xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision 502a2867ed8f1e5015efbc8078522ee2a4bcc95e)
1 
2 #include "data_bucket.h"
3 
4 /* string helpers */
5 PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val)
6 {
7 	PetscInt i;
8 
9 	*val = PETSC_FALSE;
10 	for (i=0; i<N; i++) {
11 		if (strcmp( name, gfield[i]->name) == 0 ) {
12 			*val = PETSC_TRUE;
13       PetscFunctionReturn(0);
14 		}
15 	}
16   PetscFunctionReturn(0);
17 }
18 
19 PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index)
20 {
21 	PetscInt i;
22 
23 	*index = -1;
24 	for (i=0; i<N; i++) {
25 		if (strcmp( name, gfield[i]->name ) == 0) {
26 			*index = i;
27       PetscFunctionReturn(0);
28 		}
29 	}
30   PetscFunctionReturn(0);
31 }
32 
33 #undef __FUNCT__
34 #define __FUNCT__ "DataFieldCreate"
35 PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF)
36 {
37 	DataField df;
38 
39 	df = malloc( sizeof(struct _p_DataField) );
40 	memset( df, 0, sizeof(struct _p_DataField) );
41 
42 
43 	asprintf( &df->registeration_function, "%s", registeration_function );
44 	asprintf( &df->name, "%s", name );
45 	df->atomic_size = size;
46 	df->L = L;
47 
48 	df->data = malloc( size * L ); /* allocate something so we don't have to reallocate */
49 	memset( df->data, 0, size * L );
50 
51 	*DF = df;
52   PetscFunctionReturn(0);
53 }
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "DataFieldDestroy"
57 PetscErrorCode DataFieldDestroy(DataField *DF)
58 {
59 	DataField df = *DF;
60 
61 	free(df->registeration_function);
62 	free(df->name);
63 	free(df->data);
64 	free(df);
65 
66 	*DF = NULL;
67   PetscFunctionReturn(0);
68 }
69 
70 /* data bucket */
71 #undef __FUNCT__
72 #define __FUNCT__ "DataBucketCreate"
73 PetscErrorCode DataBucketCreate(DataBucket *DB)
74 {
75 	DataBucket db;
76 
77 
78 	db = malloc( sizeof(struct _p_DataBucket) );
79 	memset( db, 0, sizeof(struct _p_DataBucket) );
80 
81 	db->finalised = PETSC_FALSE;
82 
83 	/* create empty spaces for fields */
84 	db->L         = -1;
85 	db->buffer    = 1;
86 	db->allocated = 1;
87 
88 	db->nfields   = 0;
89 	db->field     = malloc(sizeof(DataField));
90 
91 	*DB = db;
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "DataBucketDestroy"
97 PetscErrorCode DataBucketDestroy(DataBucket *DB)
98 {
99 	DataBucket db = *DB;
100 	PetscInt f;
101 	PetscErrorCode ierr;
102 
103 	/* release fields */
104 	for (f=0; f<db->nfields; f++) {
105 		ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
106 	}
107 
108 	/* this will catch the initially allocated objects in the event that no fields are registered */
109 	if (db->field != NULL) {
110 		free(db->field);
111 	}
112 
113 	free(db);
114 
115 	*DB = NULL;
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "DataBucketRegisterField"
121 PetscErrorCode DataBucketRegisterField(
122                               DataBucket db,
123                               const char registeration_function[],
124                               const char field_name[],
125                               size_t atomic_size, DataField *_gfield)
126 {
127 	PetscBool val;
128 	DataField *field,fp;
129 	PetscErrorCode ierr;
130 
131 	/* check we haven't finalised the registration of fields */
132 	/*
133    if(db->finalised==PETSC_TRUE) {
134    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
135    ERROR();
136    }
137    */
138 
139 	/* check for repeated name */
140 	StringInList( field_name, db->nfields, (const DataField*)db->field, &val );
141 	if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
142 
143 	/* create new space for data */
144 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
145 	db->field     = field;
146 
147 	/* add field */
148 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, db->allocated, &fp );CHKERRQ(ierr);
149 	db->field[ db->nfields ] = fp;
150 
151 	db->nfields++;
152 
153 	if (_gfield != NULL) {
154 		*_gfield = fp;
155 	}
156   PetscFunctionReturn(0);
157 }
158 
159 /*
160  #define DataBucketRegisterField(db,name,size,k) {\
161  char *location;\
162  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
163  _DataBucketRegisterField( (db), location, (name), (size), (k) );\
164  free(location);\
165  }
166  */
167 
168 #undef __FUNCT__
169 #define __FUNCT__ "DataBucketGetDataFieldByName"
170 PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
171 {
172 	PetscInt idx;
173 	PetscBool found;
174 
175 	StringInList(name,db->nfields,(const DataField*)db->field,&found);
176 	if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name);
177 
178 	StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);
179 
180 	*gfield = db->field[idx];
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "DataBucketQueryDataFieldByName"
186 PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
187 {
188 	*found = PETSC_FALSE;
189 	StringInList(name,db->nfields,(const DataField*)db->field,found);
190   PetscFunctionReturn(0);
191 }
192 
193 #undef __FUNCT__
194 #define __FUNCT__ "DataBucketFinalize"
195 PetscErrorCode DataBucketFinalize(DataBucket db)
196 {
197 	db->finalised = PETSC_TRUE;
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "DataFieldGetNumEntries"
203 PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum)
204 {
205 	*sum = df->L;
206   PetscFunctionReturn(0);
207 }
208 
209 #undef __FUNCT__
210 #define __FUNCT__ "DataFieldSetSize"
211 PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L)
212 {
213 	void *tmp_data;
214 
215 	if (new_L <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be <= 0");
216 
217 	if (new_L == df->L) PetscFunctionReturn(0);
218 
219 	if (new_L > df->L) {
220 
221 		tmp_data = realloc( df->data, df->atomic_size * (new_L) );
222 		df->data = tmp_data;
223 
224 		/* init new contents */
225 		memset( ( ((char*)df->data)+df->L*df->atomic_size), 0, (new_L-df->L)*df->atomic_size );
226 
227 	} else {
228 		/* reallocate pointer list, add +1 in case new_L = 0 */
229 		tmp_data = realloc( df->data, df->atomic_size * (new_L+1) );
230 		df->data = tmp_data;
231 	}
232 
233 	df->L = new_L;
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "DataFieldZeroBlock"
239 PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end)
240 {
241 	if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
242 
243 	if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
244 
245 	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);
246 
247 	memset( ( ((char*)df->data)+start*df->atomic_size), 0, (end-start)*df->atomic_size );
248   PetscFunctionReturn(0);
249 }
250 
251 /*
252  A negative buffer value will simply be ignored and the old buffer value will be used.
253  */
254 #undef __FUNCT__
255 #define __FUNCT__ "DataBucketSetSizes"
256 PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
257 {
258 	PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
259   PetscErrorCode ierr;
260 
261 	if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()");
262 
263 	current_allocated = db->allocated;
264 
265 	new_used   = L;
266 	new_unused = current_allocated - new_used;
267 	new_buffer = db->buffer;
268 	if (buffer >= 0) { /* update the buffer value */
269 		new_buffer = buffer;
270 	}
271 	new_allocated = new_used + new_buffer;
272 
273 	/* action */
274 	if (new_allocated > current_allocated) {
275 		/* increase size to new_used + new_buffer */
276 		for (f=0; f<db->nfields; f++) {
277 			ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
278 		}
279 
280 		db->L         = new_used;
281 		db->buffer    = new_buffer;
282 		db->allocated = new_used + new_buffer;
283 	}
284 	else {
285 		if (new_unused > 2 * new_buffer) {
286 
287 			/* shrink array to new_used + new_buffer */
288 			for (f=0; f<db->nfields; f++) {
289 				ierr = DataFieldSetSize( db->field[f], new_allocated );CHKERRQ(ierr);
290 			}
291 
292 			db->L         = new_used;
293 			db->buffer    = new_buffer;
294 			db->allocated = new_used + new_buffer;
295 		}
296 		else {
297 			db->L      = new_used;
298 			db->buffer = new_buffer;
299 		}
300 	}
301 
302 	/* zero all entries from db->L to db->allocated */
303 	for (f=0; f<db->nfields; f++) {
304 		DataField field = db->field[f];
305 		ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
306 	}
307   PetscFunctionReturn(0);
308 }
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "DataBucketSetInitialSizes"
312 PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
313 {
314 	PetscInt f;
315   PetscErrorCode ierr;
316 
317 	ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
318 
319 	for (f=0; f<db->nfields; f++) {
320 		DataField field = db->field[f];
321 		ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
322 	}
323   PetscFunctionReturn(0);
324 }
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "DataBucketGetSizes"
328 PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
329 {
330 	if (L) { *L = db->L; }
331 	if (buffer) { *buffer = db->buffer; }
332 	if (allocated) { *allocated = db->allocated; }
333   PetscFunctionReturn(0);
334 }
335 
336 #undef __FUNCT__
337 #define __FUNCT__ "DataBucketGetGlobalSizes"
338 PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
339 {
340 	PetscInt _L,_buffer,_allocated;
341 	PetscInt ierr;
342 
343 	_L = db->L;
344 	_buffer = db->buffer;
345 	_allocated = db->allocated;
346 
347 	if (L) {         ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm); }
348 	if (buffer) {    ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm); }
349 	if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm); }
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNCT__
354 #define __FUNCT__ "DataBucketGetGlobalSizes"
355 PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[])
356 {
357 	if (L) {      *L      = db->nfields; }
358 	if (fields) { *fields = db->field; }
359   PetscFunctionReturn(0);
360 }
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "DataFieldGetAccess"
364 PetscErrorCode DataFieldGetAccess(const DataField gfield)
365 {
366 	if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name);
367 
368 	gfield->active = PETSC_TRUE;
369   PetscFunctionReturn(0);
370 }
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "DataFieldAccessPoint"
374 PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p)
375 {
376 #ifdef DATAFIELD_POINT_ACCESS_GUARD
377 	/* debug mode */
378 	/* check point is valid */
379 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
380 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
381 
382 	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);
383 #endif
384 
385 	//*ctx_p  = (void*)( ((char*)gfield->data) + pid * gfield->atomic_size);
386 	*ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
387   PetscFunctionReturn(0);
388 }
389 
390 #undef __FUNCT__
391 #define __FUNCT__ "DataFieldAccessPointOffset"
392 PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
393 {
394 #ifdef DATAFIELD_POINT_ACCESS_GUARD
395 	/* debug mode */
396 
397 	/* check point is valid */
398 	/* 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 */
399 	if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
400 
401 	/* check point is valid */
402 	if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
403 	if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
404 
405 	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);
406 #endif
407 
408 	*ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
409   PetscFunctionReturn(0);
410 }
411 
412 #undef __FUNCT__
413 #define __FUNCT__ "DataFieldAccessPointOffset"
414 PetscErrorCode DataFieldRestoreAccess(DataField gfield)
415 {
416 	if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name );
417 
418 	gfield->active = PETSC_FALSE;
419   PetscFunctionReturn(0);
420 }
421 
422 #undef __FUNCT__
423 #define __FUNCT__ "DataFieldVerifyAccess"
424 PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size)
425 {
426 #ifdef DATAFIELD_POINT_ACCESS_GUARD
427 	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.",
428            gfield->name, gfield->atomic_size, size );
429 #endif
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "DataFieldGetAtomicSize"
435 PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size)
436 {
437   if (size) { *size = gfield->atomic_size; }
438   PetscFunctionReturn(0);
439 }
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "DataFieldGetEntries"
443 PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data)
444 {
445   if (data) {
446     *data = gfield->data;
447   }
448   PetscFunctionReturn(0);
449 }
450 
451 #undef __FUNCT__
452 #define __FUNCT__ "DataFieldRestoreEntries"
453 PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data)
454 {
455   if (data) {
456     *data = NULL;
457   }
458   PetscFunctionReturn(0);
459 }
460 
461 /* y = x */
462 #undef __FUNCT__
463 #define __FUNCT__ "DataBucketCopyPoint"
464 PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,
465                          const DataBucket yb,const PetscInt pid_y)
466 {
467 	PetscInt f;
468 	PetscErrorCode ierr;
469 
470 	for (f=0; f<xb->nfields; f++) {
471 		void *dest;
472 		void *src;
473 
474 		ierr = DataFieldGetAccess( xb->field[f] );CHKERRQ(ierr);
475 		if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f] );CHKERRQ(ierr); }
476 
477 		ierr = DataFieldAccessPoint( xb->field[f],pid_x, &src );CHKERRQ(ierr);
478 		ierr = DataFieldAccessPoint( yb->field[f],pid_y, &dest );CHKERRQ(ierr);
479 
480 		memcpy( dest, src, xb->field[f]->atomic_size );
481 
482 		ierr = DataFieldRestoreAccess( xb->field[f] );CHKERRQ(ierr);
483 		if (xb != yb) { ierr = DataFieldRestoreAccess( yb->field[f] );CHKERRQ(ierr); }
484 	}
485   PetscFunctionReturn(0);
486 }
487 
488 #undef __FUNCT__
489 #define __FUNCT__ "DataBucketCreateFromSubset"
490 PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB)
491 {
492 	PetscInt nfields;
493 	DataField *fields;
494 	DataBucketCreate(DB);
495 	PetscInt f,L,buffer,allocated,p;
496 	PetscErrorCode ierr;
497 
498 	/* copy contents of DBIn */
499 	ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
500 	ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
501 
502 	for (f=0; f<nfields; f++) {
503 		ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
504 	}
505 	ierr = DataBucketFinalize(*DB);CHKERRQ(ierr);
506 
507 	ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
508 
509 	/* now copy the desired guys from DBIn => DB */
510 	for (p=0; p<N; p++) {
511 		ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
512 	}
513   PetscFunctionReturn(0);
514 }
515 
516 // insert into an exisitng location
517 #undef __FUNCT__
518 #define __FUNCT__ "DataFieldInsertPoint"
519 PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx)
520 {
521 
522 #ifdef DATAFIELD_POINT_ACCESS_GUARD
523 	/* check point is valid */
524 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
525 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
526 #endif
527 
528   //	memcpy( (void*)((char*)field->data + index*field->atomic_size), ctx, field->atomic_size );
529 	memcpy( __DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size );
530   PetscFunctionReturn(0);
531 }
532 
533 // remove data at index - replace with last point
534 #undef __FUNCT__
535 #define __FUNCT__ "DataBucketRemovePointAtIndex"
536 PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index)
537 {
538 	PetscInt f;
539 	PetscErrorCode ierr;
540 
541 #ifdef DATAFIELD_POINT_ACCESS_GUARD
542 	/* check point is valid */
543 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
544 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
545 #endif
546 
547 	if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
548 		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 );
549 	}
550 
551 	if (index != db->L-1) { /* not last point in list */
552 		for (f=0; f<db->nfields; f++) {
553 			DataField field = db->field[f];
554 
555 			/* copy then remove */
556 			ierr = DataFieldCopyPoint( db->L-1,field, index,field );CHKERRQ(ierr);
557 
558 			//DataFieldZeroPoint(field,index);
559 		}
560 	}
561 
562 	/* decrement size */
563 	/* this will zero out an crap at the end of the list */
564 	ierr = DataBucketRemovePoint(db);CHKERRQ(ierr);
565   PetscFunctionReturn(0);
566 }
567 
568 /* copy x into y */
569 #undef __FUNCT__
570 #define __FUNCT__ "DataFieldCopyPoint"
571 PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,
572                         const PetscInt pid_y,const DataField field_y )
573 {
574 
575 #ifdef DATAFIELD_POINT_ACCESS_GUARD
576 	/* check point is valid */
577 	if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
578 	if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
579 
580 	if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
581 	if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
582 
583 	if( field_y->atomic_size != field_x->atomic_size ) {
584 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
585 	}
586 #endif
587 	/*
588    memcpy( (void*)((char*)field_y->data + pid_y*field_y->atomic_size),
589    (void*)((char*)field_x->data + pid_x*field_x->atomic_size),
590    field_x->atomic_size );
591    */
592 	memcpy(		__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
593          __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
594          field_y->atomic_size );
595   PetscFunctionReturn(0);
596 }
597 
598 
599 // zero only the datafield at this point
600 #undef __FUNCT__
601 #define __FUNCT__ "DataFieldZeroPoint"
602 PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index)
603 {
604 #ifdef DATAFIELD_POINT_ACCESS_GUARD
605 	/* check point is valid */
606 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
607 	if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
608 #endif
609 
610   //	memset( (void*)((char*)field->data + index*field->atomic_size), 0, field->atomic_size );
611 	memset( __DATATFIELD_point_access(field->data,index,field->atomic_size), 0, field->atomic_size );
612   PetscFunctionReturn(0);
613 }
614 
615 // zero ALL data for this point
616 #undef __FUNCT__
617 #define __FUNCT__ "DataBucketZeroPoint"
618 PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index)
619 {
620 	PetscInt f;
621 	PetscErrorCode ierr;
622 
623 	/* check point is valid */
624 	if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
625 	if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
626 
627 	for (f=0; f<db->nfields; f++) {
628 		DataField field = db->field[f];
629 
630 		ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr);
631 	}
632   PetscFunctionReturn(0);
633 }
634 
635 /* increment */
636 #undef __FUNCT__
637 #define __FUNCT__ "DataBucketAddPoint"
638 PetscErrorCode DataBucketAddPoint(DataBucket db)
639 {
640 	PetscErrorCode ierr;
641 
642 	ierr = DataBucketSetSizes(db,db->L+1,-1);CHKERRQ(ierr);
643   PetscFunctionReturn(0);
644 }
645 
646 /* decrement */
647 #undef __FUNCT__
648 #define __FUNCT__ "DataBucketRemovePoint"
649 PetscErrorCode DataBucketRemovePoint(DataBucket db)
650 {
651 	PetscErrorCode ierr;
652 
653 	ierr = DataBucketSetSizes(db,db->L-1,-1);CHKERRQ(ierr);
654   PetscFunctionReturn(0);
655 }
656 
657 #undef __FUNCT__
658 #define __FUNCT__ "_DataFieldViewBinary"
659 PetscErrorCode _DataFieldViewBinary(DataField field,FILE *fp)
660 {
661 	fprintf(fp,"<DataField>\n");
662 	fprintf(fp,"%d\n", field->L);
663 	fprintf(fp,"%zu\n",field->atomic_size);
664 	fprintf(fp,"%s\n", field->registeration_function);
665 	fprintf(fp,"%s\n", field->name);
666 
667 	fwrite(field->data, field->atomic_size, field->L, fp);
668   /*
669    printf("  ** wrote %zu bytes for DataField \"%s\" \n", field->atomic_size * field->L, field->name );
670    */
671 	fprintf(fp,"\n</DataField>\n");
672   PetscFunctionReturn(0);
673 }
674 
675 #undef __FUNCT__
676 #define __FUNCT__ "_DataBucketRegisterFieldFromFile"
677 PetscErrorCode _DataBucketRegisterFieldFromFile(FILE *fp,DataBucket db)
678 {
679 	PetscBool val;
680 	DataField *field;
681 
682 	DataField gfield;
683 	char dummy[100];
684 	char registeration_function[5000];
685 	char field_name[5000];
686 	PetscInt L;
687 	size_t atomic_size,strL;
688 	PetscErrorCode ierr;
689 
690 	/* check we haven't finalised the registration of fields */
691 	/*
692    if(db->finalised==PETSC_TRUE) {
693    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
694    ERROR();
695    }
696    */
697 
698 
699 	/* read file contents */
700 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
701 
702 	fscanf( fp, "%d\n",&L); //printf("read(L): %d\n", L);
703 
704 	fscanf( fp, "%zu\n",&atomic_size); //printf("read(size): %zu\n",atomic_size);
705 
706 	fgets(registeration_function,4999,fp); //printf("read(reg func): %s", registeration_function );
707 	strL = strlen(registeration_function);
708 	if (strL > 1) {
709 		registeration_function[strL-1] = 0;
710 	}
711 
712 	fgets(field_name,4999,fp); //printf("read(name): %s", field_name );
713 	strL = strlen(field_name);
714 	if (strL > 1) {
715 		field_name[strL-1] = 0;
716 	}
717 
718 #ifdef DATA_BUCKET_LOG
719 	PetscPrintf(PETSC_COMM_SELF,"  ** read L=%D; atomic_size=%zu; reg_func=\"%s\"; name=\"%s\" \n", L,atomic_size,registeration_function,field_name);
720 #endif
721 
722 
723 	/* check for repeated name */
724 	StringInList( field_name, db->nfields, (const DataField*)db->field, &val );
725 	if (val == PETSC_TRUE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot add same field twice");
726 
727 	/* create new space for data */
728 	field = realloc( db->field,     sizeof(DataField)*(db->nfields+1));
729 	db->field     = field;
730 
731 	/* add field */
732 	ierr = DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield );CHKERRQ(ierr);
733 
734 	/* copy contents of file */
735 	fread(gfield->data, gfield->atomic_size, gfield->L, fp);
736 #ifdef DATA_BUCKET_LOG
737 	PetscPrintf(PETSC_COMM_SELF,"  ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name );
738 #endif
739 	/* finish reading meta data */
740 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
741 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
742 
743 	db->field[ db->nfields ] = gfield;
744 
745 	db->nfields++;
746   PetscFunctionReturn(0);
747 }
748 
749 #undef __FUNCT__
750 #define __FUNCT__ "_DataBucketViewAscii_HeaderWrite_v00"
751 PetscErrorCode _DataBucketViewAscii_HeaderWrite_v00(FILE *fp)
752 {
753 	fprintf(fp,"<DataBucketHeader>\n");
754 	fprintf(fp,"type=DataBucket\n");
755 	fprintf(fp,"format=ascii\n");
756 	fprintf(fp,"version=0.0\n");
757 	fprintf(fp,"options=\n");
758 	fprintf(fp,"</DataBucketHeader>\n");
759   PetscFunctionReturn(0);
760 }
761 
762 #undef __FUNCT__
763 #define __FUNCT__ "_DataBucketViewAscii_HeaderRead_v00"
764 PetscErrorCode _DataBucketViewAscii_HeaderRead_v00(FILE *fp)
765 {
766 	char dummy[100];
767 	size_t strL;
768 
769 	// header open
770 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
771 
772 	// type
773 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
774 	strL = strlen(dummy);
775 	if (strL > 1) { dummy[strL-1] = 0; }
776 	if (strcmp(dummy,"type=DataBucket") != 0) {
777 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Data file doesn't contain a DataBucket type");
778 	}
779 
780 	// format
781 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
782 
783 	// version
784 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
785 	strL = strlen(dummy);
786 	if (strL > 1) { dummy[strL-1] = 0; }
787 	if (strcmp(dummy,"version=0.0") != 0) {
788 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"DataBucket file must be parsed with version=0.0 : You tried %s", dummy);
789 	}
790 
791 	// options
792 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
793 	// header close
794 	fgets(dummy,99,fp); //printf("read(header): %s", dummy );
795   PetscFunctionReturn(0);
796 }
797 
798 #undef __FUNCT__
799 #define __FUNCT__ "_DataBucketLoadFromFileBinary_SEQ"
800 PetscErrorCode _DataBucketLoadFromFileBinary_SEQ(const char filename[],DataBucket *_db)
801 {
802 	DataBucket db;
803 	FILE *fp;
804 	PetscInt L,buffer,f,nfields;
805 	PetscErrorCode ierr;
806 
807 
808 #ifdef DATA_BUCKET_LOG
809 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
810 #endif
811 
812 	/* open file */
813 	fp = fopen(filename,"rb");
814 	if (fp == NULL) {
815 		SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file with name %s", filename);
816 	}
817 
818 	/* read header */
819 	ierr = _DataBucketViewAscii_HeaderRead_v00(fp);CHKERRQ(ierr);
820 
821 	fscanf(fp,"%d\n%d\n%d\n",&L,&buffer,&nfields);
822 
823 	ierr = DataBucketCreate(&db);CHKERRQ(ierr);
824 
825 	for (f=0; f<nfields; f++) {
826 		ierr = _DataBucketRegisterFieldFromFile(fp,db);CHKERRQ(ierr);
827 	}
828 	fclose(fp);
829 
830 	ierr = DataBucketFinalize(db);CHKERRQ(ierr);
831 
832   /*
833    DataBucketSetSizes(db,L,buffer);
834    */
835 	db->L = L;
836 	db->buffer = buffer;
837 	db->allocated = L + buffer;
838 
839 	*_db = db;
840   PetscFunctionReturn(0);
841 }
842 
843 #undef __FUNCT__
844 #define __FUNCT__ "DataBucketLoadFromFile"
845 PetscErrorCode DataBucketLoadFromFile(MPI_Comm comm,const char filename[],DataBucketViewType type,DataBucket *db)
846 {
847 	PetscMPIInt nproc,rank;
848 	PetscErrorCode ierr;
849 
850 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
851 	ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
852 
853 #ifdef DATA_BUCKET_LOG
854 	PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");
855 #endif
856 	if (type == DATABUCKET_VIEW_STDOUT) {
857 
858 	} else if (type == DATABUCKET_VIEW_ASCII) {
859 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
860 	} else if (type == DATABUCKET_VIEW_BINARY) {
861 		if (nproc == 1) {
862 			ierr = _DataBucketLoadFromFileBinary_SEQ(filename,db);CHKERRQ(ierr);
863 		} else {
864 			char *name;
865 
866 			asprintf(&name,"%s_p%1.5d",filename, rank );
867 			ierr = _DataBucketLoadFromFileBinary_SEQ(name,db);CHKERRQ(ierr);
868 			free(name);
869 		}
870 	} else {
871 		SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer requested");
872 	}
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "_DataBucketViewBinary"
878 PetscErrorCode _DataBucketViewBinary(DataBucket db,const char filename[])
879 {
880 	FILE *fp = NULL;
881 	PetscInt f;
882 	PetscErrorCode ierr;
883 
884 	fp = fopen(filename,"wb");
885 	if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file with name %s", filename);
886 
887 	/* db header */
888 	ierr =_DataBucketViewAscii_HeaderWrite_v00(fp);CHKERRQ(ierr);
889 
890 	/* meta-data */
891 	fprintf(fp,"%d\n%d\n%d\n", db->L,db->buffer,db->nfields);
892 
893 	for (f=0; f<db->nfields; f++) {
894     /* load datafields */
895 		ierr = _DataFieldViewBinary(db->field[f],fp);CHKERRQ(ierr);
896 	}
897 
898 	fclose(fp);
899   PetscFunctionReturn(0);
900 }
901 
902 #undef __FUNCT__
903 #define __FUNCT__ "DataBucketView_SEQ"
904 PetscErrorCode DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type)
905 {
906   PetscErrorCode ierr;
907 
908 	switch (type) {
909 		case DATABUCKET_VIEW_STDOUT:
910 		{
911 			PetscInt f;
912 			double memory_usage_total = 0.0;
913 
914 			PetscPrintf(PETSC_COMM_SELF,"DataBucketView(SEQ): (\"%s\")\n",filename);
915 			PetscPrintf(PETSC_COMM_SELF,"  L                  = %D \n", db->L );
916 			PetscPrintf(PETSC_COMM_SELF,"  buffer             = %D \n", db->buffer );
917 			PetscPrintf(PETSC_COMM_SELF,"  allocated          = %D \n", db->allocated );
918 
919 			PetscPrintf(PETSC_COMM_SELF,"  nfields registered = %D \n", db->nfields );
920 			for (f=0; f<db->nfields; f++) {
921 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
922 
923 				PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) \n", f, db->field[f]->name, memory_usage_f  );
924 				memory_usage_total += memory_usage_f;
925 			}
926 			PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) \n", memory_usage_total );
927 		}
928 			break;
929 
930 		case DATABUCKET_VIEW_ASCII:
931 		{
932 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
933 		}
934 			break;
935 
936 		case DATABUCKET_VIEW_BINARY:
937 		{
938 			ierr = _DataBucketViewBinary(db,filename);CHKERRQ(ierr);
939 		}
940 			break;
941 
942 		case DATABUCKET_VIEW_HDF5:
943 		{
944 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No HDF5 support");
945 		}
946 			break;
947 
948 		default:
949 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
950 			break;
951 	}
952   PetscFunctionReturn(0);
953 }
954 
955 #undef __FUNCT__
956 #define __FUNCT__ "DataBucketView_MPI"
957 PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
958 {
959   PetscErrorCode ierr;
960 
961 	switch (type) {
962 		case DATABUCKET_VIEW_STDOUT:
963 		{
964 			PetscInt f;
965 			PetscInt L,buffer,allocated;
966 			double memory_usage_total,memory_usage_total_local = 0.0;
967 			PetscMPIInt rank;
968 
969 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
970 
971 			ierr = DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated);CHKERRQ(ierr);
972 
973 			for (f=0; f<db->nfields; f++) {
974 				double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
975 
976 				memory_usage_total_local += memory_usage_f;
977 			}
978 			ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
979 
980 			if (rank == 0) {
981 				PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename);
982 				PetscPrintf(comm,"  L                  = %D \n", L );
983 				PetscPrintf(comm,"  buffer (max)       = %D \n", buffer );
984 				PetscPrintf(comm,"  allocated          = %D \n", allocated );
985 
986 				PetscPrintf(comm,"  nfields registered = %D \n", db->nfields );
987 				for (f=0; f<db->nfields; f++) {
988 					double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
989 
990 					PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) : rank0\n", f, db->field[f]->name, memory_usage_f  );
991 				}
992 
993 				PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) : collective\n", memory_usage_total );
994 			}
995 
996 		}
997 			break;
998 
999 		case DATABUCKET_VIEW_ASCII:
1000 		{
1001 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying data structure");
1002 		}
1003 			break;
1004 
1005 		case DATABUCKET_VIEW_BINARY:
1006 		{
1007 			char *name;
1008 			PetscMPIInt rank;
1009 
1010 			/* create correct extension */
1011 			ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1012 			asprintf(&name,"%s_p%1.5d",filename, rank );
1013 
1014 			ierr = _DataBucketViewBinary(db,name);CHKERRQ(ierr);
1015 
1016 			free(name);
1017 		}
1018 			break;
1019 
1020 		case DATABUCKET_VIEW_HDF5:
1021 		{
1022 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5");
1023 		}
1024 			break;
1025 
1026 		default:
1027 			SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
1028 			break;
1029 	}
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 #undef __FUNCT__
1034 #define __FUNCT__ "DataBucketView"
1035 PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
1036 {
1037 	PetscMPIInt nproc;
1038   PetscErrorCode ierr;
1039 
1040 	ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
1041 	if (nproc == 1) {
1042 		ierr =DataBucketView_SEQ(db,filename,type);CHKERRQ(ierr);
1043 	} else {
1044 		ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
1045 	}
1046   PetscFunctionReturn(0);
1047 }
1048 
1049 #undef __FUNCT__
1050 #define __FUNCT__ "DataBucketDuplicateFields"
1051 PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
1052 {
1053 	DataBucket db2;
1054 	PetscInt f;
1055 	PetscErrorCode ierr;
1056 
1057 	ierr = DataBucketCreate(&db2);CHKERRQ(ierr);
1058 
1059 	/* copy contents from dbA into db2 */
1060 	for (f=0; f<dbA->nfields; f++) {
1061 		DataField field;
1062 		size_t    atomic_size;
1063 		char      *name;
1064 
1065 		field = dbA->field[f];
1066 
1067 		atomic_size = field->atomic_size;
1068 		name        = field->name;
1069 
1070 		ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
1071 	}
1072 	ierr = DataBucketFinalize(db2);CHKERRQ(ierr);
1073 	ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
1074 
1075 	/* set pointer */
1076 	*dbB = db2;
1077   PetscFunctionReturn(0);
1078 }
1079 
1080 /*
1081  Insert points from db2 into db1
1082  db1 <<== db2
1083  */
1084 #undef __FUNCT__
1085 #define __FUNCT__ "DataBucketInsertValues"
1086 PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2)
1087 {
1088 	PetscInt n_mp_points1,n_mp_points2;
1089 	PetscInt n_mp_points1_new,p;
1090 	PetscErrorCode ierr;
1091 
1092 	ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr);
1093 	ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr);
1094 
1095 	n_mp_points1_new = n_mp_points1 + n_mp_points2;
1096 	ierr = DataBucketSetSizes(db1,n_mp_points1_new,-1);CHKERRQ(ierr);
1097 
1098 	for (p=0; p<n_mp_points2; p++) {
1099 		// db1 <<== db2 //
1100 		ierr =DataBucketCopyPoint( db2,p, db1,(n_mp_points1 + p) );CHKERRQ(ierr);
1101 	}
1102   PetscFunctionReturn(0);
1103 }
1104 
1105 /* helpers for parallel send/recv */
1106 #undef __FUNCT__
1107 #define __FUNCT__ "DataBucketCreatePackedArray"
1108 PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
1109 {
1110   PetscInt       f;
1111   size_t    sizeof_marker_contents;
1112   void      *buffer;
1113 
1114   sizeof_marker_contents = 0;
1115   for (f=0; f<db->nfields; f++) {
1116     DataField df = db->field[f];
1117 
1118     sizeof_marker_contents += df->atomic_size;
1119   }
1120 
1121   buffer = malloc(sizeof_marker_contents);
1122   memset(buffer,0,sizeof_marker_contents);
1123 
1124   if (bytes) { *bytes = sizeof_marker_contents; }
1125   if (buf)   { *buf   = buffer; }
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 #undef __FUNCT__
1130 #define __FUNCT__ "DataBucketDestroyPackedArray"
1131 PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf)
1132 {
1133   if (buf) {
1134     free(*buf);
1135     *buf = NULL;
1136   }
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 #undef __FUNCT__
1141 #define __FUNCT__ "DataBucketFillPackedArray"
1142 PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
1143 {
1144   PetscInt    f;
1145   void   *data,*data_p;
1146   size_t asize,offset;
1147 
1148   offset = 0;
1149   for (f=0; f<db->nfields; f++) {
1150     DataField df = db->field[f];
1151 
1152     asize = df->atomic_size;
1153 
1154     data = (void*)( df->data );
1155     data_p = (void*)( (char*)data + index*asize );
1156 
1157     memcpy( (void*)((char*)buf + offset),  data_p,  asize);
1158     offset = offset + asize;
1159   }
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 #undef __FUNCT__
1164 #define __FUNCT__ "DataBucketInsertPackedArray"
1165 PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
1166 {
1167   PetscInt f;
1168   void *data_p;
1169   size_t offset;
1170 	PetscErrorCode ierr;
1171 
1172   offset = 0;
1173   for (f=0; f<db->nfields; f++) {
1174     DataField df = db->field[f];
1175 
1176     data_p = (void*)( (char*)data + offset );
1177 
1178     ierr = DataFieldInsertPoint(df, idx, (void*)data_p );CHKERRQ(ierr);
1179     offset = offset + df->atomic_size;
1180   }
1181   PetscFunctionReturn(0);
1182 }
1183