ASPECT
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
solution_evaluator.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2024 - 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 #ifndef _aspect_solution_evaluator_h
22 #define _aspect_solution_evaluator_h
23 
24 #include <aspect/global.h>
25 
26 
28 
29 #include <deal.II/base/array_view.h>
30 #include <deal.II/matrix_free/fe_point_evaluation.h>
31 #include <deal.II/non_matching/mapping_info.h>
32 #include <deal.II/dofs/dof_handler.h>
33 
34 namespace aspect
35 {
36  namespace internal
37  {
42  template <int dim>
44  {
45  public:
49  DynamicFEPointEvaluation(const unsigned int first_component, const unsigned int n_components)
50  : first_component (first_component),
51  n_components (n_components)
52  {}
53 
57  virtual ~DynamicFEPointEvaluation() = default;
58 
65  virtual void evaluate(const ArrayView<double> &solution_values,
66  const EvaluationFlags::EvaluationFlags flags) = 0;
67 
71  virtual
73  get_gradient(const unsigned int evaluation_point) const = 0;
74 
79  virtual
80  void
81  get_gradient(const unsigned int evaluation_point,
82  const ArrayView<Tensor<1,dim>> &gradients) const = 0;
83 
87  virtual
89  get_value(const unsigned int evaluation_point) const = 0;
90 
95  virtual
96  void
97  get_value(const unsigned int evaluation_point,
98  const ArrayView<double> &solution) const = 0;
99 
103  virtual
104  unsigned int
105  get_first_component() const final
106  {
107  return first_component;
108  }
109 
113  virtual
114  unsigned int
115  get_n_components() const final
116  {
117  return n_components;
118  }
119 
120  private:
124  unsigned int first_component;
125 
129  unsigned int n_components;
130  };
131  }
132 
138  template <int dim>
140  {
141  public:
147  SolutionEvaluator(const SimulatorAccess<dim> &simulator,
148  const UpdateFlags update_flags);
149 
154  void
155  reinit(const typename DoFHandler<dim>::active_cell_iterator &cell,
156  const ArrayView<Point<dim>> &positions);
157 
167  void
168  evaluate(const ArrayView<double> &solution_values,
169  const std::vector<EvaluationFlags::EvaluationFlags> &evaluation_flags);
170 
183  void get_solution(const unsigned int evaluation_point,
184  const ArrayView<double> &solution,
185  const std::vector<EvaluationFlags::EvaluationFlags> &evaluation_flags) const;
186 
199  void get_gradients(const unsigned int evaluation_point,
200  const ArrayView<Tensor<1, dim>> &gradients,
201  const std::vector<EvaluationFlags::EvaluationFlags> &evaluation_flags) const;
202 
207  FEPointEvaluation<dim, dim> &
208  get_velocity_or_fluid_velocity_evaluator(const bool use_fluid_velocity);
209 
213  NonMatching::MappingInfo<dim> &
214  get_mapping_info();
215 
219  unsigned int
220  n_components() const;
221 
222  private:
226  NonMatching::MappingInfo<dim> mapping_info;
227 
233  FEPointEvaluation<dim, dim> velocity;
234  std::unique_ptr<FEPointEvaluation<1, dim>> pressure;
235  FEPointEvaluation<1, dim> temperature;
236 
241  std::vector<std::unique_ptr<internal::DynamicFEPointEvaluation<dim>>> compositions;
242 
249  std::unique_ptr<FEPointEvaluation<dim, dim>> fluid_velocity;
250  std::unique_ptr<FEPointEvaluation<1, dim>> compaction_pressure;
251  std::unique_ptr<FEPointEvaluation<1, dim>> fluid_pressure;
252 
259  std::array<unsigned int, 3> melt_component_indices;
260 
266  };
267 
271  template <int dim>
272  std::unique_ptr<SolutionEvaluator<dim>>
273  construct_solution_evaluator(const SimulatorAccess<dim> &simulator_access,
274  const UpdateFlags update_flags);
275 }
276 
277 #endif
DynamicFEPointEvaluation(const unsigned int first_component, const unsigned int n_components)
virtual void evaluate(const ArrayView< double > &solution_values, const EvaluationFlags::EvaluationFlags flags)=0
virtual unsigned int get_n_components() const final
std::unique_ptr< FEPointEvaluation< 1, dim > > compaction_pressure
FEPointEvaluation< dim, dim > velocity
FEPointEvaluation< 1, dim > temperature
std::array< unsigned int, 3 > melt_component_indices
std::unique_ptr< FEPointEvaluation< dim, dim > > fluid_velocity
Definition: compat.h:59
NonMatching::MappingInfo< dim > mapping_info
virtual small_vector< Tensor< 1, dim > > get_gradient(const unsigned int evaluation_point) const =0
std::unique_ptr< FEPointEvaluation< 1, dim > > pressure
virtual unsigned int get_first_component() const final
boost::container::small_vector< T, N > small_vector
Definition: global.h:251
virtual small_vector< double > get_value(const unsigned int evaluation_point) const =0
std::unique_ptr< FEPointEvaluation< 1, dim > > fluid_pressure
const SimulatorAccess< dim > & simulator_access
std::vector< std::unique_ptr< internal::DynamicFEPointEvaluation< dim > > > compositions
std::unique_ptr< SolutionEvaluator< dim > > construct_solution_evaluator(const SimulatorAccess< dim > &simulator_access, const UpdateFlags update_flags)