Overture  Version 25
CompositeGrid.h
Go to the documentation of this file.
1 #ifndef _CompositeGrid
2 #define _CompositeGrid
3 
4 //
5 // Who to blame: Geoff Chesshire
6 //
7 
8 #include "GridCollection.h"
10 #ifdef USE_STL
11 #include "RVector.h"
12 #else
14 #include "ListOfGridCollection.h"
15 #include "ListOfCompositeGrid.h"
16 #include "ListOfIntegerArray.h"
17 #include "ListOfRealArray.h"
18 #endif // USE_STL
19 #include "TrivialArray.h"
20 #include "BoundaryAdjustment.h"
22 
23 class CompositeGrid;
24 
25 //
26 // Class for reference-counted data.
27 //
29  public GridCollectionData {
30  public:
31 
34  enum {
51  | THEinverseMap,
53  ISgivenByInterpoleePoint = INT_MIN // (sign bit) bit 31
54  };
57 
58  IntegerArray numberOfInterpolationPoints; // total number of interp. pts
59  IntegerArray numberOfImplicitInterpolationPoints; // for grids with mixed implicit/explicit
60  IntegerArray interpolationStartEndIndex; // for interp. pts ordered by interpolee grid
61 
62  // IntegerArray numberOfInterpoleePoints; // **** what is this?
63 
67 
68  IntegerArray interpolationWidth; // *wdh* move to GridCollection
70  RealArray maximumHoleCuttingDistance; // maximumHoleCuttingDistance(side,axis,grid)
71 // RealArray interpolationConditionLimit;
79 #ifdef USE_STL
80  RVector<RealDistributedArray> interpolationCoordinates;
81  RVector<IntDistributedArray> interpoleeGrid;
82  RVector<IntDistributedArray> variableInterpolationWidth;
83  RVector<IntDistributedrArray> interpoleeLocation;
84  RVector<IntDistributedArray> interpolationPoint;
85  RCVector<CompositeGrid> multigridLevel; // (overloaded)
86  RCVector<CompositeGrid> domain; // (overloaded)
87 #else
94  ListOfCompositeGrid domain; // (overloaded)
95 #endif // USE_STL
99  // IntegerArray sidesShare; // keeps track of shared sides.
100 
101  // Here are the new interpolation data arrays that are serial arrays -- these hold the
102  // interpolation data for the interpolation points on a given processor
103 
105  {
109  };
110  LocalInterpolationDataEnum localInterpolationDataState; // inidicates which grids use local arrays
111 
118  IntegerArray interpolationStartEndIndexLocal; // for interp. pts ordered by interpolee grid
119 
120  Partitioning_Type partition; // for interpolation arrays 031016
121  bool partitionInitialized; // true when the partition has been specified.
122 
123  //
124  // hybrid structured - unstructured connectivity data (should be private ??)
125  //
127 
128  //
129  // surface stitching mesh
130  //
132 
134  const Integer numberOfDimensions_ = 0,
135  const Integer numberOfComponentGrids_ = 0);
137  const CompositeGridData& x,
138  const CopyType ct = DEEP);
139  virtual ~CompositeGridData();
141  void reference(const CompositeGridData& x);
142  virtual void breakReference();
143  virtual void consistencyCheck() const;
144  virtual Integer get(
145  const GenericDataBase& db,
146  const aString& name);
147  virtual Integer put(GenericDataBase& db,
148  const aString& name,
149  int geometryToPut = -1 // by default put computedGeometry
150  ) const;
151  inline Integer update(
152  const Integer what = THEusualSuspects,
153  const Integer how = COMPUTEtheUsual)
154  { return update(*this, what, how); }
155  inline Integer update(
157  const Integer what = THEusualSuspects,
158  const Integer how = COMPUTEtheUsual)
159  { return update((GenericGridCollectionData&)x, what, how); }
160  virtual void destroy(const Integer what = NOTHING);
161  virtual Integer addRefinement(
162  const IntegerArray& range,
163  const IntegerArray& factor,
164  const Integer& level,
165  const Integer k = 0);
166  // replace refinement level level0 and higher
167  virtual Integer replaceRefinementLevels(int level0, int numberOfRefinementLevels0, IntegerArray **gridInfo );
168  virtual void deleteRefinement(const Integer& k);
169  virtual void deleteRefinementLevels(const Integer level = 0);
172  const Integer level = INTEGER_MAX)
175  const IntegerArray& factor,
176  const Integer& level,
177  const Integer k = 0);
178  virtual void makeCompleteMultigridLevels();
179  virtual void deleteMultigridCoarsening(const Integer& k);
180  virtual void deleteMultigridLevels(const Integer level = 0);
182  const Integer& k1,
183  const Integer& k2,
184  const RealArray& r,
185  const IntegerArray& interpolationStencil,
186  const intArray& useBackupRules);
188  const MappedGrid& g,
189  const Integer& k10,
190  const Integer& k20,
191  const RealArray& r,
192  const IntegerArray& interpolationStencil,
193  const intArray& useBackupRules);
194 
196  const Integer& k1,
197  const Integer& k2,
198  const realArray& r,
199  const intArray& ok,
200  const intArray& useBackupRules,
201  const Logical checkForOneSided = LogicalFalse);
202 // Logical canInterpolate(
203 // const MappedGrid& g,
204 // CompositeMask& g_mask,
205 // const Integer& k10,
206 // const Integer& k20,
207 // const RealArray& r,
208 // const LogicalArray& ok,
209 // const LogicalArray& useBackupRules,
210 // const Logical checkForOneSided);
212  const Integer& k1,
213  const Integer& k2,
214  const IntegerArray& i1,
215  const RealArray& r2,
216  const LogicalArray& ok);
217 // void adjustBoundary(
218 // const Integer& k1,
219 // const Integer& k2,
220 // const IntegerArray& i1,
221 // const RealArray& x);
222  virtual void setNumberOfGrids(const Integer& numberOfGrids_);
223  virtual void setNumberOfDimensions(const Integer& numberOfDimensions_);
224  virtual void setNumberOfDimensionsAndGrids(
225  const Integer& numberOfDimensions_,
226  const Integer& numberOfGrids_);
227  void initialize(
228  const Integer& numberOfDimensions_,
229  const Integer& numberOfGrids_);
230 
231  protected:
232 
233  // convert serial array (local) interpolation data to parallel array form
235 
236  virtual Integer update(
238  const Integer what = THEusualSuspects,
239  const Integer how = COMPUTEtheUsual);
240  virtual void referenceRefinementLevels(
242  const Integer level = INTEGER_MAX);
244  const Integer& what,
245  Integer& numberOfCollections,
246 #ifdef USE_STL
247  RCVector<CompositeGrid>& list,
248  RCVector<GridCollection>& gridCollectionList,
249  RCVector<GenericGridCollection>& genericGridCollectionList,
250 #else
251  ListOfCompositeGrid& list,
252  ListOfGridCollection& gridCollectionList,
253  ListOfGenericGridCollection& genericGridCollectionList,
254 #endif // USE_STL
255  IntegerArray& number);
256 //
257 // The following functions are declared protected here in order to disallow
258 // access to the corresponding functions in class GenericGridCollectionData.
259 //
260  inline virtual Integer addRefinement(
261  const Integer& level,
262  const Integer k = 0)
263  { return GridCollectionData::addRefinement(level, k); }
265  const Integer& level,
266  const Integer k = 0)
267  { return GridCollectionData::addMultigridCoarsening(level, k); }
268 
269  void initializePartition();
270 
271  private:
272  Integer computeGeometry(
273  const Integer& what,
274  const Integer& how);
275 //
276 // Virtual member functions used only through class ReferenceCounting:
277 //
278  private:
279  inline virtual ReferenceCounting& operator=(const ReferenceCounting& x)
280  { return operator=((CompositeGridData&)x); }
281  inline virtual void reference(const ReferenceCounting& x)
282  { reference((CompositeGridData&)x); }
283  inline virtual ReferenceCounting* virtualConstructor(
284  const CopyType ct = DEEP) const
285  { return new CompositeGridData(*this, ct); }
286  aString className;
287  public:
288  inline virtual aString getClassName() const { return className; }
289 };
290 
292  public GridCollection {
293  public:
294 // Public constants:
295 
296 // Constants to be ORed to form the first argument of update() and destroy():
297  enum {
307  };
308 
309 // Constants to be ORed to form the second argument of update():
310  enum {
312  };
313 
314 // Public data:
315 //
316 // The number of interpolation points on each component grid.
317 //
319  IntegerArray numberOfImplicitInterpolationPoints; // for grids with mixed implicit/explicit
320  IntegerArray interpolationStartEndIndex; // for interp. pts ordered by interpolee grid
321 //
322 // The number of interpolee points listed separately on each component grid.
323 //
324 // IntegerArray numberOfInterpoleePoints;
325 //
326 // The type of interpolation (toGrid, fromGrid).
327 //
329 // LogicalArray backupInterpolationIsImplicit;
330 //
331 // The width of the interpolation stencil (direction, toGrid, fromGrid).
332 //
333  IntegerArray interpolationWidth; // *wdh* moved to GridCollection
334 // IntegerArray backupInterpolationWidth;
335 //
336 // The minimum overlap for interpolation (direction, toGrid, fromGrid).
337 //
339 
340  // maximumHoleCuttingDistance(side,axis,grid) : cut hole points at most this far away from the face.
342 
343 // RealArray backupInterpolationOverlap;
344 //
345 // The maximum condition number allowed for interpolation.
346 //
347 // RealArray interpolationConditionLimit;
348 // RealArray backupInterpolationConditionLimit;
349 //
350 // The list of preferences for interpolation of each grid.
351 //
353 //
354 // Indicates which grids may interpolate from which.
355 //
357 // LogicalArray mayBackupInterpolate;
358 //
359 // Indicates which grids may cut holes in other grids.
360 //
362 //
363 // By default shared sides do not cut holes
364 //
366 //
367 // Multigrid Coarsening ratio.
368 //
370 //
371 // Multigrid Prolongation stencil width.
372 //
374 //
375 // Multigrid Restriction stencil width.
376 //
378 //
379 // Range of points interpolated from each interpolee grid.
380 //
381 // IntegerArray interpoleeGridRange;
382 #ifdef USE_STL
383 //
384 // Coordinates of interpolation point.
385 //
386  RVector<RealDistributedArray> interpolationCoordinates;
387 //
388 // Index of interpolee grid.
389 //
390  RVector<IntDistributedArray> interpoleeGrid;
391 
392 // width of interpolation (if not constant)
393  RVector<IntDistributedArray> variableInterpolationWidth;
394 //
395 // Index of the grid number of each separately-listed interpolee point.
396 //
397 // RVector<IntegerArray> interpoleePoint;
398 //
399 // Location of interpolation stencil.
400 //
401  RVector<IntDistributedArray> interpoleeLocation;
402 //
403 // Indices of interpolation point.
404 //
405  RVector<IntDistributedArray> interpolationPoint;
406 //
407 // Interpolation condition number.
408 //
409 // RVector<RealArray> interpolationCondition;
410 //
411 // Collections of grids having the same multigrid level.
412 //
413  RCVector<CompositeGrid> multigridLevel; // (overloaded)
414 
415 // Collections of grids in the same domain
416  RCVector<CompositeGrid> domain; // (overloaded)
417 #else
418 //
419 // Coordinates of interpolation point.
420 //
422 //
423 // Index of interpolee grid.
424 //
426 // width of interpolation (if not constant)
428 //
429 // Index of the grid number of each separately-listed interpolee point.
430 //
431 // ListOfIntegerArray interpoleePoint;
432 //
433 // Location of interpolation stencil.
434 //
436 //
437 // Indices of interpolation point.
438 //
440 //
441 // Interpolation condition number.
442 //
443 // ListOfRealArray interpolationCondition;
444 //
445 // Collections of grids having the same multigrid level.
446 //
448 
449 // Collections of grids in the same domain
450  ListOfCompositeGrid domain; // (overloaded)
451 
452 #endif // USE_STL
453 
454 // Data used by class Cgsh for optimization of overlap computation:
455 // RealCompositeGridFunction inverseCondition; // Inverse quality.
457  IntegerCompositeGridFunction inverseGrid; // Inverse component grid.
458 //
459 // Public member functions for access to data:
460 //
461 // The number of component grids.
462 //
467 //
468 // The tolerance for accuracy of interpolation coordinates.
469 //
470  inline Real& epsilon() { return rcData->epsilon; }
471  inline const Real& epsilon() const { return rcData->epsilon; }
472 //
473 // The type of interpolation between all pairs of grids.
474 //
477  inline const Logical& interpolationIsAllExplicit() const
481  inline const Logical& interpolationIsAllImplicit() const
483 //
484 // Public member functions:
485 //
486 // Default constructor.
487 //
488 // If numberOfDimensions_==0 (e.g., by default) then create a null
489 // CompositeGrid. Otherwise create a CompositeGrid
490 // with the given number of dimensions and number of component grids.
491 //
493  const Integer numberOfDimensions_ = 0,
494  const Integer numberOfComponentGrids_ = 0);
495 //
496 // Copy constructor. (Does deep copy by default.)
497 //
499  const CompositeGrid& x,
500  const CopyType ct = DEEP);
501 //
502 // Destructor.
503 //
504  virtual ~CompositeGrid();
505 //
506 // Assignment operator. (Does a deep copy.)
507 //
509 //
510 // Make a reference. (Does a shallow copy.)
511 //
512  void reference(const CompositeGrid& x);
513  void reference(CompositeGridData& x);
514 //
515 // Break a reference. (Replaces with a deep copy.)
516 //
517  virtual void breakReference();
518 //
519 // Change the grid to be all vertex-centered.
520 //
521  virtual void changeToAllVertexCentered();
522 //
523 // Change the grid to be all cell-centered.
524 //
525  virtual void changeToAllCellCentered();
526 
527  // reduce interpolation width of an already computed overlapping grid.
528  int changeInterpolationWidth( int width );
529 //
530 // Check that the data structure is self-consistent.
531 //
532  virtual void consistencyCheck() const;
533 //
534 // "Get" and "put" database operations.
535 //
536  virtual Integer get(
537  const GenericDataBase& db,
538  const aString& name);
539  virtual Integer put(GenericDataBase& db,
540  const aString& name,
541  int geometryToPut = -1 // by default put computedGeometry
542  ) const;
543 //
544 // Set references to reference-counted data.
545 //
546  void updateReferences(const Integer what = EVERYTHING);
547 //
548 // Update the grid.
549 //
550  inline Integer update(
551  const Integer what = THEusualSuspects,
552  const Integer how = COMPUTEtheUsual)
553  { return update(*this, what, how); }
554 //
555 // Update the grid, sharing the data of another grid.
556 //
557  inline Integer update(
558  CompositeGrid& x,
559  const Integer what = THEusualSuspects,
560  const Integer how = COMPUTEtheUsual)
561  { return update((GenericGridCollection&)x, what, how); }
562 //
563 // Destroy optional grid data.
564 //
565  virtual void destroy(const Integer what = NOTHING);
566 
567 // Add a new grid, built from a Mapping
568  virtual int add(Mapping & map);
569 
570 // Add a new grid
571  virtual int add(MappedGrid & g);
572 
573  // delete a grid:
574  virtual int deleteGrid(Integer k);
575 
576  // delete a list of grids:
577  virtual int deleteGrid(const IntegerArray & gridsToDelete );
578 
579  // delete a list of grids (use this when you also know the list of grids to save).
580  virtual int deleteGrid(const IntegerArray & gridsToDelete, const IntegerArray & gridsToSave );
581 
582 //
583 // Add a refinement grid to the collection.
584 //
585  virtual Integer addRefinement(
586  const IntegerArray& range, // The indexRange of the refinement grid.
587  const IntegerArray& factor, // The refinement factor w.r.t. level-1.
588  const Integer& level, // The refinement level number.
589  const Integer k = 0); // The index of an ancestor of the refinement.
591  const IntegerArray& range,
592  const Integer& factor,
593  const Integer& level,
594  const Integer k = 0) {
595  IntegerArray factors(3); factors = factor;
596  return addRefinement(range, factors, level, k);
597  }
598  // replace refinement level level0 and higher
599  virtual Integer replaceRefinementLevels(int level0, int numberOfRefinementLevels0, IntegerArray **gridInfo );
600 //
601 // Delete all multigrid levels of refinement grid k.
602 //
603  virtual void deleteRefinement(const Integer& k);
604 //
605 // Delete all grids with refinement level greater than the given level.
606 //
607  virtual void deleteRefinementLevels(const Integer level = 0);
608 //
609 // Reference x[i] for refinementLevelNumber(i) <= level.
610 // Delete all other grids.
611 //
613  CompositeGrid& x,
614  const Integer level = INTEGER_MAX)
616 
617 
618 // Return the number of possible multigrid levels supported by the current numbers of grid lines
619 // on the different grids.
621 
622 //
623 // Add a multigrid coarsening of grid k.
624 //
626  const IntegerArray& factor, // The coarsening factor w.r.t level-1
627  const Integer& level, // The multigrid level number.
628  const Integer k = 0); // The index of the corresponding grid
629  // at any finer multigrid level.
630 //
631 // Add multigrid coarsenings of grids in order to complete the multigrid levels.
632 //
633  virtual void makeCompleteMultigridLevels();
634 //
635 // Delete grid k, a multigrid coarsening, and all of its multigrid coarsenings.
636 //
637  virtual void deleteMultigridCoarsening(const Integer& k);
638 //
639 // Delete all of the grids with multigrid level greater than the given level.
640 //
641  virtual void deleteMultigridLevels(const Integer level = 0);
642 
643  // Save the grid in a file
644  int saveGridToAFile(const aString & gridFileName, const aString & gridName );
645 
646 //
647 // Set the number of grids.
648 //
649  virtual void setNumberOfGrids(const Integer& numberOfGrids_);
650 //
651 // Set the number of dimensions.
652 //
653  virtual void setNumberOfDimensions(const Integer& numberOfDimensions_);
654 //
655 // Set the number of dimensions and grids.
656 //
657  virtual void setNumberOfDimensionsAndGrids(
658  const Integer& numberOfDimensions_,
659  const Integer& numberOfGrids_);
660 
661  // return size of this object
662  virtual real sizeOf(FILE *file = NULL, bool returnSizeOfReference = false ) const;
663 
664  // assign overlap parameters such as interpolationWidth etc.
665  int setOverlapParameters();
666  int setOverlapParameters(CompositeGrid & cg); // base parameters on this CompositeGrid.
667 
668 //
669 // Specify the set of processes over which GridFunctions are distributed.
670 // We now support only the specification of a contiguous range of process IDs.
671 //
672  void specifyProcesses(const Range& range);
673 
674 //
675 // Initialize the CompositeGrid with the given number of dimensions and grids.
676 // These grids have their gridNumbers, baseGridNumbers and componentGridNumbers
677 // set to [0, ..., numberOfGrids_-1], and their refinementLevelNumbers and
678 // multigridLevelNumbers set to zero.
679 //
680  virtual void initialize(
681  const Integer& numberOfDimensions_,
682  const Integer& numberOfGrids_);
683 //
684 // Return the bounds on the cube of points used for interpolation.
685 //
686 // Input:
687 // k1 Index of the grid containing the interpolation points.
688 // k2 Index of the grid containing the interpolee points.
689 // r(0:2,p1:p2) Interpolation coordinates of the interpolation points.
690 // The dimensions of this array determine p1 and p2.
691 // useBackupRules(p1:p2)
692 // Must be LogicalFalse unless backup interpolation
693 // rules should be used to determine the stencil.
694 // Output:
695 // interpolationStencil(0:1,0:2,p1:p2)
696 // lower and upper index bounds of the stencils.
697 //
699  const Integer& k1,
700  const Integer& k2,
701  const RealArray& r,
702  const IntegerArray& interpolationStencil,
703  const intArray& useBackupRules) {
704  rcData->getInterpolationStencil(k1, k2, r, interpolationStencil,
705  useBackupRules);
706  }
708  const MappedGrid& g, // The unrefined grid corresponding to grid[k20].
709  const Integer& k1,
710  const Integer& k2,
711  const RealArray& r,
712  const IntegerArray& interpolationStencil,
713  const intArray& useBackupRules) {
714  rcData->getInterpolationStencil(g, k1, k2, r, interpolationStencil,
715  useBackupRules);
716  }
717 //
718 // Determine whether points on grid k1 at r in the coordinates of grids k2
719 // can be interpolated from grids k2.
720 //
722  const Integer& k1,
723  const Integer& k2,
724  const realArray& r, // (0:numberOfDimensions-1,p1:p2)
725  const intArray& ok, // (p1:p2)
726  const intArray& useBackupRules, // (p1:p2)
727  const Logical checkForOneSided = LogicalFalse) {
728  return rcData->canInterpolate(k1, k2, r, ok,
729  useBackupRules, checkForOneSided);
730  }
731 // inline Logical canInterpolate(
732 // const MappedGrid& g, // The unrefined grid corresponding to grid[k20].
733 // CompositeMask& g_mask, // Masks on k2 (possibly) and its siblings.
734 // const Integer& k1,
735 // const Integer& k2,
736 // const RealArray& r, // (0:numberOfDimensions-1,p1:p2)
737 // const LogicalArray& ok, // (p1:p2)
738 // const LogicalArray& useBackupRules, // (p1:p2)
739 // const Logical checkForOneSided = LogicalFalse) {
740 // return rcData->canInterpolate(g, g_mask, k1, k2, r, ok,
741 // useBackupRules, checkForOneSided);
742 // }
743 //
744 // Check if these boundary discretization points of this grid
745 // lie at least epsilon inside the parameter space of grid g2.
746 //
747 // void isInteriorBoundaryPoint(
748 // const Integer& k1,
749 // const Integer& k2,
750 // const IntegerArray& i1,
751 // const RealArray& r2,
752 // const LogicalArray& ok)
753 // { rcData->isInteriorBoundaryPoint(k1, k2, i1, r2, ok); }
754 //
755 // Adjust the position of interpolation points to take
756 // into account mismatch between shared boundaries.
757 //
758 // inline void adjustBoundary(
759 // const Integer& k1,
760 // const Integer& k2,
761 // const IntegerArray& i1,
762 // const RealArray& x)
763 // { rcData->adjustBoundary(k1, k2, i1, x); }
764 
765  // Here is the master grid.
767 //
768 // Pointer to reference-counted data.
769 //
772  inline CompositeGridData* operator->() { return rcData; }
773  inline const CompositeGridData* operator->() const { return rcData; }
774  inline CompositeGridData& operator*() { return *rcData; }
775  inline const CompositeGridData& operator*() const { return *rcData; }
776 
777  //
778  // set the hybridConnectivity for this CompositeGrid
779  //
780  void setHybridConnectivity(const int grid,
781  intArray * gridIndex2UVertex_,
782  intArray & uVertex2GridIndex_,
783  intArray * gridVertex2UVertex_, // could be build from gridIndex2UVertex ?
784  intArray & boundaryFaceMapping_);
785 
786  //
787  // get the hybrid connectivity for this composite grid
788  //
790 
791  void setSurfaceStitching( UnstructuredMapping *stitching );
792 
794 
796 
797  private:
798  virtual Integer update(
800  const Integer what = THEusualSuspects,
801  const Integer how = COMPUTEtheUsual);
802  virtual void referenceRefinementLevels(
804  const Integer level = INTEGER_MAX);
805  protected:
806 
807  CompositeGrid *master; // overloaded
808 
809 //
810 // The following functions are declared protected here in order to disallow
811 // access to the corresponding functions in class GenericGridCollectionData.
812 //
813  inline virtual Integer addRefinement(
814  const Integer& level,
815  const Integer k = 0)
816  { return GridCollection::addRefinement(level, k); }
818  const Integer& level,
819  const Integer k = 0)
820  { return GridCollection::addMultigridCoarsening(level, k); }
821  inline virtual void initialize(
822  const Integer& numberOfGrids_)
823  { GridCollection::initialize(numberOfGrids_); }
824 //
825 // Virtual member functions used only through class ReferenceCounting:
826 //
827  private:
828  inline virtual ReferenceCounting& operator=(const ReferenceCounting& x)
829  { return operator=((CompositeGrid&)x); }
830  inline virtual void reference(const ReferenceCounting& x)
831  { reference((CompositeGrid&)x); }
832  inline virtual ReferenceCounting* virtualConstructor(
833  const CopyType ct = DEEP) const
834  { return new CompositeGrid(*this, ct); }
835  aString className;
836  public:
837  inline virtual aString getClassName() const { return className; }
838 };
839 //
840 // Stream output operator.
841 //
842 ostream& operator<<(ostream& s, const CompositeGrid& g);
843 
844 #endif // _CompositeGrid