Overture  Version 25
Interpolant.h
Go to the documentation of this file.
1 #ifndef INTERPOLANT_H
2 #define INTERPOLANT_H "Interpolant.h"
3 
4 //============================================================================
5 // Interpolant Class
6 // Interpolate and Swap Periodic Edges
7 // -----------------------------------
8 // Supply implicit and explicit interpolation operators for
9 // realCompositeGridFunctions. By default the interpolation
10 // operator will swap periodic edges.
11 //
12 // Notes:
13 // o A CompositeGrid holds a pointer to an Interpolant. This pointer is set by an
14 // Interpolant when the constructor Interpolant(CompositeGrid &) or the member function
15 // updateToMatchGrid( CompositeGrid & ) is called. Grid functions look for
16 // the pointer in the CompositeGrid when their interpolate function is called.
17 // o The first time that an Interpolant is created the pointer in the CompositeGrid
18 // will be set. Subsequently, when new Interpolant's are created they will
19 // simply reference the existing one. Thus the CompositeGrid will point to an
20 // Interpolant (and grid functions will know how to interpolate) as long as at
21 // least one Interpolant remains in scope.
22 // o The explicit interpolation is defined in Interpolant.C
23 // o The implicit interpolation is peformed by Oges - the overlapping grid
24 // equation solver.
25 //
26 // Usage:
27 // CompositeGrid cg(...);
28 // realCompositeGridFunction u(cg);
29 // ...
30 // Interpolant interpolant(cg);
31 // ...
32 // interpolant.interpolate( u ) ;
33 // ..or..
34 // u.interpolate(); // implicitly knows about the interpolant through the CompositeGrid
35 //
36 //============================================================================
37 
38 #include "CompositeGrid.h"
39 #include "CompositeGridFunction.h"
40 // include "MultigridCompositeGridFunction.h"
41 #include "ListOfRealArray.h"
42 #include "ListOfDoubleArray.h"
43 
44 class Oges; // forward declaration
45 class InterpolateRefinements; // forward declaration
46 class ParallelOverlappingGridInterpolator; // forward declaration
47 
49 {
50  public:
51 
52  Interpolant();
55  Interpolant(const Interpolant & interpolant, const CopyType copyType=DEEP );
56  virtual ~Interpolant();
57  Interpolant & operator= ( const Interpolant & interpolant );
58  void reference( const Interpolant & interpolant );
59  virtual void breakReference();
60 
61  void updateToMatchGrid(CompositeGrid & cg, int refinementLevel=0 );
62 
64  const Range & C0 = nullRange, // optionally specify components to interpolate
65  const Range & C1 = nullRange,
66  const Range & C2 = nullRange );
67 
69  const Range & C0 = nullRange,
70  const Range & C1 = nullRange,
71  const Range & C2 = nullRange );
72 
73  int interpolate( int gridToInterpolate, // only interpolate this grid.
75  const Range & C0 = nullRange, // optionally specify components to interpolate
76  const Range & C1 = nullRange,
77  const Range & C2 = nullRange );
78 
80  const IntegerArray & gridsToInterpolate, // specify which grids to interpolate
81  const Range & C0 = nullRange, // optionally specify components to interpolate
82  const Range & C1 = nullRange,
83  const Range & C2 = nullRange );
84 
86  const IntegerArray & gridsToInterpolate, // specify which grids to interpolate
87  const IntegerArray & gridsToInterpolateFrom, // specify which grids to interpolate from
88  const Range & C0 = nullRange, // optionally specify components to interpolate
89  const Range & C1 = nullRange,
90  const Range & C2 = nullRange );
91 
92  int interpolate( realArray & ui, // save results here
93  int gridToInterpolate, // only interpolate values on this grid that
94  int interpoleeGrid, // interpolate from this grid.
96  const Range & C0 = nullRange, // optionally specify components to interpolate
97  const Range & C1 = nullRange,
98  const Range & C2 = nullRange );
99 
100  bool interpolationIsExplicit() const;
101  bool interpolationIsImplicit() const;
102 
104  {
107  optimizedC, // use C style loops
108  numberOfInterpolationMethods // counts number in this list
109  };
110 
112 
113  int setMaximumNumberOfIterations(int maximumNumberOfIterations); // for iterative interpolation
114 
116  {
117  precomputeAllCoefficients, // requires w^d coefficients per interp pt (w=width of interp stencil)
118  precomputeSomeCoefficients, // requires w*d coefficients per interp pt (d=dimension, 1,2, or 3)
119  precomputeNoCoefficients // requires d coefficinets per interp point
120  };
121  // For wider interpolation stencils we may want to use less storage
123 
125  {
129  };
130 
134 
136  {
140  };
141 
142  int setInterpolationOption(InterpolationOptionEnum option, bool trueOrFalse );
144 
145  // supply an AMR interpolation object:
147  int setMaximumRefinementLevelToInterpolate(int maxLevelToInterpolate );
149 
150  // return size of this object
151  virtual real sizeOf(FILE *file = NULL ) const;
152 
153  int static printStatistics( FILE *file= stdout ); // statistics for all Interpolants
154 
155  int printMyStatistics( FILE *file= stdout ); // statistics for this Interpolant too
156 
157  int interpolateRefinementLevel( const int refinementLevel,
159  const Range & C0 = nullRange, // optionally specify components to interpolate
160  const Range & C1 = nullRange,
161  const Range & C2 = nullRange );
162 
163  BoundaryConditionParameters bcParams; // used to apply BC's for AMR grids.
164 
165  // for testing interpolation and amr interpolation
166  static int testInterpolation( CompositeGrid & cg, int problemType );
167 
168  int debug;
169 
170  protected:
171 
172  int getComponentRanges(const Range & C0, const Range & C1, const Range & C2, Range C[4],
174 
175  private:
176 
177  InterpolationMethodEnum interpolationMethod;
178  ImplicitInterpolationMethodEnum implicitInterpolationMethod;
179  ExplicitInterpolationStorageOptionEnum explicitInterpolationStorageOption;
180  int *useVariableWidthInterpolation;
181 
182  int maximumNumberOfIterations;
183 
184  real tolerance; // tolerance for implicit interpolation by iteration.
185  bool explicitInterpolation;
186  bool interpolationIsInitialized;
187  bool initializeParallelInterpolator;
188 
189  bool interpolateRefinementBoundaries; // if true, interpolate all refinement boundaries
190  bool interpolateHidden; // if true, interpolate hidden coarse grid points from higher level refinemnts
191  bool interpolateOverlappingRefinementBoundaries;
192  bool interpRefinementsWasNewed;
193  int maximumRefinementLevelToInterpolate; // only interpolate levels <= this level
194 
195  ListOfRealDistributedArray coeff; // coeff's for explicit interpolation
196  IntegerArray width;
197  CompositeGrid cg;
198 
199  // realCompositeGridFunction v; // holds components, no need to reference count
200  static real timeForExplicitInterpolation;
201  static real timeForImplicitInterpolation;
202  static real timeForIterativeImplicitInterpolation;
203  static real timeForInitializeInterpolation;
204  static real timeForAMRInterpolation;
205  static real timeForAMRCoarseFromFine;
206  static real timeForAMRExtrapolateRefinementBoundaries;
207  static real timeForAMRExtrapolateAll;
208  static real timeForAMRExtrapInterpolationNeighbours;
209  static real timeForAMRRefinementBoundaries;
210 
211 
212  int updateForAdaptiveGrid;
213  InterpolateRefinements *interpRefinements;
214 
215  void initialize();
216  int initializeInterpolation();
217  int initializeExplicitInterpolation();
218 
219  int internalInterpolate( realCompositeGridFunction & u,
220  const Range C[],
221  const IntegerArray & gridToInterpolate = Overture::nullIntArray(),
222  const IntegerArray & gridsToInterpolateFrom = Overture::nullIntArray());
223 
224  int explicitInterpolate(realCompositeGridFunction & u,
225  const Range C[],
226  const IntegerArray & gridsToInterpolate = Overture::nullIntArray(),
227  const IntegerArray & gridsToInterpolateFrom = Overture::nullIntArray() ) const;
228 
229  int implicitInterpolateByIteration(realCompositeGridFunction & u,
230  const Range C[],
231  const IntegerArray & gridToInterpolate = Overture::nullIntArray(),
232  const IntegerArray & gridsToInterpolateFrom = Overture::nullIntArray() ) const;
233 
234  virtual ReferenceCounting& operator=( const ReferenceCounting & x)
235  { return operator=( *(Interpolant*) & x ); }
236  virtual void reference( const ReferenceCounting & x)
237  { reference( (Interpolant &) x ); }
238  virtual ReferenceCounting* virtualConstructor(const CopyType ct = DEEP) const
239  { return ::new Interpolant(*this, ct); }
240 
241 
242  // holds data to be reference counted (not A++ since they are ref. counted)
243  class RCData : public ReferenceCounting
244  {
245  public:
246  friend class GridCollectionFunction;
247  RCData();
248  ~RCData();
249  RCData & operator=( const RCData & rcdata );
250 
251  Oges *implicitInterpolant;
252  ParallelOverlappingGridInterpolator *parallelInterpolator; // holds a pointer to the parallel interpolator
253 
254  private:
255 
256  // These are used by list's of ReferenceCounting objects
257  virtual void reference( const ReferenceCounting & mgf )
258  { RCData::reference( (RCData&) mgf ); }
259  virtual ReferenceCounting & operator=( const ReferenceCounting & mgf )
260  { return RCData::operator=( (RCData&) mgf ); }
261  virtual ReferenceCounting* virtualConstructor( const CopyType = DEEP ) const
262  { return ::new RCData(); }
263  };
264 
265 public: // TEMP
266  RCData *rcData;
267 };
268 
269 
270 
271 #endif