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  using namespace dealii;
32 
33  namespace internal
34  {
43  template <int dim>
45  {
46  public:
51  virtual
52  bool
53  need_material_properties() const;
54 
59  virtual
60  void
61  create_additional_material_model_outputs (const unsigned int n_points,
63 
68  virtual
69  void
70  setup(const unsigned int q_points);
71 
79  virtual
80  void
81  operator()(const MaterialModel::MaterialModelInputs<dim> &in,
83  const FEValues<dim> &fe_values,
84  const LinearAlgebra::BlockVector &solution,
85  std::vector<double> &output) = 0;
86 
90  virtual
91  ~FunctorBase();
92  };
93  }
94 
109  template <int dim>
110  class LateralAveraging : public SimulatorAccess<dim>
111  {
112  public:
129  std::vector<std::vector<double>>
130  compute_lateral_averages(const unsigned int n_slices,
131  const std::vector<std::string> &property_names) const;
132 
155  std::vector<std::vector<double>>
156  compute_lateral_averages(const std::vector<double> &depth_bounds,
157  const std::vector<std::string> &property_names) const;
158 
189  std::vector<std::vector<double>>
190  compute_lateral_averages(const std::vector<double> &depth_bounds,
191  std::vector<std::unique_ptr<internal::FunctorBase<dim>>> &functors) const;
192 
202  void
203  get_temperature_averages(std::vector<double> &values) const;
204 
217  void
218  get_composition_averages(const unsigned int composition_index,
219  std::vector<double> &values) const;
220 
228  void
229  get_viscosity_averages(std::vector<double> &values) const;
230 
238  void
239  get_log_viscosity_averages(std::vector<double> &values) const;
240 
248  void
249  get_velocity_magnitude_averages(std::vector<double> &values) const;
250 
258  void
259  get_sinking_velocity_averages(std::vector<double> &values) const;
260 
268  void
269  get_rising_velocity_averages(std::vector<double> &values) const;
270 
278  void
279  get_Vs_averages(std::vector<double> &values) const;
280 
288  void
289  get_Vp_averages(std::vector<double> &values) const;
290 
299  void
300  get_vertical_heat_flux_averages(std::vector<double> &values) const;
301 
311  void
312  get_vertical_mass_flux_averages(std::vector<double> &values) const;
313  };
314 }
315 
316 
317 #endif
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:261
Definition: compat.h:59
Definition: compat.h:42