ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2024 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 
22 #ifndef _aspect_mesh_refinement_interface_h
23 #define _aspect_mesh_refinement_interface_h
24 
25 #include <aspect/global.h>
26 #include <aspect/plugins.h>
28 
29 #include <memory>
30 #include <deal.II/base/table_handler.h>
31 #include <deal.II/base/parameter_handler.h>
32 #include <deal.II/distributed/tria.h>
33 
34 
35 namespace aspect
36 {
37  template <int dim> class Simulator;
38  template <int dim> class SimulatorAccess;
39 
40 
47  namespace MeshRefinement
48  {
49 
71  template <int dim>
73  {
74  public:
84  virtual
85  void
86  execute (Vector<float> &error_indicators) const;
87 
98  virtual
99  void
100  tag_additional_cells () const;
101  };
102 
103 
104 
105 
106 
107 
114  template <int dim>
115  class Manager : public Plugins::ManagerBase<Interface<dim>>, public SimulatorAccess<dim>
116  {
117  public:
124  virtual
125  void
126  execute (Vector<float> &error_indicators) const;
127 
133  virtual
134  void
135  tag_additional_cells () const;
136 
141  static
142  void
143  declare_parameters (ParameterHandler &prm);
144 
150  void
151  parse_parameters (ParameterHandler &prm) override;
152 
167  template <typename MeshRefinementType,
168  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshRefinementType>::value>>
169  DEAL_II_DEPRECATED
170  bool
171  has_matching_mesh_refinement_strategy () const;
172 
189  template <typename MeshRefinementType,
190  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshRefinementType>::value>>
191  DEAL_II_DEPRECATED
192  const MeshRefinementType &
193  get_matching_mesh_refinement_strategy () const;
194 
212  static
213  void
214  register_mesh_refinement_criterion (const std::string &name,
215  const std::string &description,
216  void (*declare_parameters_function) (ParameterHandler &),
217  std::unique_ptr<Interface<dim>> (*factory_function) ());
218 
228  static
229  void
230  write_plugin_graph (std::ostream &output_stream);
231 
235  DeclException1 (ExcMeshRefinementNameNotFound,
236  std::string,
237  << "Could not find entry <"
238  << arg1
239  << "> among the names of registered mesh refinement objects.");
240  private:
246  { plus, max };
247 
252 
258 
263  std::vector<double> scaling_factors;
264  };
265 
266 
267 
268  template <int dim>
269  template <typename MeshRefinementType, typename>
270  inline
271  bool
273  {
274  return this->template has_matching_active_plugin<MeshRefinementType>();
275  }
276 
277 
278 
279  template <int dim>
280  template <typename MeshRefinementType, typename>
281  inline
282  const MeshRefinementType &
284  {
285  return this->template get_matching_active_plugin<MeshRefinementType>();
286  }
287 
288 
289 
297 #define ASPECT_REGISTER_MESH_REFINEMENT_CRITERION(classname,name,description) \
298  template class classname<2>; \
299  template class classname<3>; \
300  namespace ASPECT_REGISTER_MESH_REFINEMENT_CRITERION_ ## classname \
301  { \
302  aspect::internal::Plugins::RegisterHelper<aspect::MeshRefinement::Interface<2>,classname<2>> \
303  dummy_ ## classname ## _2d (&aspect::MeshRefinement::Manager<2>::register_mesh_refinement_criterion, \
304  name, description); \
305  aspect::internal::Plugins::RegisterHelper<aspect::MeshRefinement::Interface<3>,classname<3>> \
306  dummy_ ## classname ## _3d (&aspect::MeshRefinement::Manager<3>::register_mesh_refinement_criterion, \
307  name, description); \
308  }
309  }
310 }
311 
312 
313 #endif
virtual void parse_parameters(ParameterHandler &prm)
MergeOperation merge_operation
Definition: interface.h:251
virtual void tag_additional_cells() const
virtual void execute(Vector< float > &error_indicators) const
DEAL_II_DEPRECATED bool has_matching_mesh_refinement_strategy() const
void write_plugin_graph(std::ostream &output_stream)
DEAL_II_DEPRECATED const MeshRefinementType & get_matching_mesh_refinement_strategy() const
std::vector< double > scaling_factors
Definition: interface.h:263
Definition: compat.h:59
static void declare_parameters(ParameterHandler &prm)
DeclException1(ProbabilityFunctionNegative, Point< dim >,<< "Your probability density function in the particle generator " "returned a negative probability density for the following position: "<< arg1<< ". Please check your function expression.")