ASPECT
introspection.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 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 
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 #include <map>
35 
36 namespace aspect
37 {
38  using namespace dealii;
39 
44  template <int dim>
45  std::vector<VariableDeclaration<dim>>
46  construct_default_variables (const Parameters<dim> &parameters);
47 
48 
55  {
59  enum Type
60  {
61  chemical_composition = 0,
62  stress = 1,
63  strain = 2,
64  grain_size = 3,
65  porosity = 4,
66  density = 5,
67  entropy = 6,
68  generic = 7,
69  unspecified = 8
70  } type;
71 
75  constexpr static unsigned int n_types = 9;
76 
81  static
82  Type
83  parse_type(const std::string &input)
84  {
85  if (input == "chemical composition")
87  else if (input == "stress")
89  else if (input == "strain")
91  else if (input == "grain size")
93  else if (input == "porosity")
95  else if (input == "density")
97  else if (input == "entropy")
99  else if (input == "generic")
101  else if (input == "unspecified")
103  else
104  AssertThrow(false, ExcNotImplemented());
105 
107  }
108  };
109 
127  template <int dim>
129  {
130  public:
134  Introspection (const std::vector<VariableDeclaration<dim>> &variables,
135  const Parameters<dim> &parameters);
136 
140  ~Introspection ();
141 
142 
143 
156  const unsigned int n_components;
157 
161  const unsigned int n_compositional_fields;
162 
168 
174 
180  {
181  std::array<unsigned int, dim> velocities;
182  unsigned int pressure;
183  unsigned int temperature;
184  std::vector<unsigned int> compositional_fields;
185  };
191 
197  const unsigned int n_blocks;
198 
204  {
205  unsigned int velocities;
206  unsigned int pressure;
207  unsigned int temperature;
208  std::vector<unsigned int> compositional_fields;
209 
217  std::vector<unsigned int> compositional_field_sparsity_pattern;
218  };
219 
225 
230  struct Extractors
231  {
232  Extractors (const ComponentIndices &component_indices);
233 
235  const FEValuesExtractors::Scalar pressure;
236  const FEValuesExtractors::Scalar temperature;
237  const std::vector<FEValuesExtractors::Scalar> compositional_fields;
238  };
239 
245 
256  {
257  unsigned int velocities;
258  unsigned int pressure;
259  unsigned int temperature;
260  std::vector<unsigned int> compositional_fields;
261  };
262 
268 
277  {
278  unsigned int max_degree;
279  unsigned int velocities;
280  unsigned int temperature;
281  std::vector<unsigned int> compositional_fields;
283  };
284 
290 
310  struct Quadratures
311  {
312  Quadrature<dim> velocities;
313  Quadrature<dim> pressure;
314  Quadrature<dim> temperature;
315  Quadrature<dim> compositional_field_max;
316  std::vector<Quadrature<dim>> compositional_fields;
317  Quadrature<dim> system;
318  };
319 
325 
335  {
336  Quadrature<dim-1> velocities;
337  Quadrature<dim-1> pressure;
338  Quadrature<dim-1> temperature;
339  Quadrature<dim-1> compositional_fields;
340  Quadrature<dim-1> system;
341  };
342 
348 
355  {
357 
361  ComponentMask velocities;
362 
366  ComponentMask pressure;
367 
371  ComponentMask temperature;
372 
379  std::vector<ComponentMask> compositional_fields;
380 
385  ComponentMask compositions;
386  };
387 
394 
408  std::vector<types::global_dof_index> system_dofs_per_block;
409 
415  struct IndexSets
416  {
424 
432  std::vector<IndexSet> system_partitioning;
433 
442  std::vector<IndexSet> system_relevant_partitioning;
443 
449  std::vector<IndexSet> stokes_partitioning;
450 
458 
464 
470  };
471 
478 
484 
491 
508  const std::vector<unsigned int> &
509  get_composition_base_element_indices() const;
510 
518  const std::vector<unsigned int> &
519  get_compositional_field_indices_with_base_element(const unsigned int base_element_index) const;
520 
529  unsigned int
530  compositional_index_for_name (const std::string &name) const;
531 
538  std::string
539  name_for_compositional_index (const unsigned int index) const;
540 
544  const std::vector<std::string> &
545  get_composition_names () const;
546 
551  const std::vector<CompositionalFieldDescription> &
552  get_composition_descriptions () const;
553 
561  const std::vector<std::string> &
562  chemical_composition_field_names () const;
563 
571  const std::vector<unsigned int> &
572  chemical_composition_field_indices () const;
573 
581  unsigned int
582  n_chemical_composition_fields () const;
583 
592  bool
593  composition_type_exists (const CompositionalFieldDescription::Type &type) const;
594 
604  unsigned int
605  find_composition_type (const CompositionalFieldDescription::Type &type) const;
606 
615  bool
616  compositional_name_exists (const std::string &name) const;
617 
622  const std::vector<unsigned int> &
623  get_indices_for_fields_of_type (const CompositionalFieldDescription::Type &type) const;
624 
629  const std::vector<std::string> &
630  get_names_for_fields_of_type (const CompositionalFieldDescription::Type &type) const;
631 
636  unsigned int
637  get_number_of_fields_of_type (const CompositionalFieldDescription::Type &type) const;
638 
646  bool
647  is_stokes_component (const unsigned int component_index) const;
648 
656  bool
657  is_composition_component (const unsigned int component_index) const;
658 
659  private:
664  std::vector<std::string> composition_names;
665 
671  std::vector<CompositionalFieldDescription> composition_descriptions;
672 
678  std::vector<std::vector<std::string>> composition_names_for_type;
679 
685  std::vector<std::vector<unsigned int>> composition_indices_for_type;
686 
691  std::vector<unsigned int> composition_base_element_indices;
692 
698  std::map<unsigned int, std::vector<unsigned int>> compositional_field_indices_with_base_element;
699  };
700 }
701 
702 
703 #endif
std::vector< IndexSet > system_relevant_partitioning
const std::vector< bool > use_discontinuous_composition_discretization
std::map< unsigned int, std::vector< unsigned int > > compositional_field_indices_with_base_element
const FEValuesExtractors::Vector velocities
Quadrature< dim > compositional_field_max
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
std::vector< std::vector< unsigned int > > composition_indices_for_type
std::vector< Quadrature< dim > > compositional_fields
::TrilinosWrappers::MPI::Vector Vector
Definition: global.h:256
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
std::vector< std::vector< std::string > > composition_names_for_type
const unsigned int n_compositional_fields
const ComponentIndices component_indices
std::vector< unsigned int > compositional_fields
Definition: compat.h:59
std::vector< unsigned int > compositional_fields
const bool use_discontinuous_temperature_discretization
std::array< unsigned int, dim > velocities
const FaceQuadratures face_quadratures
const Quadratures quadratures
const BlockIndices block_indices
std::vector< unsigned int > compositional_field_sparsity_pattern
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
Definition: compat.h:42
std::vector< unsigned int > composition_base_element_indices
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:83