Overture  Version 25
OgmgParameters.h
Go to the documentation of this file.
1 #ifndef OGMG_PARAMETERS_H
2 #define OGMG_PARAMETERS_H
3 
4 // This class holds parameters used by Ogmg.
5 // You can save commonly used parameter values in an object of this class.
6 // You can interactively update parameter values.
7 // To make these parameters known to a particular Ogmg object, use the
8 // Ogmg::setOgmgParameters function.
9 
10 #include "Overture.h"
11 #include "OgesParameters.h"
12 
13 // forward declarations:
14 class Ogmg;
16 class DialogData;
17 class LoadBalancer;
18 
19 
21 {
22  public:
23 
24  enum OptionEnum // ** remember to change the documentation if you change this list.
25  {
26  THEnumberOfCycles, // number of times to iterate on each level (=1 -> V-cycle, 2=W-cycle)
27  THEnumberOfPreSmooths, // (minimum) number of pre-smooths (level)
28  THEnumberOfPostSmooths, // (minimum) number of postsmooths (level)
29  THEnumberOfSmooths, // will go away
30  THEnumberOfSubSmooths, // number of sub smooths (grid,level)
31  THEsmootherType, // type of smooth (grid,level)
32  THEsmoothingRateCutoff, // continue smoothing until smoothing rate is bigger than this
33  THEuseDirectSolverOnCoarseGrid, // if true use Oges, if false use a 'smoother' on the coarse grid.
47  };
48 
50  {
51  cycleTypeC, // normal "gamma" cycle, gamma=1: V, gamma=2 : W
52  cycleTypeF, // F cycle
53  cycleAdaptive // adaptive cycle (not implemented yet)
54  };
55 
56 
57  // Here are smoothers:
59  {
60  Jacobi=0,
62  redBlack, // red-black Gauss-Seidel
74  };
75 
76 
77  // Convergence can be determined in different ways: *wdh* 110310
79  {
80  residualConverged, // l2Norm(residual) < residualTolerance*l2Nnorm(f) + absoluteTolerance
81  errorEstimateConverged, // max(error estimate) < errorTolerance
82  residualConvergedOldWay // max(residual) < relativeTolerance*numberOfGridPoints (old way - keep as an option for now)
83  };
84 
85 
86  // these are the numbers defining values in the cg.boundaryCondition(side,axis) array
87  // **These should match those in OgesParameters.h**
89  {
92  mixed=3,
102  };
103 
105  {
110  };
111 
112 
113  enum
114  {
115  allGrids=-99,
116  allLevels=-100,
118  };
119 
120  OgmgParameters();
121  OgmgParameters(CompositeGrid & cg); // apply parameters to this grid.
122  ~OgmgParameters();
123 
124  virtual OgmgParameters& operator=(const OgmgParameters& par);
125 
126  // Automatically choose good parameters (smoothers etc.) for a grid
127  int chooseGoodMultigridParameters(CompositeGrid & cg, int maxLevels=useLevelsInGrid, int robustnessLevel=1 );
128 
129  // print out current values of parameters
130  int display(FILE *file = stdout);
131 
132  int get( OptionEnum option, int & value ) const;
133  int get( OptionEnum option, real & value ) const;
134 
135  int get( const GenericDataBase & dir, const aString & name);
136 
137  int put( GenericDataBase & dir, const aString & name) const;
138 
139  int set( CompositeGrid & cg); // apply parameters to this grid.
140 
141  // aString getSolverName();
142 
143  // ----------------Functions to set parameters -----------------------------
144  int set( OptionEnum option, int value=0 );
145  int set( OptionEnum option, float value );
146  int set( OptionEnum option, double value );
147 
148  // Set the absolute tolerance for the residual based convergence criteria.
150 
151  // Set the tolerance for the error, iterate until the norm of the estimated error is
152  // less than "errorTolerance"
154 
155  int setMaximumNumberOfIterations( const int max );
156 
157  // set the mean value of the solution for a singular problem
158  int setMeanValueForSingularProblem( const real meanValue );
159 
160  // set option concerning the null vector for singular problems
161  int setNullVectorOption( NullVectorOptionsEnum option, const aString & fileName );
162 
163  // specify the number of iterations per level (1=V cycle, 2=W cycle)
164  int setNumberOfCycles( const int & number, const int & level=allLevels );
165 
166  // Set the number of smooths on a given level
167  int setNumberOfSmooths(const int numberOfPreSmooths, const int numberOfPostSmooths, const int level);
168 
169  int setNumberOfSubSmooths( const int & numberOfSmooths, const int & grid, const int & level=allLevels);
170 
171  int setParameters( const Ogmg & ogmg); // set all parameters equal to those values found in ogmg.
172 
173  // Indicate if the problem is singular
174  int setProblemIsSingular( const bool trueOrFalse=TRUE );
175 
176  // Set the relative tolerance for the residual based convergence criteria (scaled by the l2Norm(f))
178 
179  int setSmootherType(const SmootherTypeEnum & smoother,
180  const int & grid=allGrids,
181  const int & level=allLevels );
182 
183  int updateToMatchGrid( CompositeGrid & cg, int maxLevels=useLevelsInGrid );
184  int update( GenericGraphicsInterface & gi, CompositeGrid & cg ); // update parameters interactively
185 
186  protected:
187 
188  int buildOptionsDialog( DialogData & dialog );
189 
190  int numberOfMultigridLevels() const; // this is how many levels the parameters thinks exists
191  int numberOfComponentGrids() const; // this is how many grids the parameters thinks exists
192 
193 
194  void init();
196 
197  int set( OptionEnum option, int value, real rvalue );
198  int get( OptionEnum option, int & value, real & rvalue ) const;
199 
201 
202  IntegerArray numberOfCycles; // number of times to iterate on each level (=1 -> V-cycle, 2=W-cycle)
203  IntegerArray numberOfSmooths; // (minimum) number of smooths (level)
204  IntegerArray numberOfSubSmooths; // number of sub smooths (grid,level)
205  IntegerArray smootherType; // type of smooth (grid,level)
206 
207  // Lower and upper bounds for the adaptive smoothing algorithm:
210 
211 
212  real omegaJacobi; // relaxation parameter for Jacobi smoother
213  real omegaGaussSeidel; // relaxation parameter for Gauss-Seidel smoother
214  real omegaRedBlack; // relaxation parameter for red-black smoother
217  real variableOmegaScaleFactor; // for scaling automatically chosen omega
218  bool useLocallyOptimalOmega; // use locally optimal omega on non-uniform grids
219  bool useLocallyOptimalLineOmega; // use locally optimal omega on non-uniform grids (line solvers)
222  bool useNewRedBlackSmoother; // this one works in parallel
223 
224  real smoothingRateCutoff; // continue smoothing until smoothing rate is bigger than this
225  bool useDirectSolverOnCoarseGrid; // if false use a 'smoother' on the coarse grid.
226  int numberOfIterationsOnCoarseGrid; // if iterating on the coarse grid
227  bool useDirectSolverForOneLevel; // if true call direct solver if there is only 1 level, otherwise use smoother
228 
231  int subSmoothReferenceGrid; // reference grid (uses 1 sub-smooth) for variable sub-smooths
232 
233  bool useFullMultigrid; // use a full MG cycle (i.e. start solution process from the coarse grid)
234  int minimumNumberOfInitialSmooths; // smooth at least this many times on the first iteration.
235 
242  OgesParameters *nullVectorParameters; // Oges parameters when solving for the null vector
243 
244 
246 
249  bool useErrorEstimateForConvergence; // if true check the error estimate in the convergence test
250 
251  OgesParameters ogesParameters; // for the coarse grid solver
253 
254  OgesParameters *ogesSmoothParameters; // for Oges solvers used as smoothers
256 
257  LoadBalancer *loadBalancer; // for parallel load balancing
258 
260 
263 
264  // we can read/save the multigrid composite grid instead of generating it:
267 
269  bool useNewAutoSubSmooth; // use defects computed earlier for determining defect ratios
270 
273  bool useOptimizedVersion; // if true use the new optimized smoothers etc.
274  bool decoupleCoarseGridEquations; // if true set interpolation points to zero on all coarser levels.
275 
277 
278  int fineToCoarseTransferWidth; // 1=injection, 3= full weighting
279  int coarseToFineTransferWidth; // 2=linear interpolation, 4=cubic interpolation
280 
281 
282  // boundary conditions are : extrapolation, equation, combination
283  // extrapolation: usually means a dirichlet BC on the boundary
284  // equation : usually means a neumann BC on the ghost line
285 
286  // possible conditions on the boundary or ghost line are
287  // imposeDirichlet : impose a dirichlet BC
288  // imposeExtrapolation : explicitly impose an extrapolation equation
289  // injection :
290  // partialWeighting : full weighting in the tangential directions
291  // halfWeighting : full weighting but exclude ghost points.
292  // lumpedPartialWeighting : lump partial weighting coefficients in tangential directions.
293  // imposeNeumann
295  {
303  };
304 
305  // These hold the default options for the boundary and ghostline
306  BoundaryAveragingOption boundaryAveragingOption[2]; // [0] : for extrapolation BC, [1] : for equation BC
308 
310  {
314  };
316 
317  int useSymmetryForNeumannOnLowerLevels; // replace Neumann BC by even symmetry for l>0
318  int useSymmetryForDirichletOnLowerLevels; // replace Neumann BC by even symmetry for l>0
322  int solveEquationWithBoundaryConditions; // for 4th order, solve PDE on boundary with BC's
323 
325  {
326  useUnknown=-1, // bogus value
331  };
332 
339 
346 
348 
350 
351  int numberOfBoundaryLayersToSmooth; // extra smoothing at boundaries
354 
355  bool combineSmoothsWithIBS; // if true then merge regular smooths with IBS smooths
356  int numberOfInterpolationLayersToSmooth; // for smoothing points near interpolation points
357  int numberOfInterpolationSmoothIterations; // for smoothing interpolation neighbours
359  int numberOfIBSIterations; // global iterations of interp. boundary smoothing
360 
361  int gridOrderingForSmooth; // 0= 1...ng 1= ng...1 2=alternate
363  IntegerArray totalNumberOfSmoothsPerLevel; // counts smooths per level
364  IntegerArray totalNumberOfSubSmooths; // counts sub-smooths per grid
365 
366  int coarseGridInterpolationWidth; // -1 = use default
367 
368  int alternateSmoothingDirections; // 0=no, 1=yes
369 
371 
373 
375 
376  bool allowInterpolationFromGhostPoints; // for coarse level interpolation points that are outside a grid
377  bool allowExtrapolationOfInterpolationPoints; // for coarse level interpolation points that are outside a grid
378 
379  bool saveGridCheckFile; // save the coarse grid check file
380 
381  int chooseGoodParametersOption; // choose good parameters has different levels of robustness, 0=OFF
382 
383  friend class Ogmg;
384 };
385 
386 #endif