1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_AMG_PARAMETERS_HH
4 #define DUNE_AMG_PARAMETERS_HH
5 
6 #include <cstddef>
7 
8 namespace Dune
9 {
10   namespace Amg
11   {
12     /**
13      * @addtogroup ISTL_PAAMG
14      *
15      * @{
16      */
17     /** @file
18      * @author Markus Blatt
19      * @brief Parameter classes for customizing AMG
20      *
21      * All parameters of the AMG can be set by using the class Parameter, which
22      * can be provided to CoarsenCriterion via its constructor.
23      */
24 
25     /**
26      * @brief Parameters needed to check whether a node depends on another.
27      */
28     class DependencyParameters
29     {
30     public:
31       /** @brief Constructor */
DependencyParameters()32       DependencyParameters()
33         : alpha_(1.0/3.0), beta_(1.0E-5)
34       {}
35 
36       /**
37        * @brief Set threshold for marking nodes as isolated.
38        * The default value is 1.0E-5.
39        */
setBeta(double b)40       void setBeta(double b)
41       {
42         beta_ = b;
43       }
44 
45       /**
46        * @brief Get the threshold for marking nodes as isolated.
47        * The default value is 1.0E-5.
48        * @return beta
49        */
beta() const50       double beta() const
51       {
52         return beta_;
53       }
54 
55       /**
56        * @brief Set the scaling value for marking connections as strong.
57        * Default value is 1/3
58        */
setAlpha(double a)59       void setAlpha(double a)
60       {
61         alpha_ = a;
62       }
63 
64       /**
65        * @brief Get the scaling value for marking connections as strong.
66        * Default value is 1/3
67        */
alpha() const68       double alpha() const
69       {
70         return alpha_;
71       }
72 
73     private:
74       double alpha_, beta_;
75     };
76 
77     /**
78      * @brief Parameters needed for the aggregation process
79      */
80     class AggregationParameters :
81       public DependencyParameters
82     {
83     public:
84       /**
85        * @brief Constructor.
86        *
87        * The parameters will be initialized with default values suitable
88        * for 2D isotropic problems.
89        *
90        * If that does not fit your needs either use setDefaultValuesIsotropic
91        * setDefaultValuesAnisotropic or setup the values by hand
92        */
AggregationParameters()93       AggregationParameters()
94         : maxDistance_(2), minAggregateSize_(4), maxAggregateSize_(6),
95           connectivity_(15), skipiso_(false)
96       {}
97 
98       /**
99        * @brief Sets reasonable default values for an isotropic problem.
100        *
101        * Reasonable means that we should end up with cube aggregates of
102        * diameter 2.
103        *
104        * @param dim The dimension of the problem.
105        * @param diameter The preferred diameter for the aggregation.
106        */
setDefaultValuesIsotropic(std::size_t dim,std::size_t diameter=2)107       void setDefaultValuesIsotropic(std::size_t dim, std::size_t diameter=2)
108       {
109         maxDistance_=diameter-1;
110         std::size_t csize=1;
111 
112         for(; dim>0; dim--) {
113           csize*=diameter;
114           maxDistance_+=diameter-1;
115         }
116         minAggregateSize_=csize;
117         maxAggregateSize_=static_cast<std::size_t>(csize*1.5);
118       }
119 
120       /**
121        * @brief Sets reasonable default values for an anisotropic problem.
122        *
123        * Reasonable means that we should end up with cube aggregates with
124        * sides of diameter 2 and sides in one dimension that are longer
125        * (e.g. for 3D: 2x2x3).
126        *
127        * @param dim The dimension of the problem.
128        * @param diameter The preferred diameter for the aggregation.
129        */
setDefaultValuesAnisotropic(std::size_t dim,std::size_t diameter=2)130       void setDefaultValuesAnisotropic(std::size_t dim,std::size_t diameter=2)
131       {
132         setDefaultValuesIsotropic(dim, diameter);
133         maxDistance_+=dim-1;
134       }
135       /**
136        * @brief Get the maximal distance allowed between two nodes in a aggregate.
137        *
138        * The distance between two nodes in a aggregate is the minimal number of edges
139        * it takes to travel from one node to the other without leaving the aggregate.
140        * @return The maximum distance allowed.
141        */
maxDistance() const142       std::size_t maxDistance() const { return maxDistance_;}
143 
144       /**
145        * @brief Set the maximal distance allowed between two nodes in a aggregate.
146        *
147        * The distance between two nodes in a aggregate is the minimal number of edges
148        * it takes to travel from one node to the other without leaving the aggregate.
149        * The default value is 2.
150        * @param distance The maximum distance allowed.
151        */
setMaxDistance(std::size_t distance)152       void setMaxDistance(std::size_t distance) { maxDistance_ = distance;}
153 
154       /**
155        * @brief Whether isolated aggregates will not be represented on
156        * the coarse level.
157        * @return True if these aggregates will be skipped.
158        */
skipIsolated() const159       bool skipIsolated() const
160       {
161         return skipiso_;
162       }
163 
164       /**
165        * @brief Set whether isolated aggregates will not be represented on
166        * the coarse level.
167        * @param skip True if these aggregates will be skipped.
168        */
setSkipIsolated(bool skip)169       void setSkipIsolated(bool skip)
170       {
171         skipiso_=skip;
172       }
173 
174       /**
175        * @brief Get the minimum number of nodes a aggregate has to consist of.
176        * @return The minimum number of nodes.
177        */
minAggregateSize() const178       std::size_t minAggregateSize() const { return minAggregateSize_;}
179 
180       /**
181        * @brief Set the minimum number of nodes a aggregate has to consist of.
182        *
183        * the default value is 4.
184        * @return The minimum number of nodes.
185        */
setMinAggregateSize(std::size_t size)186       void setMinAggregateSize(std::size_t size){ minAggregateSize_=size;}
187 
188       /**
189        * @brief Get the maximum number of nodes a aggregate is allowed to have.
190        * @return The maximum number of nodes.
191        */
maxAggregateSize() const192       std::size_t maxAggregateSize() const { return maxAggregateSize_;}
193 
194       /**
195        * @brief Set the maximum number of nodes a aggregate is allowed to have.
196        *
197        * The default values is 6.
198        * @param size The maximum number of nodes.
199        */
setMaxAggregateSize(std::size_t size)200       void setMaxAggregateSize(std::size_t size){ maxAggregateSize_ = size;}
201 
202       /**
203        * @brief Get the maximum number of connections a aggregate is allowed to have.
204        *
205        * This limit exists to achieve sparsity of the coarse matrix. the default value is 15.
206        *
207        * @return The maximum number of connections a aggregate is allowed to have.
208        */
maxConnectivity() const209       std::size_t maxConnectivity() const { return connectivity_;}
210 
211       /**
212        * @brief Set the maximum number of connections a aggregate is allowed to have.
213        *
214        * This limit exists to achieve sparsity of the coarse matrix. the default value is 15.
215        *
216        * @param connectivity The maximum number of connections a aggregate is allowed to have.
217        */
setMaxConnectivity(std::size_t connectivity)218       void setMaxConnectivity(std::size_t connectivity){ connectivity_ = connectivity;}
219 
220     private:
221       std::size_t maxDistance_, minAggregateSize_, maxAggregateSize_, connectivity_;
222       bool skipiso_;
223 
224     };
225 
226 
227     /**
228      * @brief Identifiers for the different accumulation modes.
229      */
230     enum AccumulationMode {
231       /**
232        * @brief No data accumulution.
233        *
234        * The coarse level data will be distributed to all processes.
235        */
236       noAccu = 0,
237       /**
238        * @brief Accumulate data to one process at once
239        *
240        * Once no further coarsening is possible all data will be accumulated to one process
241        */
242       atOnceAccu=1,
243       /**
244        * @brief Successively accumulate to fewer processes.
245        */
246       successiveAccu=2
247     };
248 
249 
250 
251 
252     /**
253      * @brief Parameters for the complete coarsening process.
254      */
255     class CoarseningParameters : public AggregationParameters
256     {
257     public:
258       /**
259        * @brief Set the maximum number of levels allowed in the hierarchy.
260        */
setMaxLevel(int l)261       void setMaxLevel(int l)
262       {
263         maxLevel_ = l;
264       }
265       /**
266        * @brief Get the maximum number of levels allowed in the hierarchy.
267        */
maxLevel() const268       int maxLevel() const
269       {
270         return maxLevel_;
271       }
272 
273       /**
274        * @brief Set the maximum number of unknowns allowed on the coarsest level.
275        */
setCoarsenTarget(int nodes)276       void setCoarsenTarget(int nodes)
277       {
278         coarsenTarget_ = nodes;
279       }
280 
281       /**
282        * @brief Get the maximum number of unknowns allowed on the coarsest level.
283        */
coarsenTarget() const284       int coarsenTarget() const
285       {
286         return coarsenTarget_;
287       }
288 
289       /**
290        * @brief Set the minimum coarsening rate to be achieved in each coarsening.
291        *
292        * The default value is 1.2
293        */
setMinCoarsenRate(double rate)294       void setMinCoarsenRate(double rate)
295       {
296         minCoarsenRate_ = rate;
297       }
298 
299       /**
300        * @brief Get the minimum coarsening rate to be achieved.
301        */
minCoarsenRate() const302       double minCoarsenRate() const
303       {
304         return minCoarsenRate_;
305       }
306 
307       /**
308        * @brief Whether the data should be accumulated on fewer processes on coarser levels.
309        */
accumulate() const310       AccumulationMode accumulate() const
311       {
312         return accumulate_;
313       }
314       /**
315        * @brief Set whether he data should be accumulated on fewer processes on coarser levels.
316        */
setAccumulate(AccumulationMode accu)317       void setAccumulate(AccumulationMode accu)
318       {
319         accumulate_=accu;
320       }
321 
setAccumulate(bool accu)322       void setAccumulate(bool accu){
323         accumulate_=accu ? successiveAccu : noAccu;
324       }
325       /**
326        * @brief Set the damping factor for the prolongation.
327        *
328        * @param d The new damping factor.
329        */
setProlongationDampingFactor(double d)330       void setProlongationDampingFactor(double d)
331       {
332         dampingFactor_ = d;
333       }
334 
335       /**
336        * @brief Get the damping factor for the prolongation.
337        *
338        * @return d The damping factor.
339        */
getProlongationDampingFactor() const340       double getProlongationDampingFactor() const
341       {
342         return dampingFactor_;
343       }
344       /**
345        * @brief Constructor
346        * @param maxLevel The maximum number of levels allowed in the matrix hierarchy (default: 100).
347        * @param coarsenTarget If the number of nodes in the matrix is below this threshold the
348        * coarsening will stop (default: 1000).
349        * @param minCoarsenRate If the coarsening rate falls below this threshold the
350        * coarsening will stop (default: 1.2)
351        * @param prolongDamp The damping factor to apply to the prolongated update (default: 1.6)
352        * @param accumulate Whether to accumulate the data onto fewer processors on coarser levels.
353        */
CoarseningParameters(int maxLevel=100,int coarsenTarget=1000,double minCoarsenRate=1.2,double prolongDamp=1.6,AccumulationMode accumulate=successiveAccu)354       CoarseningParameters(int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2,
355                            double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu)
356         : maxLevel_(maxLevel), coarsenTarget_(coarsenTarget), minCoarsenRate_(minCoarsenRate),
357           dampingFactor_(prolongDamp), accumulate_( accumulate)
358       {}
359 
360     private:
361       /**
362        * @brief The maximum number of levels allowed in the hierarchy.
363        */
364       int maxLevel_;
365       /**
366        * @brief The maximum number of unknowns allowed on the coarsest level.
367        */
368       int coarsenTarget_;
369       /**
370        * @brief The minimum coarsening rate to be achieved.
371        */
372       double minCoarsenRate_;
373       /**
374        * @brief The damping factor to apply to the prologated correction.
375        */
376       double dampingFactor_;
377       /**
378        * @brief Whether the data should be agglomerated to fewer processor on
379        * coarser levels.
380        */
381       AccumulationMode accumulate_;
382     };
383 
384     /**
385      * @brief All parameters for AMG.
386      *
387      * Instances of this class can be provided to CoarsenCriterion via its
388      * constructor.
389      */
390     class Parameters : public CoarseningParameters
391     {
392     public:
393       /**
394        * @brief Set the debugging level.
395        *
396        * @param level If 0 no debugging output will be generated.
397        * @warning In parallel the level has to be consistent over all procceses.
398        */
setDebugLevel(int level)399       void setDebugLevel(int level)
400       {
401         debugLevel_ = level;
402       }
403 
404       /**
405        * @brief Get the debugging Level.
406        *
407        * @return 0 if no debugging output will be generated.
408        */
debugLevel() const409       int debugLevel() const
410       {
411         return debugLevel_;
412       }
413 
414       /**
415        * @brief Set the number of presmoothing steps to apply
416        * @param steps The number of steps:
417        */
setNoPreSmoothSteps(std::size_t steps)418       void setNoPreSmoothSteps(std::size_t steps)
419       {
420         preSmoothSteps_=steps;
421       }
422       /**
423        * @brief Get the number of presmoothing steps to apply
424        * @return The number of steps:
425        */
getNoPreSmoothSteps() const426       std::size_t getNoPreSmoothSteps() const
427       {
428         return preSmoothSteps_;
429       }
430 
431       /**
432        * @brief Set the number of postsmoothing steps to apply
433        * @param steps The number of steps:
434        */
setNoPostSmoothSteps(std::size_t steps)435       void setNoPostSmoothSteps(std::size_t steps)
436       {
437         postSmoothSteps_=steps;
438       }
439       /**
440        * @brief Get the number of postsmoothing steps to apply
441        * @return The number of steps:
442        */
getNoPostSmoothSteps() const443       std::size_t getNoPostSmoothSteps() const
444       {
445         return postSmoothSteps_;
446       }
447 
448       /**
449        * @brief Set the value of gamma; 1 for V-cycle, 2 for W-cycle
450        */
setGamma(std::size_t gamma)451       void setGamma(std::size_t gamma)
452       {
453         gamma_=gamma;
454       }
455       /**
456        * @brief Get the value of gamma; 1 for V-cycle, 2 for W-cycle
457        */
getGamma() const458       std::size_t getGamma() const
459       {
460         return gamma_;
461       }
462 
463       /**
464        * @brief Set whether to use additive multigrid.
465        * @param additive True if multigrid should be additive.
466        */
setAdditive(bool additive)467       void setAdditive(bool additive)
468       {
469         additive_=additive;
470       }
471 
472       /**
473        * @brief Get whether to use additive multigrid.
474        * @return True if multigrid should be additive.
475        */
getAdditive() const476       bool getAdditive() const
477       {
478         return additive_;
479       }
480 
481       /**
482        * @brief Constructor
483        * @param maxLevel The maximum number of levels allowed in the matrix hierarchy (default: 100).
484        * @param coarsenTarget If the number of nodes in the matrix is below this threshold the
485        * coarsening will stop (default: 1000).
486        * @param minCoarsenRate If the coarsening rate falls below this threshold the
487        * coarsening will stop (default: 1.2)
488        * @param prolongDamp The damping factor to apply to the prolongated update (default: 1.6)
489        * @param accumulate Whether to accumulate the data onto fewer processors on coarser levels.
490        */
Parameters(int maxLevel=100,int coarsenTarget=1000,double minCoarsenRate=1.2,double prolongDamp=1.6,AccumulationMode accumulate=successiveAccu)491       Parameters(int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2,
492                  double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu)
493         : CoarseningParameters(maxLevel, coarsenTarget, minCoarsenRate, prolongDamp, accumulate)
494           , debugLevel_(2), preSmoothSteps_(2), postSmoothSteps_(2), gamma_(1),
495           additive_(false)
496       {}
497     private:
498       int debugLevel_;
499       std::size_t preSmoothSteps_;
500       std::size_t postSmoothSteps_;
501       std::size_t gamma_;
502       bool additive_;
503     };
504 
505   } //namespace AMG
506 
507 } //namespace Dune
508 #endif
509