Overture  Version 25
GridCollection.h
Go to the documentation of this file.
1 #ifndef _GridCollection
2 #define _GridCollection
3 
4 //
5 // Who to blame: Geoff Chesshire
6 //
7 
10 #include "MappedGrid.h"
12 #ifndef USE_STL
13 #include "ListOfMappedGrid.h"
14 #include "ListOfGridCollection.h"
15 #endif // USE_STL
16 #include "GCMath.h"
17 
18 class AMR_RefinementLevelInfo;
19 class Interpolant;
20 
21 class GridCollection;
23 
24 //
25 // Class for reference-counted data.
26 //
29  public:
30  enum {
50  // THEminMaxEdgeLength = MappedGridData::THEminMaxEdgeLength,
73  };
78 #ifdef USE_STL
79  RCVector<MappedGrid> grid; // (overloaded)
80  RCVector<GridCollection> baseGrid; // (overloaded)
81  RCVector<GridCollection> refinementLevel; // (overloaded)
82  RCVector<GridCollection> componentGrid; // (overloaded)
83  RCVector<GridCollection> multigridLevel; // (overloaded)
84  RCVector<GridCollection> domain; // (overloaded)
85 #else
86  ListOfMappedGrid grid; // (overloaded)
87  ListOfGridCollection baseGrid; // (overloaded)
91  ListOfGridCollection domain; // (overloaded)
92 #endif // USE_STL
94 
97 
98 // ********
99 // // here we move the stuff that used to be in the CompositeGrid Class
100 // IntegerArray *numberOfInterpolationPoints;
101 // IntegerArray *interpolationWidth;
102 
103  // arrays needed for hybrid grids:
104 // IntegerArray *numberOfHybridInterfaceNodes; // numberOfHybridNodes(grid)
105 // IntegerArray *hybridInterfaceGridIndex; // hybridInterfaceGridIndex(numberOfHybridBoundaryNodes,4)
106 // IntegerArray *hybridGlobalID; // hybridGlobalID[grid](numberOfHybridNodes(grid))
107 // ********
108 
110  const Integer numberOfDimensions_ = 0,
111  const Integer numberOfGrids_ = 0);
113  const GridCollectionData& x,
114  const CopyType ct = DEEP);
115  virtual ~GridCollectionData();
117  inline MappedGrid& operator[](const Integer& i) { return grid[i]; }
118  inline const MappedGrid& operator[](const Integer& i) const
119  { return grid[i]; }
120  void reference(const GridCollectionData& x);
121  virtual void breakReference();
122  virtual void consistencyCheck() const;
123  virtual Integer get(
124  const GenericDataBase& db,
125  const aString& name);
126  virtual Integer put(GenericDataBase& db,
127  const aString& name,
128  int geometryToPut = -1 // by default put computedGeometry
129  ) const;
130  inline Integer update(
131  const Integer what = THEusualSuspects,
132  const Integer how = COMPUTEtheUsual)
133  { return update(*this, what, how); }
134  inline Integer update(
136  const Integer what = THEusualSuspects,
137  const Integer how = COMPUTEtheUsual)
138  { return update((GenericGridCollectionData&)x, what, how); }
139 
141 
142  virtual void destroy(const Integer what = NOTHING);
143  virtual Integer addRefinement(
144  const IntegerArray& range,
145  const IntegerArray& factor,
146  const Integer& level,
147  const Integer k = 0);
148 
149  virtual Integer replaceRefinementLevels(int level0, int numberOfRefinementLevels0, IntegerArray **gridInfo );
150 
151  virtual void deleteRefinement(const Integer& k);
152  virtual void deleteRefinementLevels(const Integer level = 0);
155  const Integer level = INTEGER_MAX)
158  const IntegerArray& factor,
159  const Integer& level,
160  const Integer k = 0);
161  virtual void deleteMultigridCoarsening(const Integer& k);
162  virtual void deleteMultigridLevels(const Integer level = 0);
163  virtual void setNumberOfGrids(const Integer& numberOfGrids_);
164  virtual void setNumberOfDimensions(const Integer& numberOfDimensions_);
165  virtual void setNumberOfDimensionsAndGrids(
166  const Integer& numberOfDimensions_,
167  const Integer& numberOfGrids_);
168  void initialize(
169  const Integer& numberOfDimensions_,
170  const Integer& numberOfGrids_);
171 
172  protected:
173 
174  virtual Integer update(
176  const Integer what = THEusualSuspects,
177  const Integer how = COMPUTEtheUsual);
178  virtual void referenceRefinementLevels(
180  const Integer level = INTEGER_MAX);
182  const Integer& what,
183  Integer& numberOfCollections,
184 #ifdef USE_STL
185  RCVector<GridCollection>& list,
186  RCVector<GenericGridCollection>& genericList,
187 #else
188  ListOfGridCollection& list,
189  ListOfGenericGridCollection& genericList,
190 #endif // USE_STL
191  IntegerArray& collection);
192 //
193  int updateRefinementGrid( int n, int b, int p,
194  const IntegerArray& range,
195  const IntegerArray& factor,
196  const Integer& level );
197 
198 // The following functions are declared protected here in order to disallow
199 // access to the corresponding functions in class GenericGridCollectionData.
200 //
201  inline virtual Integer addRefinement(
202  const Integer& level,
203  const Integer k = 0) {
204  if (&level || &k); // Avoid compiler warnings.
205  cerr << "virtual void GridCollectionData::addRefinement(const Integer& level, const Integer k) must not be called!"
206  << "It must have been called illegally through the base class GenericGridCollectionData."
207  << endl;
208  abort();
209  return -1;
210  }
212  const Integer& level,
213  const Integer k = 0) {
214  if (&level || &k); // Avoid compiler warnings.
215  cerr << "virtual void GridCollectionData::addMultigridCoarsening(const Integer& level, const Integer k) must not be called!"
216  << "It must have been called illegally through the base class GenericGridCollectionData."
217  << endl;
218  abort();
219  return -1;
220  }
221 //
222 // Virtual member functions used only through class ReferenceCounting:
223 //
224  private:
225  inline virtual ReferenceCounting& operator=(const ReferenceCounting& x)
226  { return operator=((GridCollectionData&)x); }
227  inline virtual void reference(const ReferenceCounting& x)
228  { reference((GridCollectionData&)x); }
229  inline virtual ReferenceCounting* virtualConstructor(
230  const CopyType ct = DEEP) const
231  { return new GridCollectionData(*this, ct); }
232  aString className;
233  public:
234  inline virtual aString getClassName() const { return className; }
235 };
236 
238  public GenericGridCollection {
239  public:
240 // Public constants:
241 
242 // Constants to be ORed to form the first argument of update() and destroy():
243  enum {
264  // THEminMaxEdgeLength = GridCollectionData::THEminMaxEdgeLength,
268  };
269 
270 // Constants to be ORed to form the second argument of update():
271  enum {
277  };
278 
279 // Mask bits to access data in mask.
280  enum {
292  };
293 
294 // Public data:
295 //
299 
301 
302 #ifdef USE_STL
303 // Lists of component grids.
304  RCVector<MappedGrid> grid; // (overloaded)
305 
306 // Collections of grids having the same base grid.
307  RCVector<GridCollection> baseGrid; // (overloaded)
308 
309 // Collections of grids having the same refinement level.
310  RCVector<GridCollection> refinementLevel; // (overloaded)
311 
312 // Collections of multigrid coarsenings of the same component grid
313  RCVector<GridCollection> componentGrid; // (overloaded)
314 
315 // Collections of grids having the same multigrid level.
316  RCVector<GridCollection> multigridLevel; // (overloaded)
317 
318 // Collections of grids in the same domain
319  RCVector<GridCollection> domain; // (overloaded)
320 #else
321 // Lists of component grids.
322  ListOfMappedGrid grid; // (overloaded)
323 
324 // Collections of grids having the same base grid.
326 
327 // Collections of grids having the same refinement level.
329 
330 // Collections of multigrid coarsenings of the same component grid
332 
333 // Collections of grids having the same multigrid level.
335 
336 // Collections of grids in the same domain
337  ListOfGridCollection domain; // (overloaded)
338 #endif // USE_STL
339 
340  AMR_RefinementLevelInfo* refinementLevelInfo;
341 //
342 // Public member functions for access to data:
343 //
344 // The number of dimensions of the MappedGrids in grid.
345 //
347  { return rcData->numberOfDimensions; }
348  inline const Integer& numberOfDimensions() const
349  { return rcData->numberOfDimensions; }
350 //
351 // Public member functions:
352 //
353 // Default constructor.
354 //
355 // Create a GridCollection with the given number of dimensions
356 // and number of component grids.
357 //
359  const Integer numberOfDimensions_ = 0,
360  const Integer numberOfGrids_ = 0);
361 //
362 // Copy constructor. (Does deep copy by default.)
363 //
365  const GridCollection& x,
366  const CopyType ct = DEEP);
367 //
368 // Destructor.
369 //
370  virtual ~GridCollection();
371 //
372 // Assignment operator. (Does a deep copy.)
373 //
375 //
376 // Get a reference to a component grid using C or Fortran indexing.
377 //
378  inline MappedGrid& operator[](const Integer& i) { return grid[i]; }
379  inline const MappedGrid& operator[](const Integer& i) const
380  { return grid[i]; }
381 //
382 // Make a reference. (Does a shallow copy.)
383 //
384  void reference(const GridCollection& x);
385  void reference(GridCollectionData& x);
386 //
387 // Break a reference. (Replaces with a deep copy.)
388 //
389  virtual void breakReference();
390 //
391 // Change the grid to be all vertex-centered.
392 //
393  virtual void changeToAllVertexCentered();
394 //
395 // Change the grid to be all cell-centered.
396 //
397  virtual void changeToAllCellCentered();
398 //
399 // Check that the data structure is self-consistent.
400 //
401  virtual void consistencyCheck() const;
402 
403 // display parallel distribution
404  void displayDistribution(const aString & label, FILE *file=stdout );
405 
406  int numberOfGridPoints() const; // return the number of grid points
407 
408  // return a pointer to the Interpolant.
410 //
411 // "Get" and "put" database operations.
412 //
413  virtual Integer get(
414  const GenericDataBase& db,
415  const aString& name);
416  virtual Integer put(GenericDataBase& db,
417  const aString& name,
418  int geometryToPut = -1 // by default put computedGeometry
419  ) const;
420 //
421 // Set references to reference-counted data.
422 //
423  void updateReferences(const Integer what = EVERYTHING);
424 //
425 // Update the grid.
426 //
427  inline Integer update(
428  const Integer what = THEusualSuspects,
429  const Integer how = COMPUTEtheUsual)
430  { return update(*this, what, how); }
431 //
432 // Update the grid, sharing the data of another grid.
433 //
434  inline Integer update(
435  GridCollection& x,
436  const Integer what = THEusualSuspects,
437  const Integer how = COMPUTEtheUsual)
438  { return update((GenericGridCollection&)x, what, how); }
439 
440  // get the list of parent/child/siblings for AMR
442 
443  // update the parent/child/sibling info for AMR
444  virtual void updateParentChildSiblingInfo();
445 
446  // Call this function when the grid collection has changed
448 
449 //
450 // Destroy optional grid data.
451 //
452  virtual void destroy(const Integer what = NOTHING);
453 
454 // Add a new grid.
455  virtual int add(MappedGrid & g);
456 // Add a new grid, built from a Mapping
457  virtual int add(Mapping & map);
458 
459  // delete a grid:
460  virtual int deleteGrid(Integer k);
461 
462  // delete a list of grids:
463  virtual int deleteGrid(const IntegerArray & gridsToDelete );
464 
465  // delete a list of grids (use this when you also know the list of grids to save).
466  virtual int deleteGrid(const IntegerArray & gridsToDelete, const IntegerArray & gridsToSave );
467 
468 //
469 // Add a refinement grid to the collection.
470 //
471  virtual Integer addRefinement(
472  const IntegerArray& range, // The indexRange of the refinement grid.
473  const IntegerArray& factor, // The refinement factor w.r.t. level-1.
474  const Integer& level, // The refinement level number.
475  const Integer k = 0); // The index of an ancestor of the refinement.
477  const IntegerArray& range,
478  const Integer& factor,
479  const Integer& level,
480  const Integer k = 0) {
481  IntegerArray factors(3); factors = factor;
482  return addRefinement(range, factors, level, k);
483  }
484 
485  // replace refinement level level0 and higher
486  virtual Integer replaceRefinementLevels(int level0, int numberOfRefinementLevels0, IntegerArray **gridInfo );
487 
488  //
489 // Delete all multigrid levels of refinement grid k.
490 //
491  virtual void deleteRefinement(const Integer& k);
492 //
493 // Delete all grids with refinement level greater than the given level.
494 //
495  virtual void deleteRefinementLevels(const Integer level = 0);
496 //
497 // Reference x[i] for refinementLevelNumber(i) <= level.
498 // Delete all other grids.
499 //
501  GridCollection& x,
502  const Integer level = INTEGER_MAX)
504 //
505 // Add a multigrid coarsening of grid k.
506 //
508  const IntegerArray& factor, // The coarsening factor w.r.t level-1
509  const Integer& level, // The multigrid level number.
510  const Integer k = 0); // The index of the corresponding grid
511  // at any finer multigrid level.
512 //
513 // Delete grid k, a multigrid coarsening, and all of its multigrid coarsenings.
514 //
515  virtual void deleteMultigridCoarsening(const Integer& k);
516 //
517 // Delete all of the grids with multigrid level greater than the given level.
518 //
519  virtual void deleteMultigridLevels(const Integer level = 0);
520 
521 // Specify the grids that belong to a domain
522  virtual void addToDomain( int domain, const IntegerArray & gridList );
523 
524 //
525 // Set the number of grids.
526 //
527  virtual void setNumberOfGrids(const Integer& numberOfGrids_);
528 //
529 // Set the number of dimensions.
530 //
531  virtual void setNumberOfDimensions(const Integer& numberOfDimensions_);
532 //
533 // Set the number of dimensions and grids.
534 //
535  virtual void setNumberOfDimensionsAndGrids(
536  const Integer& numberOfDimensions_,
537  const Integer& numberOfGrids_);
538 
539 
540  // return size of this object
541  virtual real sizeOf(FILE *file = NULL, bool returnSizeOfReference = false ) const;
542 
543  int setMaskAtRefinements();
544 
545 //
546 // Specify the set of processes over which GridFunctions are distributed.
547 // We now support only the specification of a contiguous range of process IDs.
548 //
549  void specifyProcesses(const Range& range);
550 
551 
552  // numberOfHybridInterfaceNodes(grid)
553 // inline IntegerArray& numberOfHybridInterfaceNodes(int grid){ return *rcData->numberOfHybridInterfaceNodes; }
554 // inline const IntegerArray& numberOfHybridInterfaceNodes(int grid) const{ return *rcData->numberOfHybridInterfaceNodes; }
555 //
556 // // hybridInterfaceGridIndex(numberOfHybridBoundaryNodes,4)
557 // inline IntegerArray& hybridInterfaceGridIndex(){ return *rcData->hybridInterfaceGridIndex; }
558 // inline const IntegerArray& hybridInterfaceGridIndex() const{ return *rcData->hybridInterfaceGridIndex; }
559 //
560 // inline IntegerArray& hybridGlobalID(){ return *rcData->hybridGlobalID; }
561 // inline const IntegerArray& hybridGlobalID() const{ return *rcData->hybridGlobalID; }
562 
563 
564 //
565 // Initialize the GridCollection with the given number of dimensions and grids.
566 // These grids have their gridNumbers, baseGridNumbers and componentGridNumbers
567 // set to [0, ..., numberOfGrids_-1], and their refinementLevelNumbers and
568 // multigridLevelNumbers set to zero.
569 //
570  virtual void initialize(
571  const Integer& numberOfDimensions_,
572  const Integer& numberOfGrids_);
573 
574  // Here is the master grid.
576 //
577 // Pointer to reference-counted data.
578 //
581  inline GridCollectionData* operator->() { return rcData; }
582  inline const GridCollectionData* operator->() const { return rcData; }
583  inline GridCollectionData& operator*() { return *rcData; }
584  inline const GridCollectionData& operator*() const { return *rcData; }
585 
586  private:
587  virtual Integer update(
589  const Integer what = THEusualSuspects,
590  const Integer how = COMPUTEtheUsual);
591  virtual void referenceRefinementLevels(
593  const Integer level = INTEGER_MAX);
594 
595  protected:
596 //
597 // The following functions are declared protected here in order to disallow
598 // access to the corresponding functions in class GenericGridCollection.
599 //
600  inline virtual Integer addRefinement(
601  const Integer& level,
602  const Integer k = 0) {
603  if (&level || &k); // Avoid compiler warnings.
604  cerr << "virtual void GridCollection::addRefinement(const Integer& level, const Integer k) must not be called!"
605  << "It must have been called illegally through the base class GenericGridCollection."
606  << endl;
607  abort();
608  return -1;
609  }
611  const Integer& level,
612  const Integer k = 0) {
613  if (&level || &k); // Avoid compiler warnings.
614  cerr << "virtual void GridCollection::addMultigridCoarsening(const Integer& level, const Integer k) must not be called!"
615  << "It must have been called illegally through the base class GenericGridCollection."
616  << endl;
617  abort();
618  return -1;
619  }
620  inline virtual void initialize(
621  const Integer& numberOfGrids_) {
622  if (&numberOfGrids_); // Avoid compiler warnings.
623  cerr << "virtual void GridCollection::initialize(const Integer& numberOfGrids_) must not be called!"
624  << "It must have been called illegally through the base class GenericGridCollection."
625  << endl;
626  abort();
627  }
628 
629  protected:
631 
632 //
633 // Virtual member functions used only through class ReferenceCounting:
634 //
635  private:
636  inline virtual ReferenceCounting& operator=(const ReferenceCounting& x)
637  { return operator=((GridCollection&)x); }
638  inline virtual void reference(const ReferenceCounting& x)
639  { reference((GridCollection&)x); }
640  inline virtual ReferenceCounting* virtualConstructor(
641  const CopyType ct = DEEP) const
642  { return new GridCollection(*this, ct); }
643  aString className;
644  public:
645  inline virtual aString getClassName() const { return className; }
646 };
647 //
648 // Stream output operator.
649 //
650 ostream& operator<<(ostream& s, const GridCollection& g);
651 
652 #endif // _GridCollection