ASPECT
lateral_averaging.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2023 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_lateral_averaging_h
23 #define _aspect_lateral_averaging_h
24 
26 
27 #include <deal.II/fe/fe_values.h>
28 
29 namespace aspect
30 {
31  namespace internal
32  {
41  template <int dim>
43  {
44  public:
49  virtual
50  bool
52 
57  virtual
58  void
59  create_additional_material_model_outputs (const unsigned int n_points,
61 
66  virtual
67  void
68  setup(const unsigned int q_points);
69 
77  virtual
78  void
81  const FEValues<dim> &fe_values,
82  const LinearAlgebra::BlockVector &solution,
83  std::vector<double> &output) = 0;
84 
88  virtual
89  ~FunctorBase();
90  };
91  }
92 
107  template <int dim>
108  class LateralAveraging : public SimulatorAccess<dim>
109  {
110  public:
127  std::vector<std::vector<double>>
128  compute_lateral_averages(const unsigned int n_slices,
129  const std::vector<std::string> &property_names) const;
130 
153  std::vector<std::vector<double>>
154  compute_lateral_averages(const std::vector<double> &depth_bounds,
155  const std::vector<std::string> &property_names) const;
156 
187  std::vector<std::vector<double>>
188  compute_lateral_averages(const std::vector<double> &depth_bounds,
189  std::vector<std::unique_ptr<internal::FunctorBase<dim>>> &functors) const;
190 
200  void
201  get_temperature_averages(std::vector<double> &values) const;
202 
215  void
216  get_composition_averages(const unsigned int composition_index,
217  std::vector<double> &values) const;
218 
226  void
227  get_viscosity_averages(std::vector<double> &values) const;
228 
236  void
237  get_log_viscosity_averages(std::vector<double> &values) const;
238 
246  void
247  get_velocity_magnitude_averages(std::vector<double> &values) const;
248 
256  void
257  get_sinking_velocity_averages(std::vector<double> &values) const;
258 
266  void
267  get_rising_velocity_averages(std::vector<double> &values) const;
268 
276  void
277  get_Vs_averages(std::vector<double> &values) const;
278 
286  void
287  get_Vp_averages(std::vector<double> &values) const;
288 
297  void
298  get_vertical_heat_flux_averages(std::vector<double> &values) const;
299 
309  void
310  get_vertical_mass_flux_averages(std::vector<double> &values) const;
311  };
312 }
313 
314 
315 #endif
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:268
virtual void create_additional_material_model_outputs(const unsigned int n_points, MaterialModel::MaterialModelOutputs< dim > &outputs) const
virtual void operator()(const MaterialModel::MaterialModelInputs< dim > &in, const MaterialModel::MaterialModelOutputs< dim > &out, const FEValues< dim > &fe_values, const LinearAlgebra::BlockVector &solution, std::vector< double > &output)=0
virtual bool need_material_properties() const
Definition: compat.h:59
virtual void setup(const unsigned int q_points)