xref: /petsc/src/sys/utils/segbuffer.c (revision 47452e7be4b12ea2a938340a38eb49cf6e2cf13a)
10f453b92SJed Brown #include <petscsys.h>
20f453b92SJed Brown 
30f453b92SJed Brown /* Segmented (extendable) array implementation */
40f453b92SJed Brown struct _n_PetscSegBuffer {
50f453b92SJed Brown   PetscInt unitbytes;
60f453b92SJed Brown   PetscInt alloc;
70f453b92SJed Brown   PetscInt used;
80f453b92SJed Brown   PetscInt tailused;
90f453b92SJed Brown   PetscSegBuffer tail;
100f453b92SJed Brown   union {                       /* Dummy types to ensure alignment */
110f453b92SJed Brown     PetscReal dummy_real;
120f453b92SJed Brown     PetscInt  dummy_int;
130f453b92SJed Brown     char      array[1];
140f453b92SJed Brown   } u;
150f453b92SJed Brown };
160f453b92SJed Brown 
170f453b92SJed Brown #undef __FUNCT__
180f453b92SJed Brown #define __FUNCT__ "PetscSegBufferAlloc_Private"
190f453b92SJed Brown static PetscErrorCode PetscSegBufferAlloc_Private(PetscSegBuffer *seg,PetscInt count)
200f453b92SJed Brown {
210f453b92SJed Brown   PetscErrorCode ierr;
220f453b92SJed Brown   PetscSegBuffer newseg,s;
230f453b92SJed Brown   PetscInt       alloc;
240f453b92SJed Brown 
250f453b92SJed Brown   PetscFunctionBegin;
260f453b92SJed Brown   s = *seg;
270f453b92SJed Brown   /* Grow at least fast enough to hold next item, like Fibonacci otherwise (up to 1MB chunks) */
280f453b92SJed Brown   alloc = PetscMax(s->used+count,PetscMin(1000000/s->unitbytes+1,s->alloc+s->tailused));
290f453b92SJed Brown   ierr  = PetscMalloc(offsetof(struct _n_PetscSegBuffer,u)+alloc*s->unitbytes,&newseg);CHKERRQ(ierr);
300f453b92SJed Brown   ierr  = PetscMemzero(newseg,offsetof(struct _n_PetscSegBuffer,u));CHKERRQ(ierr);
310f453b92SJed Brown 
320f453b92SJed Brown   newseg->unitbytes = s->unitbytes;
330f453b92SJed Brown   newseg->tailused  = s->used + s->tailused;
340f453b92SJed Brown   newseg->tail      = s;
350f453b92SJed Brown   newseg->alloc     = alloc;
360f453b92SJed Brown   *seg              = newseg;
370f453b92SJed Brown   PetscFunctionReturn(0);
380f453b92SJed Brown }
390f453b92SJed Brown 
400f453b92SJed Brown #undef __FUNCT__
410f453b92SJed Brown #define __FUNCT__ "PetscSegBufferCreate"
420f453b92SJed Brown /*@C
430f453b92SJed Brown    PetscSegBufferCreate - create segmented buffer
440f453b92SJed Brown 
450f453b92SJed Brown    Not Collective
460f453b92SJed Brown 
470f453b92SJed Brown    Input Arguments:
480f453b92SJed Brown +  unitbytes - number of bytes that each entry will contain
490f453b92SJed Brown -  expected - expected/typical number of entries
500f453b92SJed Brown 
510f453b92SJed Brown    Output Argument:
520f453b92SJed Brown .  seg - segmented buffer object
530f453b92SJed Brown 
540f453b92SJed Brown    Level: developer
550f453b92SJed Brown 
56*47452e7bSJed Brown .seealso: PetscSegBufferGet(), PetscSegBufferExtractAlloc(), PetscSegBufferExtractTo(), PetscSegBufferExtractInPlace(), PetscSegBufferDestroy()
570f453b92SJed Brown @*/
580f453b92SJed Brown PetscErrorCode PetscSegBufferCreate(PetscInt unitbytes,PetscInt expected,PetscSegBuffer *seg)
590f453b92SJed Brown {
600f453b92SJed Brown   PetscErrorCode ierr;
610f453b92SJed Brown 
620f453b92SJed Brown   PetscFunctionBegin;
630f453b92SJed Brown   ierr = PetscMalloc(offsetof(struct _n_PetscSegBuffer,u)+expected*unitbytes,seg);CHKERRQ(ierr);
640f453b92SJed Brown   ierr = PetscMemzero(*seg,offsetof(struct _n_PetscSegBuffer,u));CHKERRQ(ierr);
650f453b92SJed Brown 
660f453b92SJed Brown   (*seg)->unitbytes = unitbytes;
670f453b92SJed Brown   (*seg)->alloc     = expected;
680f453b92SJed Brown   PetscFunctionReturn(0);
690f453b92SJed Brown }
700f453b92SJed Brown 
710f453b92SJed Brown #undef __FUNCT__
720f453b92SJed Brown #define __FUNCT__ "PetscSegBufferGet"
730f453b92SJed Brown /*@C
740f453b92SJed Brown    PetscSegBufferGet - get new buffer space from a segmented buffer
750f453b92SJed Brown 
760f453b92SJed Brown    Not Collective
770f453b92SJed Brown 
780f453b92SJed Brown    Input Arguments:
790f453b92SJed Brown +  seg - address of segmented buffer
800f453b92SJed Brown -  count - number of entries needed
810f453b92SJed Brown 
820f453b92SJed Brown    Output Argument:
830f453b92SJed Brown .  buf - address of new buffer for contiguous data
840f453b92SJed Brown 
850f453b92SJed Brown    Level: developer
860f453b92SJed Brown 
87*47452e7bSJed Brown .seealso: PetscSegBufferCreate(), PetscSegBufferExtractAlloc(), PetscSegBufferExtractTo(), PetscSegBufferExtractInPlace(), PetscSegBufferDestroy()
880f453b92SJed Brown @*/
890f453b92SJed Brown PetscErrorCode PetscSegBufferGet(PetscSegBuffer *seg,PetscInt count,void *buf)
900f453b92SJed Brown {
910f453b92SJed Brown   PetscErrorCode ierr;
920f453b92SJed Brown   PetscSegBuffer s;
930f453b92SJed Brown 
940f453b92SJed Brown   PetscFunctionBegin;
950f453b92SJed Brown   s = *seg;
960f453b92SJed Brown   if (PetscUnlikely(s->used + count > s->alloc)) {ierr = PetscSegBufferAlloc_Private(seg,count);CHKERRQ(ierr);}
970f453b92SJed Brown   s = *seg;
980f453b92SJed Brown   *(char**)buf = &s->u.array[s->used*s->unitbytes];
990f453b92SJed Brown   s->used += count;
1000f453b92SJed Brown   PetscFunctionReturn(0);
1010f453b92SJed Brown }
1020f453b92SJed Brown 
1030f453b92SJed Brown #undef __FUNCT__
1040f453b92SJed Brown #define __FUNCT__ "PetscSegBufferDestroy"
1050f453b92SJed Brown /*@C
1060f453b92SJed Brown    PetscSegBufferDestroy - destroy segmented buffer
1070f453b92SJed Brown 
1080f453b92SJed Brown    Not Collective
1090f453b92SJed Brown 
1100f453b92SJed Brown    Input Arguments:
1110f453b92SJed Brown .  seg - address of segmented buffer object
1120f453b92SJed Brown 
1130f453b92SJed Brown    Level: developer
1140f453b92SJed Brown 
1150f453b92SJed Brown .seealso: PetscSegBufferCreate()
1160f453b92SJed Brown @*/
1170f453b92SJed Brown PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *seg)
1180f453b92SJed Brown {
1190f453b92SJed Brown   PetscErrorCode ierr;
1200f453b92SJed Brown   PetscSegBuffer s;
1210f453b92SJed Brown 
1220f453b92SJed Brown   PetscFunctionBegin;
1230f453b92SJed Brown   for (s=*seg; s;) {
1240f453b92SJed Brown     PetscSegBuffer tail = s->tail;
1250f453b92SJed Brown     ierr = PetscFree(s);CHKERRQ(ierr);
1260f453b92SJed Brown     s = tail;
1270f453b92SJed Brown   }
1280f453b92SJed Brown   *seg = NULL;
1290f453b92SJed Brown   PetscFunctionReturn(0);
1300f453b92SJed Brown }
1310f453b92SJed Brown 
1320f453b92SJed Brown #undef __FUNCT__
133*47452e7bSJed Brown #define __FUNCT__ "PetscSegBufferExtractTo"
1340f453b92SJed Brown /*@C
135*47452e7bSJed Brown    PetscSegBufferExtractTo - extract contiguous data to provided buffer and reset segmented buffer
136*47452e7bSJed Brown 
137*47452e7bSJed Brown    Not Collective
138*47452e7bSJed Brown 
139*47452e7bSJed Brown    Input Argument:
140*47452e7bSJed Brown +  seg - segmented buffer
141*47452e7bSJed Brown -  contig - allocated buffer to hold contiguous data
142*47452e7bSJed Brown 
143*47452e7bSJed Brown    Level: developer
144*47452e7bSJed Brown 
145*47452e7bSJed Brown .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferDestroy(), PetscSegBufferExtractAlloc(), PetscSegBufferExtractInPlace()
146*47452e7bSJed Brown @*/
147*47452e7bSJed Brown PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer *seg,void *contig)
148*47452e7bSJed Brown {
149*47452e7bSJed Brown   PetscErrorCode ierr;
150*47452e7bSJed Brown   PetscInt       unitbytes;
151*47452e7bSJed Brown   PetscSegBuffer s,t;
152*47452e7bSJed Brown   char           *ptr;
153*47452e7bSJed Brown 
154*47452e7bSJed Brown   PetscFunctionBegin;
155*47452e7bSJed Brown   s = *seg;
156*47452e7bSJed Brown 
157*47452e7bSJed Brown   unitbytes = s->unitbytes;
158*47452e7bSJed Brown 
159*47452e7bSJed Brown   ptr  = ((char*)contig) + s->tailused*unitbytes;
160*47452e7bSJed Brown   ierr = PetscMemcpy(ptr,s->u.array,s->used*unitbytes);CHKERRQ(ierr);
161*47452e7bSJed Brown   for (t=s->tail; t;) {
162*47452e7bSJed Brown     PetscSegBuffer tail = t->tail;
163*47452e7bSJed Brown     ptr -= t->used*unitbytes;
164*47452e7bSJed Brown     ierr = PetscMemcpy(ptr,t->u.array,t->used*unitbytes);CHKERRQ(ierr);
165*47452e7bSJed Brown     ierr = PetscFree(t);CHKERRQ(ierr);
166*47452e7bSJed Brown     t    = tail;
167*47452e7bSJed Brown   }
168*47452e7bSJed Brown   if (ptr != contig) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Tail count does not match");
169*47452e7bSJed Brown   s->used             = 0;
170*47452e7bSJed Brown   s->tailused         = 0;
171*47452e7bSJed Brown   s->tail             = NULL;
172*47452e7bSJed Brown   PetscFunctionReturn(0);
173*47452e7bSJed Brown }
174*47452e7bSJed Brown 
175*47452e7bSJed Brown #undef __FUNCT__
176*47452e7bSJed Brown #define __FUNCT__ "PetscSegBufferExtractAlloc"
177*47452e7bSJed Brown /*@C
178*47452e7bSJed Brown    PetscSegBufferExtractAlloc - extract contiguous data to new allocation and reset segmented buffer
1790f453b92SJed Brown 
1800f453b92SJed Brown    Not Collective
1810f453b92SJed Brown 
1820f453b92SJed Brown    Input Argument:
1830f453b92SJed Brown .  seg - segmented buffer
1840f453b92SJed Brown 
1850f453b92SJed Brown    Output Argument:
1860f453b92SJed Brown .  contiguous - address of new array containing contiguous data, caller frees with PetscFree()
1870f453b92SJed Brown 
1880f453b92SJed Brown    Level: developer
1890f453b92SJed Brown 
190*47452e7bSJed Brown    Developer Notes: 'seg' argument is a pointer so that implementation could reallocate, though this is not currently done
191*47452e7bSJed Brown 
192*47452e7bSJed Brown .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferDestroy(), PetscSegBufferExtractTo(), PetscSegBufferExtractInPlace()
1930f453b92SJed Brown @*/
194*47452e7bSJed Brown PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer *seg,void *contiguous)
1950f453b92SJed Brown {
1960f453b92SJed Brown   PetscErrorCode ierr;
197*47452e7bSJed Brown   PetscSegBuffer s;
198*47452e7bSJed Brown   void           *contig;
1990f453b92SJed Brown 
2000f453b92SJed Brown   PetscFunctionBegin;
2010f453b92SJed Brown   s = *seg;
2020f453b92SJed Brown 
203*47452e7bSJed Brown   ierr = PetscMalloc((s->used+s->tailused)*s->unitbytes,&contig);CHKERRQ(ierr);
204*47452e7bSJed Brown   ierr = PetscSegBufferExtractTo(seg,contig);CHKERRQ(ierr);
205*47452e7bSJed Brown   *(void**)contiguous = contig;
206*47452e7bSJed Brown   PetscFunctionReturn(0);
2070f453b92SJed Brown }
208*47452e7bSJed Brown 
209*47452e7bSJed Brown #undef __FUNCT__
210*47452e7bSJed Brown #define __FUNCT__ "PetscSegBufferExtractInPlace"
211*47452e7bSJed Brown /*@C
212*47452e7bSJed Brown    PetscSegBufferExtractInPlace - extract in-place contiguous representation of data and reset segmented buffer for reuse
213*47452e7bSJed Brown 
214*47452e7bSJed Brown    Collective
215*47452e7bSJed Brown 
216*47452e7bSJed Brown    Input Arguments:
217*47452e7bSJed Brown .  seg - segmented buffer object
218*47452e7bSJed Brown 
219*47452e7bSJed Brown    Output Arguments:
220*47452e7bSJed Brown .  contig - address of pointer to contiguous memory
221*47452e7bSJed Brown 
222*47452e7bSJed Brown    Level: developer
223*47452e7bSJed Brown 
224*47452e7bSJed Brown .seealso: PetscSegBufferExtractAlloc(), PetscSegBufferExtractTo()
225*47452e7bSJed Brown @*/
226*47452e7bSJed Brown PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer *seg,void *contig)
227*47452e7bSJed Brown {
228*47452e7bSJed Brown   PetscErrorCode ierr;
229*47452e7bSJed Brown 
230*47452e7bSJed Brown   PetscFunctionBegin;
231*47452e7bSJed Brown   if (!(*seg)->tail) {
232*47452e7bSJed Brown     *(char**)contig = (*seg)->u.array;
233*47452e7bSJed Brown   } else {
234*47452e7bSJed Brown     PetscSegBuffer s = *seg,newseg;
235*47452e7bSJed Brown 
236*47452e7bSJed Brown     ierr = PetscSegBufferCreate(s->unitbytes,s->used+s->tailused,&newseg);CHKERRQ(ierr);
237*47452e7bSJed Brown     ierr = PetscSegBufferExtractTo(seg,newseg->u.array);CHKERRQ(ierr);
238*47452e7bSJed Brown     ierr = PetscSegBufferDestroy(seg);CHKERRQ(ierr);
239*47452e7bSJed Brown     *seg = newseg;
240*47452e7bSJed Brown     *(void**)contig = newseg->u.array;
241*47452e7bSJed Brown   }
2420f453b92SJed Brown   PetscFunctionReturn(0);
2430f453b92SJed Brown }
244