ASPECT
introspection.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_introspection_h
23 #define _aspect_introspection_h
24 
25 #include <deal.II/base/index_set.h>
26 #include <deal.II/base/quadrature.h>
27 #include <deal.II/fe/component_mask.h>
28 #include <deal.II/fe/fe_values_extractors.h>
29 #include <deal.II/fe/fe.h>
30 
32 #include <aspect/parameters.h>
33 
34 
35 namespace aspect
36 {
37  using namespace dealii;
38 
43  template <int dim>
44  std::vector<VariableDeclaration<dim>>
45  construct_default_variables (const Parameters<dim> &parameters);
46 
47 
54  {
58  enum Type
59  {
67  generic,
68  unspecified
69  } type;
70 
75  static
76  Type
77  parse_type(const std::string &input)
78  {
79  if (input == "chemical composition")
81  else if (input == "stress")
83  else if (input == "strain")
85  else if (input == "grain size")
87  else if (input == "porosity")
89  else if (input == "density")
91  else if (input == "entropy")
93  else if (input == "generic")
95  else if (input == "unspecified")
97  else
98  AssertThrow(false, ExcNotImplemented());
99 
101  }
102  };
103 
121  template <int dim>
123  {
124  public:
128  Introspection (const std::vector<VariableDeclaration<dim>> &variables,
129  const Parameters<dim> &parameters);
130 
134  ~Introspection ();
135 
136 
137 
150  const unsigned int n_components;
151 
155  const unsigned int n_compositional_fields;
156 
162 
168 
174  {
175  unsigned int velocities[dim];
176  unsigned int pressure;
177  unsigned int temperature;
178  std::vector<unsigned int> compositional_fields;
179  };
185 
191  const unsigned int n_blocks;
192 
198  {
199  unsigned int velocities;
200  unsigned int pressure;
201  unsigned int temperature;
202  std::vector<unsigned int> compositional_fields;
203  };
209 
214  struct Extractors
215  {
216  Extractors (const ComponentIndices &component_indices);
217 
219  const FEValuesExtractors::Scalar pressure;
220  const FEValuesExtractors::Scalar temperature;
221  const std::vector<FEValuesExtractors::Scalar> compositional_fields;
222  };
228 
240  {
241  unsigned int velocities;
242  unsigned int pressure;
243  unsigned int temperature;
244  unsigned int compositional_fields;
245  };
251 
260  {
261  unsigned int velocities;
262  unsigned int temperature;
263  unsigned int compositional_fields;
264  };
270 
290  struct Quadratures
291  {
292  Quadrature<dim> velocities;
293  Quadrature<dim> pressure;
294  Quadrature<dim> temperature;
295  Quadrature<dim> compositional_fields;
296  Quadrature<dim> system;
297  };
303 
313  {
314  Quadrature<dim-1> velocities;
315  Quadrature<dim-1> pressure;
316  Quadrature<dim-1> temperature;
317  Quadrature<dim-1> compositional_fields;
318  Quadrature<dim-1> system;
319  };
325 
332  {
334 
335  ComponentMask velocities;
336  ComponentMask pressure;
337  ComponentMask temperature;
338  std::vector<ComponentMask> compositional_fields;
339  };
346 
360  std::vector<types::global_dof_index> system_dofs_per_block;
361 
367  struct IndexSets
368  {
376 
384  std::vector<IndexSet> system_partitioning;
385 
394  std::vector<IndexSet> system_relevant_partitioning;
395 
401  std::vector<IndexSet> stokes_partitioning;
402 
410 
416 
422  };
429 
435 
442 
455  unsigned int
456  compositional_index_for_name (const std::string &name) const;
457 
464  std::string
465  name_for_compositional_index (const unsigned int index) const;
466 
470  const std::vector<std::string> &
471  get_composition_names () const;
472 
477  const std::vector<CompositionalFieldDescription> &
478  get_composition_descriptions () const;
479 
484  const std::vector<std::string> &
485  chemical_composition_field_names () const;
486 
491  const std::vector<unsigned int> &
492  chemical_composition_field_indices () const;
493 
498  unsigned int
499  n_chemical_composition_fields () const;
500 
509  bool
510  composition_type_exists (const CompositionalFieldDescription::Type &type) const;
511 
512 
522  unsigned int
523  find_composition_type (const CompositionalFieldDescription::Type &type) const;
524 
525 
534  bool
535  compositional_name_exists (const std::string &name) const;
536 
541  const std::vector<unsigned int>
542  get_indices_for_fields_of_type (const CompositionalFieldDescription::Type &type) const;
543 
548  const std::vector<std::string>
549  get_names_for_fields_of_type (const CompositionalFieldDescription::Type &type) const;
550 
555  unsigned int
556  get_number_of_fields_of_type (const CompositionalFieldDescription::Type &type) const;
557 
565  bool
566  is_stokes_component (const unsigned int component_index) const;
567 
568  private:
573  std::vector<std::string> composition_names;
574 
580  std::vector<CompositionalFieldDescription> composition_descriptions;
581 
587  std::vector<std::string> chemical_composition_names;
588 
594  std::vector<unsigned int> chemical_composition_indices;
595 
601  const unsigned int n_chemical_compositions;
602 
603  };
604 }
605 
606 
607 #endif
std::vector< IndexSet > system_relevant_partitioning
std::vector< unsigned int > chemical_composition_indices
const FEValuesExtractors::Vector velocities
std::vector< VariableDeclaration< dim > > construct_default_variables(const Parameters< dim > &parameters)
Parameters< dim >::AdvectionFieldMethod::Kind temperature_method
std::vector< unsigned int > compositional_fields
std::vector< unsigned int > compositional_fields
const unsigned int n_blocks
const bool use_discontinuous_composition_discretization
::TrilinosWrappers::MPI::Vector Vector
Definition: global.h:233
const BaseElements base_elements
const PolynomialDegree polynomial_degree
const FEValuesExtractors::Scalar temperature
const unsigned int n_components
const Extractors extractors
std::vector< IndexSet > system_partitioning
std::vector< typename Parameters< dim >::AdvectionFieldMethod::Kind > compositional_field_methods
std::vector< CompositionalFieldDescription > composition_descriptions
Quadrature< dim > compositional_fields
const unsigned int n_compositional_fields
const ComponentIndices component_indices
const bool use_discontinuous_temperature_discretization
const FaceQuadratures face_quadratures
const unsigned int n_chemical_compositions
const Quadratures quadratures
const BlockIndices block_indices
const std::vector< FEValuesExtractors::Scalar > compositional_fields
const FEValuesExtractors::Scalar pressure
const ComponentMasks component_masks
std::vector< types::global_dof_index > system_dofs_per_block
std::vector< std::string > chemical_composition_names
Definition: compat.h:42
std::vector< std::string > composition_names
std::vector< ComponentMask > compositional_fields
std::vector< IndexSet > stokes_partitioning
static Type parse_type(const std::string &input)
Definition: introspection.h:77