ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2022 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  using namespace dealii;
38 
39  template <int dim> class Simulator;
40  template <int dim> class SimulatorAccess;
41 
42 
49  namespace MeshRefinement
50  {
51 
73  template <int dim>
74  class Interface
75  {
76  public:
81  virtual ~Interface () = default;
82 
88  virtual void initialize ();
89 
99  virtual
100  void
101  update ();
102 
112  virtual
113  void
114  execute (Vector<float> &error_indicators) const;
115 
126  virtual
127  void
128  tag_additional_cells () const;
129 
143  static
144  void
145  declare_parameters (ParameterHandler &prm);
146 
153  virtual
154  void
155  parse_parameters (ParameterHandler &prm);
156  };
157 
158 
159 
160 
161 
162 
169  template <int dim>
170  class Manager : public ::aspect::SimulatorAccess<dim>
171  {
172  public:
177  ~Manager () override;
178 
185  virtual
186  void
187  update ();
188 
195  virtual
196  void
197  execute (Vector<float> &error_indicators) const;
198 
204  virtual
205  void
206  tag_additional_cells () const;
207 
212  static
213  void
214  declare_parameters (ParameterHandler &prm);
215 
221  void
222  parse_parameters (ParameterHandler &prm);
223 
233  template <typename MeshRefinementType,
234  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshRefinementType>::value>>
235  bool
236  has_matching_mesh_refinement_strategy () const;
237 
249  template <typename MeshRefinementType,
250  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshRefinementType>::value>>
251  const MeshRefinementType &
252  get_matching_mesh_refinement_strategy () const;
253 
271  static
272  void
273  register_mesh_refinement_criterion (const std::string &name,
274  const std::string &description,
275  void (*declare_parameters_function) (ParameterHandler &),
276  std::unique_ptr<Interface<dim>> (*factory_function) ());
277 
287  static
288  void
289  write_plugin_graph (std::ostream &output_stream);
290 
294  DeclException1 (ExcMeshRefinementNameNotFound,
295  std::string,
296  << "Could not find entry <"
297  << arg1
298  << "> among the names of registered mesh refinement objects.");
299  private:
305  { plus, max };
306 
311 
317 
322  std::vector<double> scaling_factors;
323 
328  std::list<std::unique_ptr<Interface<dim>>> mesh_refinement_objects;
329  };
330 
331 
332 
333  template <int dim>
334  template <typename MeshRefinementType, typename>
335  inline
336  bool
338  {
339  for (typename std::list<std::unique_ptr<Interface<dim>>>::const_iterator
340  p = mesh_refinement_objects.begin();
341  p != mesh_refinement_objects.end(); ++p)
342  if (Plugins::plugin_type_matches<MeshRefinementType>(*(*p)))
343  return true;
344 
345  return false;
346  }
347 
348 
349 
350  template <int dim>
351  template <typename MeshRefinementType, typename>
352  inline
353  const MeshRefinementType &
355  {
356  AssertThrow(has_matching_mesh_refinement_strategy<MeshRefinementType> (),
357  ExcMessage("You asked MeshRefinement::Manager::get_matching_mesh_refinement_strategy() for a "
358  "mesh refinement strategy of type <" + boost::core::demangle(typeid(MeshRefinementType).name()) + "> "
359  "that could not be found in the current model. Activate this "
360  "mesh refinement strategy in the input file."));
361 
362  for (typename std::list<std::unique_ptr<Interface<dim>>>::const_iterator
363  p = mesh_refinement_objects.begin();
364  p != mesh_refinement_objects.end(); ++p)
365  if (Plugins::plugin_type_matches<MeshRefinementType>(*(*p)))
366  return Plugins::get_plugin_as_type<MeshRefinementType>(*(*p));
367 
368  // We will never get here, because we had the Assert above. Just to avoid warnings.
369  typename std::list<std::unique_ptr<Interface<dim>>>::const_iterator mesh_refinement_strategy;
370  return Plugins::get_plugin_as_type<MeshRefinementType>(*(*mesh_refinement_strategy));
371  }
372 
373 
374 
382 #define ASPECT_REGISTER_MESH_REFINEMENT_CRITERION(classname,name,description) \
383  template class classname<2>; \
384  template class classname<3>; \
385  namespace ASPECT_REGISTER_MESH_REFINEMENT_CRITERION_ ## classname \
386  { \
387  aspect::internal::Plugins::RegisterHelper<aspect::MeshRefinement::Interface<2>,classname<2>> \
388  dummy_ ## classname ## _2d (&aspect::MeshRefinement::Manager<2>::register_mesh_refinement_criterion, \
389  name, description); \
390  aspect::internal::Plugins::RegisterHelper<aspect::MeshRefinement::Interface<3>,classname<3>> \
391  dummy_ ## classname ## _3d (&aspect::MeshRefinement::Manager<3>::register_mesh_refinement_criterion, \
392  name, description); \
393  }
394  }
395 }
396 
397 
398 #endif
MergeOperation merge_operation
Definition: interface.h:310
void write_plugin_graph(std::ostream &output_stream)
std::vector< double > scaling_factors
Definition: interface.h:322
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:42
std::list< std::unique_ptr< Interface< dim > > > mesh_refinement_objects
Definition: interface.h:328
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.")