DynListI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 template<class T, Foam::label staticSize>
28 {
29  return dataPtr_;
30 }
31 
32 template<class T, Foam::label staticSize>
33 inline const T* Foam::DynList<T, staticSize>::data() const
34 {
35  return dataPtr_;
36 }
37 
38 template<class T, Foam::label staticSize>
40 {
41  checkAllocation();
42 
43  if( s > staticSize )
44  {
45  if( s > nAllocated_ )
46  {
47  //- allocates enough space for the elements
48  T* newData = new T[s];
49 
50  for(label i=0;i<nextFree_;++i)
51  newData[i] = this->operator[](i);
52 
53  if( nAllocated_ > staticSize )
54  delete [] dataPtr_;
55 
56  dataPtr_ = newData;
57  nAllocated_ = s;
58  }
59  else if( s < nAllocated_ )
60  {
61  //- shrinks the list
62  T* newData = new T[s];
63 
64  for(label i=0;i<s;++i)
65  newData[i] = this->operator[](i);
66 
67  delete [] dataPtr_;
68 
69  dataPtr_ = newData;
70  nAllocated_ = s;
71  }
72  }
73  else
74  {
75  if( nAllocated_ > staticSize )
76  {
77  //- delete dynamically allocated data
78  for(label i=0;i<s;++i)
79  staticData_[i] = dataPtr_[i];
80 
81  delete [] dataPtr_;
82  }
83 
84  dataPtr_ = staticData_;
85  nAllocated_ = staticSize;
86  }
87 }
88 
89 template<class T, Foam::label staticSize>
91 {
92  if( (i < 0) || (i >= nextFree_) )
93  {
95  (
96  "void Foam::DynList<T, label, Offset>::"
97  "checkIndex(const label i) const"
98  ) << "Index " << i << " is not in range " << 0
99  << " and " << nextFree_ << abort(FatalError);
100  }
101 }
102 
103 template<class T, Foam::label staticSize>
105 {
106  if( nextFree_ > nAllocated_ )
108  (
109  "template<class T, Foam::label staticSize> "
110  "inline void Foam::DynList<T, staticSize>::"
111  "checkAllocation() const"
112  ) << "nextFree_ is out of scope 0 " << " and " << nAllocated_
113  << abort(FatalError);
114 }
115 
116 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117 
118 //- Construct null
119 template<class T, Foam::label staticSize>
121 :
122  dataPtr_(NULL),
123  nAllocated_(0),
124  staticData_(),
125  nextFree_(0)
126 {
127  setSize(0);
128 
129  # ifdef DEBUG
130  checkAllocation();
131  # endif
132 }
133 
134 
135 template<class T, Foam::label staticSize>
137 :
138  dataPtr_(NULL),
139  nAllocated_(0),
140  staticData_(),
141  nextFree_(0)
142 {
143  setSize(s);
144 
145  # ifdef DEBUG
146  checkAllocation();
147  # endif
148 }
149 
150 template<class T, Foam::label staticSize>
152 :
153  dataPtr_(NULL),
154  nAllocated_(0),
155  staticData_(),
156  nextFree_(0)
157 {
158  setSize(s);
159 
160  for(label i=0;i<s;++i)
161  this->operator[](i) = val;
162 
163  # ifdef DEBUG
164  checkAllocation();
165  # endif
166 }
167 
168 template<class T, Foam::label staticSize>
170 :
171  dataPtr_(NULL),
172  nAllocated_(0),
173  staticData_(),
174  nextFree_(0)
175 {
176  setSize(ul.size());
177 
178  forAll(ul, i)
179  this->operator[](i) = ul[i];
180 
181  # ifdef DEBUG
182  checkAllocation();
183  # endif
184 }
185 
186 template<class T, Foam::label staticSize>
187 template<class ListType>
188 inline Foam::DynList<T, staticSize>::DynList(const ListType& l)
189 :
190  dataPtr_(NULL),
191  nAllocated_(0),
192  staticData_(),
193  nextFree_(0)
194 {
195  setSize(l.size());
196  for(label i=0;i<nextFree_;++i)
197  this->operator[](i) = l[i];
198 
199  # ifdef DEBUG
200  checkAllocation();
201  # endif
202 }
203 
204 //- Copy construct
205 template<class T, Foam::label staticSize>
207 (
208  const DynList<T, staticSize>& dl
209 )
210 :
211  dataPtr_(NULL),
212  nAllocated_(0),
213  staticData_(),
214  nextFree_(0)
215 {
216  setSize(dl.size());
217  for(label i=0;i<nextFree_;++i)
218  this->operator[](i) = dl[i];
219 
220  # ifdef DEBUG
221  checkAllocation();
222  # endif
223 }
224 
225 template<class T, Foam::label staticSize>
227 {
228  allocateSize(0);
229 }
230 
231 
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233 
234 template<class T, Foam::label staticSize>
236 {
237  # ifdef DEBUG
238  checkAllocation();
239  # endif
240 
241  return nextFree_;
242 }
243 
244 template<class T, Foam::label staticSize>
246 {
247  # ifdef DEBUG
248  checkAllocation();
249  # endif
250 
251  if( !contiguous<T>() )
252  {
253  FatalErrorIn("DynList<T>::byteSize()")
254  << "Cannot return the binary size of a list of "
255  "non-primitive elements"
256  << abort(FatalError);
257  }
258 
259  return nextFree_*sizeof(T);
260 }
261 
262 template<class T, Foam::label staticSize>
264 {
265  # ifdef DEBUG
266  checkAllocation();
267  # endif
268 
269  allocateSize(s);
270  nextFree_ = s;
271 
272  # ifdef DEBUG
273  checkAllocation();
274  # endif
275 }
276 
277 
278 template<class T, Foam::label staticSize>
280 {
281  # ifdef DEBUG
282  checkAllocation();
283  # endif
284 
285  nextFree_ = 0;
286 }
287 
288 
289 template<class T, Foam::label staticSize>
291 {
292  # ifdef DEBUG
293  checkAllocation();
294  # endif
295 
296  allocateSize(nextFree_);
297 
298  # ifdef DEBUG
299  checkAllocation();
300  # endif
301 }
302 
303 template<class T, Foam::label staticSize>
305 {
306  # ifdef DEBUG
307  checkAllocation();
308  # endif
309 
310  if( nextFree_ >= nAllocated_ )
311  {
312  const label newSize = 2*nAllocated_+2;
313  allocateSize(newSize);
314  }
315 
316  # ifdef DEBUG
317  checkAllocation();
318  # endif
319 
320  this->operator[](nextFree_++) = e;
321 }
322 
323 template<class T, Foam::label staticSize>
325 {
326  # ifdef DEBUG
327  checkAllocation();
328  # endif
329 
330  if( !contains(e) )
331  append(e);
332 
333  # ifdef DEBUG
334  checkAllocation();
335  # endif
336 }
337 
338 template<class T, Foam::label staticSize>
339 inline bool Foam::DynList<T, staticSize>::contains(const T& e) const
340 {
341  # ifdef DEBUG
342  checkAllocation();
343  # endif
344 
345  for(label i=0;i<nextFree_;++i)
346  {
347  if( this->operator[](i) == e )
348  return true;
349  }
350 
351  return false;
352 }
353 
354 template<class T, Foam::label staticSize>
356 (
357  const T& e
358 ) const
359 {
360  # ifdef DEBUG
361  checkAllocation();
362  # endif
363 
364  for(label i=0;i<nextFree_;++i)
365  {
366  if( this->operator[](i) == e )
367  return i;
368  }
369 
370  return -1;
371 }
372 
373 template<class T, Foam::label staticSize>
375 {
376  # ifdef DEBUG
377  checkAllocation();
378  # endif
379 
380  return this->operator[](nextFree_-1);
381 }
382 
383 template<class T, Foam::label staticSize>
385 {
386  # ifdef DEBUG
387  checkAllocation();
388  # endif
389 
390  if( nextFree_ == 0 )
391  {
393  (
394  "void Foam::DynList<T, staticSize>::remove()"
395  ) << "List is empty" << abort(FatalError);
396  }
397 
398  T el = this->operator[](nextFree_-1);
399  --nextFree_;
400  return el;
401 }
402 
403 template<class T, Foam::label staticSize>
405 {
406  # ifdef DEBUG
407  checkAllocation();
408  # endif
409 
410  if( nextFree_ == 0 )
411  {
413  (
414  "void Foam::DynList<T, staticSize>::remove()"
415  ) << "List is empty" << abort(FatalError);
416  }
417 
418  T el = this->operator[](i);
419  this->operator[](i) = this->operator[](nextFree_-1);
420  --nextFree_;
421 
422  # ifdef DEBUG
423  checkAllocation();
424  # endif
425 
426  return el;
427 }
428 
429 template<class T, Foam::label staticSize>
431 {
432  # ifdef DEBUG
433  checkAllocation();
434  # endif
435 
436  return this->operator()(i);
437 }
438 
439 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
440 
441 template<class T, Foam::label staticSize>
443 {
444  # ifdef DEBUG
445  checkAllocation();
446  # endif
447 
448  nextFree_ = Foam::max(nextFree_, i + 1);
449 
450  if( nextFree_ >= nAllocated_ )
451  {
452  allocateSize(2 * nextFree_+1);
453  }
454 
455  # ifdef DEBUG
456  checkAllocation();
457  # endif
458 
459  return this->operator[](i);
460 }
461 
462 template<class T, Foam::label staticSize>
464 {
465  # ifdef FULLDEBUG
466  checkAllocation();
467  checkIndex(i);
468  # endif
469 
470  return dataPtr_[i];
471 }
472 
473 template<class T, Foam::label staticSize>
475 {
476  # ifdef FULLDEBUG
477  checkAllocation();
478  checkIndex(i);
479  # endif
480 
481  return dataPtr_[i];
482 }
483 
484 template<class T, Foam::label staticSize>
486 (
487  const label index,
488  const label offset
489 ) const
490 {
491  return (index + offset) % nextFree_;
492 }
493 
494 template<class T, Foam::label staticSize>
496 (
497  const label index,
498  const label offset
499 ) const
500 {
501  return (index + nextFree_ - offset) % nextFree_;
502 }
503 
504 template<class T, Foam::label staticSize>
506 (
507  const label index,
508  const label offset
509 ) const
510 {
511  return operator[](fcIndex(index, offset));
512 }
513 
514 template<class T, Foam::label staticSize>
516 (
517  const label index,
518  const label offset
519 ) const
520 {
521  return operator[](rcIndex(index, offset));
522 }
523 
524 template<class T, Foam::label staticSize>
526 {
527  # ifdef DEBUG
528  checkAllocation();
529  # endif
530 
531  for(label i=0;i<nextFree_;++i)
532  operator[](i) = t;
533 }
534 
535 template<class T, Foam::label staticSize>
537 (
538  const DynList<T, staticSize>& dl
539 )
540 {
541  # ifdef DEBUG
542  checkAllocation();
543  # endif
544 
545  allocateSize(dl.size());
546  nextFree_ = dl.size();
547 
548  # ifdef DEBUG
549  checkAllocation();
550  # endif
551 
552  for(label i=0;i<nextFree_;++i)
553  this->operator[](i) = dl[i];
554 }
555 
556 template<class T, Foam::label staticSize>
557 template<class ListType>
558 inline void Foam::DynList<T, staticSize>::operator=(const ListType& l)
559 {
560  # ifdef DEBUG
561  checkAllocation();
562  # endif
563 
564  allocateSize(l.size());
565  nextFree_ = l.size();
566 
567  # ifdef DEBUG
568  checkAllocation();
569  # endif
570 
571  for(label i=0;i<nextFree_;++i)
572  this->operator[](i) = l[i];
573 }
574 
575 template<class T, Foam::label staticSize>
577 (
578  const DynList<T, staticSize>& DL
579 ) const
580 {
581  if( nextFree_ != DL.nextFree_ )
582  return false;
583 
584  forAll(DL, i)
585  if( this->operator[](i) != DL[i] )
586  return false;
587 
588  return true;
589 }
590 
591 template<class T, Foam::label staticSize>
593 (
594  const DynList<T, staticSize>& DL
595 ) const
596 {
597  return !operator==(DL);
598 }
599 
600 
601 // ************************************************************************* //
setSize
points setSize(newPointi)
Foam::DynList::newElmt
T & newElmt(const label)
return a refence to the element. Resize the list if necessary
Definition: DynListI.H:430
Foam::DynList::~DynList
~DynList()
Definition: DynListI.H:226
Foam::DynList::data
T * data()
access to the data pointer
Definition: DynListI.H:27
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynList::operator()
T & operator()(const label)
Definition: DynListI.H:442
Foam::DynList::removeElement
T removeElement(const label i)
Definition: DynListI.H:404
Foam::DynList::fcIndex
label fcIndex(const label index, const label offset=1) const
return forward and reverse circular indices
Definition: DynListI.H:486
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::DynList::byteSize
label byteSize() const
Definition: DynListI.H:245
Foam::DynList::removeLastElement
T removeLastElement()
Return and remove the last element.
Definition: DynListI.H:384
Foam::DynList::shrink
void shrink()
Shrink the List<T> to the number of elements used.
Definition: DynListI.H:290
Foam::DynList::contains
bool contains(const T &e) const
check if the element is in the list (takes linear time)
Definition: DynListI.H:339
Foam::DynList::DynList
DynList()
Construct null.
Definition: DynListI.H:120
Foam::DynList::rcElement
const T & rcElement(const label index, const label offset=1) const
Definition: DynListI.H:516
Foam::DynList::allocateSize
void allocateSize(const label)
allocate list size
Definition: DynListI.H:39
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::DynList::operator=
void operator=(const T &)
Assignment of all entries to the given value.
Definition: DynListI.H:525
Foam::DynList::appendIfNotIn
void appendIfNotIn(const T &e)
Definition: DynListI.H:324
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::DynList
Definition: DynList.H:53
Foam::DynList::checkIndex
void checkIndex(const label) const
check if index is inside the scope (used for debugging only)
Definition: DynListI.H:90
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::DynList::rcIndex
label rcIndex(const label index, const label offset=1) const
Definition: DynListI.H:496
Foam::DynList::containsAtPosition
label containsAtPosition(const T &e) const
Definition: DynListI.H:356
T
const volScalarField & T
Definition: createFields.H:25
Foam::DynList::lastElement
const T & lastElement() const
return a const reference to the last element
Definition: DynListI.H:374
Foam::DynList::setSize
void setSize(const label)
Reset size of List.
Definition: DynListI.H:263
Foam::UList< T >
Foam::DynList::checkAllocation
void checkAllocation() const
check if nAllocated_ is greater or equal to nextFree_
Definition: DynListI.H:104
Foam::DynList::fcElement
const T & fcElement(const label index, const label offset=1) const
return forward and reverse circular elements
Definition: DynListI.H:506
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::DynList::size
label size() const
Definition: DynListI.H:235
Foam::UList::size
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::DynList::operator[]
const T & operator[](const label) const
return access to an element
Definition: DynListI.H:463
Foam::DynList::clear
void clear()
Clear the list, i.e. set next free to zero.
Definition: DynListI.H:279
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304