deal2lkit: A ToolKit library for Deal.II
parsed_grid_refinement.h
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2015 by the deal2lkit authors
4 //
5 // This file is part of the deal2lkit library.
6 //
7 // The deal2lkit library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal2lkit distribution.
13 //
14 //-----------------------------------------------------------
15 
16 #ifndef _d2k_parsed_grid_refinement_h
17 #define _d2k_parsed_grid_refinement_h
18 
19 #include <deal.II/grid/tria.h>
21 
22 #ifdef DEAL_II_WITH_MPI
23 #ifdef DEAL_II_WITH_P4EST
26 #endif
27 #endif
28 
29 
30 #include <deal2lkit/config.h>
32 
33 using namespace dealii;
34 
35 
36 D2K_NAMESPACE_OPEN
41 {
42 public:
46  ParsedGridRefinement(const std::string &name="",
47  const std::string &strategy="fraction",
48  const double &top_parameter=.3,
49  const double &bottom_parameter=.1,
50  const unsigned int &max_cells=0,
51  const unsigned int &order=2);
52 
56  virtual void declare_parameters(ParameterHandler &prm);
57 
58 
68  template<int dim, class Vector , int spacedim>
69  void mark_cells(const Vector &criteria,
70  Triangulation< dim, spacedim > &tria) const;
71 
72 #ifdef DEAL_II_WITH_MPI
73 #ifdef DEAL_II_WITH_P4EST
74 
85  template<int dim, class Vector , int spacedim>
86  void mark_cells(const Vector &criteria,
88 #endif
89 #endif
90 
91 private:
95  std::string strategy;
96  double top_parameter;
98  unsigned int max_cells;
99  unsigned int order;
100 };
101 
102 
103 
104 #ifdef DEAL_II_WITH_MPI
105 #ifdef DEAL_II_WITH_P4EST
106 template<int dim, class Vector , int spacedim>
107 void ParsedGridRefinement::mark_cells(const Vector &criteria,
109 {
110  if (strategy == "number")
112  criteria,
113  top_parameter, bottom_parameter,
114  max_cells ? max_cells :
115  std::numeric_limits< unsigned int >::max());
116  else if (strategy == "fraction")
118  criteria,
119  top_parameter, bottom_parameter);
120  else
121  Assert(false, ExcInternalError());
122 
123 }
124 #endif
125 #endif
126 
127 template<int dim, class Vector , int spacedim>
129  Triangulation< dim, spacedim > &tria) const
130 {
131  if (strategy == "number")
133  criteria,
134  top_parameter, bottom_parameter,
135  max_cells ? max_cells :
136  std::numeric_limits< unsigned int >::max());
137  else if (strategy == "fraction")
139  criteria,
140  top_parameter, bottom_parameter,
141  max_cells ? max_cells :
142  std::numeric_limits< unsigned int >::max());
143  // This one does not seem to work properly
144  // else if(strategy == "optimize")
145  // GridRefinement::refine_and_coarsen_optimize (tria, order);
146  else
147  Assert(false, ExcInternalError());
148 }
149 
150 D2K_NAMESPACE_CLOSE
151 
152 
153 #endif
154 
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const VectorType &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
A parameter acceptor base class.
A wrapper for refinement strategies.
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
#define Assert(cond, exc)
std::string strategy
Default expression of this function.
void refine_and_coarsen_fixed_fraction(parallel::distributed::Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error)
void refine_and_coarsen_fixed_number(parallel::distributed::Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
static::ExceptionBase & ExcInternalError()
void mark_cells(const Vector &criteria, Triangulation< dim, spacedim > &tria) const
Mark cells a the triangulation for refinement or coarsening, according to the given strategy applied ...