ASPECT
interface.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 #ifndef _aspect_material_model_interface_h
22 #define _aspect_material_model_interface_h
23 
24 #include <aspect/global.h>
25 #include <aspect/plugins.h>
27 
28 #include <deal.II/base/point.h>
29 #include <deal.II/base/quadrature.h>
30 #include <deal.II/base/symmetric_tensor.h>
31 #include <deal.II/base/parameter_handler.h>
32 #include <deal.II/dofs/dof_handler.h>
33 #include <deal.II/dofs/dof_accessor.h>
34 #include <deal.II/fe/mapping.h>
35 #include <deal.II/fe/fe_values.h>
36 #include <deal.II/fe/component_mask.h>
37 #include <deal.II/numerics/data_postprocessor.h>
38 #include <deal.II/base/signaling_nan.h>
39 
40 namespace aspect
41 {
42  template <int dim>
43  struct Introspection;
44 
45 
53  namespace MaterialModel
54  {
55  using namespace dealii;
56 
61  namespace NonlinearDependence
62  {
85  {
87 
88  none = 1,
90  pressure = 4,
93 
95  };
96 
97 
102  const Dependence d2)
103  {
104  return Dependence(static_cast<int>(d1) | static_cast<int>(d2));
105  }
106 
108  const Dependence d2)
109  {
110  d1 = (d1 | d2);
111  return d1;
112  }
113 
119  {
125 
131 
137 
143 
149 
155  ModelDependence ();
156  };
157 
166  bool
167  identifies_single_variable(const Dependence dependence);
168  }
169 
174  namespace MaterialProperties
175  {
185  enum Property
186  {
188 
189  none = 1,
191  density = 4,
201 
204  specific_heat |
209  viscosity |
214  };
215 
220  inline Property operator | (const Property d1,
221  const Property d2)
222  {
223  return Property(static_cast<int>(d1) | static_cast<int>(d2));
224  }
225 
227  const Property d2)
228  {
229  return (d1 | d2);
230  }
231  }
232 
233 
234  // Forward declaration:
235  template <int dim>
237 
244  template <int dim>
246  {
247  public:
258  MaterialModelInputs(const unsigned int n_points,
259  const unsigned int n_comp);
260 
276  MaterialModelInputs(const DataPostprocessorInputs::Vector<dim> &input_data,
277  const Introspection<dim> &introspection,
278  const bool compute_strain_rate = true);
279 
280 
300  MaterialModelInputs(const FEValuesBase<dim,dim> &fe_values,
301  const typename DoFHandler<dim>::active_cell_iterator &cell,
302  const Introspection<dim> &introspection,
303  const LinearAlgebra::BlockVector &solution_vector,
304  const bool compute_strain_rate = true);
305 
320 
324  MaterialModelInputs (MaterialModelInputs &&) noexcept = default;
325 
330  MaterialModelInputs &operator= (const MaterialModelInputs &source) = delete;
331 
335  MaterialModelInputs &operator= (MaterialModelInputs &&) = default;
336 
342  void reinit(const FEValuesBase<dim,dim> &fe_values,
343  const typename DoFHandler<dim>::active_cell_iterator &cell,
344  const Introspection<dim> &introspection,
345  const LinearAlgebra::BlockVector &solution_vector,
346  const bool compute_strain_rate = true);
347 
352  unsigned int n_evaluation_points() const;
353 
360  bool requests_property(const MaterialProperties::Property &property) const;
361 
366  std::vector<Point<dim>> position;
367 
371  std::vector<double> temperature;
372 
376  std::vector<double> pressure;
377 
382  std::vector<Tensor<1,dim>> pressure_gradient;
383 
391  std::vector<Tensor<1,dim>> velocity;
392 
398  std::vector<std::vector<double>> composition;
399 
412  std::vector<SymmetricTensor<2,dim>> strain_rate;
413 
429  typename DoFHandler<dim>::active_cell_iterator current_cell;
430 
439 
449  template <class AdditionalInputType>
450  AdditionalInputType *get_additional_input();
451 
456  template <class AdditionalInputType>
457  const AdditionalInputType *get_additional_input() const;
458 
464  std::vector<std::unique_ptr<AdditionalMaterialInputs<dim>>> additional_inputs;
465  };
466 
467 
468  // Forward declaration:
469  template <int dim>
471 
472 
479  template <int dim>
481  {
482  public:
493  MaterialModelOutputs (const unsigned int n_points,
494  const unsigned int n_comp);
495 
510 
514  MaterialModelOutputs (MaterialModelOutputs &&) noexcept = default;
515 
520  MaterialModelOutputs &operator= (const MaterialModelOutputs &source) = delete;
521 
525  MaterialModelOutputs &operator= (MaterialModelOutputs &&) noexcept = default;
526 
531  unsigned int n_evaluation_points() const;
532 
536  std::vector<double> viscosities;
537 
541  std::vector<double> densities;
542 
547  std::vector<double> thermal_expansion_coefficients;
548 
552  std::vector<double> specific_heat;
553 
557  std::vector<double> thermal_conductivities;
558 
563  std::vector<double> compressibilities;
564 
570  std::vector<double> entropy_derivative_pressure;
571 
578  std::vector<double> entropy_derivative_temperature;
579 
616  std::vector<std::vector<double>> reaction_terms;
617 
623  std::vector<std::unique_ptr<AdditionalMaterialOutputs<dim>>> additional_outputs;
624 
632  template <class AdditionalOutputType>
633  AdditionalOutputType *get_additional_output();
634 
638  template <class AdditionalOutputType>
639  const AdditionalOutputType *get_additional_output() const;
640 
645  void move_additional_outputs_from(MaterialModelOutputs<dim> &other);
646  };
647 
648 
649 
664  namespace MaterialAveraging
665  {
705  {
717  };
718 
719 
726  std::string get_averaging_operation_names ();
727 
733  AveragingOperation parse_averaging_operation_name (const std::string &s);
734 
740  template <int dim>
741  void average (const AveragingOperation operation,
742  const typename DoFHandler<dim>::active_cell_iterator &cell,
743  const Quadrature<dim> &quadrature_formula,
744  const Mapping<dim> &mapping,
745  const MaterialProperties::Property &requested_properties,
746  MaterialModelOutputs<dim> &values_out);
747 
753  void average_property (const AveragingOperation operation,
754  const FullMatrix<double> &projection_matrix,
755  const FullMatrix<double> &expansion_matrix,
756  std::vector<double> &values_out);
757  }
758 
759 
778  template <int dim>
780  {
781  public:
786  virtual ~AdditionalMaterialInputs() = default;
787 
792  virtual void
793  fill (const LinearAlgebra::BlockVector &solution,
794  const FEValuesBase<dim> &fe_values,
795  const Introspection<dim> &introspection) = 0;
796  };
797 
798 
825  template <int dim>
827  {
828  public:
833  virtual ~AdditionalMaterialOutputs() = default;
834 
835  virtual void average (const MaterialAveraging::AveragingOperation /*operation*/,
836  const FullMatrix<double> &/*projection_matrix*/,
837  const FullMatrix<double> &/*expansion_matrix*/)
838  {}
839  };
840 
841 
863  template <int dim>
865  {
866  public:
875  NamedAdditionalMaterialOutputs(const std::vector<std::string> &output_names);
876 
887  NamedAdditionalMaterialOutputs(const std::vector<std::string> &output_names,
888  const unsigned int n_points);
889 
893  ~NamedAdditionalMaterialOutputs() override;
894 
899  const std::vector<std::string> &get_names() const;
900 
905  virtual std::vector<double> get_nth_output(const unsigned int idx) const;
906 
908  const FullMatrix<double> &/*projection_matrix*/,
909  const FullMatrix<double> &/*expansion_matrix*/) override
910  {}
911 
912 
917  std::vector<std::vector<double>> output_values;
918 
919  private:
920  const std::vector<std::string> names;
921  };
922 
923 
929  template <int dim>
931  {
932  public:
933  SeismicAdditionalOutputs(const unsigned int n_points);
934 
935  std::vector<double> get_nth_output(const unsigned int idx) const override;
936 
942  std::vector<double> vs;
943 
949  std::vector<double> vp;
950  };
951 
952 
973  template <int dim>
975  {
976  public:
977  ReactionRateOutputs (const unsigned int n_points,
978  const unsigned int n_comp);
979 
980  std::vector<double> get_nth_output(const unsigned int idx) const override;
981 
989  std::vector<std::vector<double>> reaction_rates;
990  };
991 
992 
993 
999  template <int dim>
1001  {
1002  public:
1003  PhaseOutputs(const unsigned int n_points);
1004  };
1005 
1006 
1007 
1026  template <int dim>
1028  {
1029  public:
1030  PrescribedFieldOutputs (const unsigned int n_points,
1031  const unsigned int n_comp);
1032 
1033  std::vector<double> get_nth_output(const unsigned int idx) const override;
1034 
1042  std::vector<std::vector<double>> prescribed_field_outputs;
1043  };
1044 
1058  template <int dim>
1060  {
1061  public:
1062  PrescribedTemperatureOutputs (const unsigned int n_points);
1063 
1064  std::vector<double> get_nth_output(const unsigned int idx) const override;
1065 
1073  std::vector<double> prescribed_temperature_outputs;
1074  };
1075 
1082  template <int dim>
1084  {
1085  public:
1086  AdditionalMaterialOutputsStokesRHS(const unsigned int n_points)
1087  : rhs_u(n_points), rhs_p(n_points), rhs_melt_pc(n_points)
1088  {}
1089 
1091  = default;
1092 
1094  const FullMatrix<double> &/*projection_matrix*/,
1095  const FullMatrix<double> &/*expansion_matrix*/) override
1096  {
1097  // TODO: not implemented
1098  }
1099 
1105  std::vector<Tensor<1,dim>> rhs_u;
1106 
1111  std::vector<double> rhs_p;
1112 
1117  std::vector<double> rhs_melt_pc;
1118  };
1119 
1120 
1121 
1137  template <int dim>
1139  {
1140  public:
1144  explicit PrescribedPlasticDilation (const unsigned int n_points);
1145 
1149  std::vector<double> get_nth_output(const unsigned int idx) const override;
1150 
1155  std::vector<double> dilation;
1156  };
1157 
1158 
1159 
1166  template <int dim>
1168  {
1169  public:
1170  ElasticOutputs(const unsigned int n_points)
1171  : elastic_force(n_points, numbers::signaling_nan<SymmetricTensor<2,dim>>())
1172  , viscoelastic_strain_rate(n_points, numbers::signaling_nan<SymmetricTensor<2,dim>>())
1173  {}
1174 
1175  ~ElasticOutputs() override
1176  = default;
1177 
1179  const FullMatrix<double> &/*projection_matrix*/,
1180  const FullMatrix<double> &/*expansion_matrix*/) override
1181  {
1182  AssertThrow(operation == MaterialAveraging::AveragingOperation::none,ExcNotImplemented());
1183  return;
1184  }
1185 
1191  std::vector<SymmetricTensor<2,dim>> elastic_force;
1192 
1197  std::vector<SymmetricTensor<2,dim>> viscoelastic_strain_rate;
1198  };
1199 
1200 
1201 
1212  template <int dim>
1214  {
1215  public:
1216  EnthalpyOutputs(const unsigned int n_points)
1217  : enthalpies_of_fusion(n_points, numbers::signaling_nan<double>())
1218  {}
1219 
1220  virtual ~EnthalpyOutputs()
1221  = default;
1222 
1224  const FullMatrix<double> &/*projection_matrix*/,
1225  const FullMatrix<double> &/*expansion_matrix*/) override
1226  {
1227  AssertThrow(operation == MaterialAveraging::AveragingOperation::none,ExcNotImplemented());
1228  return;
1229  }
1230 
1236  std::vector<double> enthalpies_of_fusion;
1237  };
1238 
1239 
1240 
1256  template <int dim>
1257  class Interface : public Plugins::InterfaceBase
1258  {
1259  public:
1267 
1275 
1288  get_model_dependence () const;
1289 
1299  virtual bool is_compressible () const = 0;
1308  virtual
1309  void evaluate (const MaterialModel::MaterialModelInputs<dim> &in,
1311 
1317  virtual
1318  void
1319  create_additional_named_outputs (MaterialModelOutputs &outputs) const;
1320 
1321 
1329  virtual
1330  void
1331  fill_additional_material_model_inputs(MaterialModel::MaterialModelInputs<dim> &input,
1332  const LinearAlgebra::BlockVector &solution,
1333  const FEValuesBase<dim> &fe_values,
1334  const Introspection<dim> &introspection) const;
1335 
1336  protected:
1353  };
1354 
1370  template <int dim>
1371  void
1372  register_material_model (const std::string &name,
1373  const std::string &description,
1374  void (*declare_parameters_function) (ParameterHandler &),
1375  std::unique_ptr<Interface<dim>> (*factory_function) ());
1376 
1387  template <int dim>
1388  std::unique_ptr<Interface<dim>>
1389  create_material_model (const std::string &model_name);
1390 
1391 
1403  template <int dim>
1404  std::unique_ptr<Interface<dim>>
1405  create_material_model (ParameterHandler &prm);
1406 
1407 
1414  template <int dim>
1415  std::string
1417 
1418 
1424  template <int dim>
1425  void
1426  declare_parameters (ParameterHandler &prm);
1427 
1428 
1429 
1439  template <int dim>
1440  void
1441  write_plugin_graph (std::ostream &output_stream);
1442 
1443 
1444 
1445 // --------------------- template function definitions ----------------------------------
1446 
1447  template <int dim>
1448  template <class AdditionalInputType>
1450  {
1451  for (unsigned int i=0; i<additional_inputs.size(); ++i)
1452  {
1453  AdditionalInputType *result = dynamic_cast<AdditionalInputType *> (additional_inputs[i].get());
1454  if (result)
1455  return result;
1456  }
1457  return nullptr;
1458  }
1459 
1460 
1461  template <int dim>
1462  template <class AdditionalInputType>
1463  const AdditionalInputType *MaterialModelInputs<dim>::get_additional_input() const
1464  {
1465  for (unsigned int i=0; i<additional_inputs.size(); ++i)
1466  {
1467  const AdditionalInputType *result = dynamic_cast<const AdditionalInputType *> (additional_inputs[i].get());
1468  if (result)
1469  return result;
1470  }
1471  return nullptr;
1472  }
1473 
1474 
1475  template <int dim>
1476  template <class AdditionalOutputType>
1478  {
1479  for (unsigned int i=0; i<additional_outputs.size(); ++i)
1480  {
1481  AdditionalOutputType *result = dynamic_cast<AdditionalOutputType *> (additional_outputs[i].get());
1482  if (result)
1483  return result;
1484  }
1485  return nullptr;
1486  }
1487 
1488 
1489  template <int dim>
1490  template <class AdditionalOutputType>
1491  const AdditionalOutputType *MaterialModelOutputs<dim>::get_additional_output() const
1492  {
1493  for (unsigned int i=0; i<additional_outputs.size(); ++i)
1494  {
1495  const AdditionalOutputType *result = dynamic_cast<const AdditionalOutputType *> (additional_outputs[i].get());
1496  if (result)
1497  return result;
1498  }
1499  return nullptr;
1500  }
1501 
1502 
1503  template <int dim>
1505  {
1506  Assert(this->additional_outputs.empty(), ExcMessage("Destination of move needs to be empty!"));
1507  this->additional_outputs = std::move(other.additional_outputs);
1508  }
1509 
1510 
1518 #define ASPECT_REGISTER_MATERIAL_MODEL(classname,name,description) \
1519  template class classname<2>; \
1520  template class classname<3>; \
1521  namespace ASPECT_REGISTER_MATERIAL_MODEL_ ## classname \
1522  { \
1523  aspect::internal::Plugins::RegisterHelper<aspect::MaterialModel::Interface<2>,classname<2>> \
1524  dummy_ ## classname ## _2d (&aspect::MaterialModel::register_material_model<2>, \
1525  name, description); \
1526  aspect::internal::Plugins::RegisterHelper<aspect::MaterialModel::Interface<3>,classname<3>> \
1527  dummy_ ## classname ## _3d (&aspect::MaterialModel::register_material_model<3>, \
1528  name, description); \
1529  }
1530  }
1531 }
1532 
1533 
1534 #endif
void write_plugin_graph(std::ostream &output_stream)
std::vector< double > entropy_derivative_pressure
Definition: interface.h:570
std::vector< double > entropy_derivative_temperature
Definition: interface.h:578
void average_property(const AveragingOperation operation, const FullMatrix< double > &projection_matrix, const FullMatrix< double > &expansion_matrix, std::vector< double > &values_out)
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:239
std::vector< std::vector< double > > composition
Definition: interface.h:398
std::vector< double > thermal_expansion_coefficients
Definition: interface.h:547
virtual void average(const MaterialAveraging::AveragingOperation, const FullMatrix< double > &, const FullMatrix< double > &)
Definition: interface.h:835
std::vector< Point< dim > > position
Definition: interface.h:366
std::vector< std::vector< double > > prescribed_field_outputs
Definition: interface.h:1042
AdditionalInputType * get_additional_input()
Definition: interface.h:1449
void register_material_model(const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< Interface< dim >>(*factory_function)())
std::vector< std::unique_ptr< AdditionalMaterialOutputs< dim > > > additional_outputs
Definition: interface.h:623
std::vector< double > thermal_conductivities
Definition: interface.h:557
ElasticOutputs(const unsigned int n_points)
Definition: interface.h:1170
std::vector< std::vector< double > > reaction_rates
Definition: interface.h:989
void declare_parameters(ParameterHandler &prm)
bool identifies_single_variable(const Dependence dependence)
void average(const MaterialAveraging::AveragingOperation, const FullMatrix< double > &, const FullMatrix< double > &) override
Definition: interface.h:1093
std::vector< SymmetricTensor< 2, dim > > strain_rate
Definition: interface.h:412
std::vector< std::vector< double > > reaction_terms
Definition: interface.h:616
AdditionalOutputType * get_additional_output()
Definition: interface.h:1477
const std::vector< std::string > names
Definition: interface.h:920
Property operator|=(Property &d1, const Property d2)
Definition: interface.h:226
void average(const MaterialAveraging::AveragingOperation operation, const FullMatrix< double > &, const FullMatrix< double > &) override
Definition: interface.h:1223
std::vector< Tensor< 1, dim > > velocity
Definition: interface.h:391
EnthalpyOutputs(const unsigned int n_points)
Definition: interface.h:1216
AveragingOperation parse_averaging_operation_name(const std::string &s)
Property operator|(const Property d1, const Property d2)
Definition: interface.h:220
AdditionalMaterialOutputsStokesRHS(const unsigned int n_points)
Definition: interface.h:1086
void move_additional_outputs_from(MaterialModelOutputs< dim > &other)
Definition: interface.h:1504
std::vector< double > enthalpies_of_fusion
Definition: interface.h:1236
std::vector< SymmetricTensor< 2, dim > > elastic_force
Definition: interface.h:1191
std::vector< std::vector< double > > output_values
Definition: interface.h:917
std::vector< std::unique_ptr< AdditionalMaterialInputs< dim > > > additional_inputs
Definition: interface.h:464
void average(const MaterialAveraging::AveragingOperation operation, const FullMatrix< double > &, const FullMatrix< double > &) override
Definition: interface.h:1178
Definition: compat.h:59
DoFHandler< dim >::active_cell_iterator current_cell
Definition: interface.h:429
NonlinearDependence::ModelDependence model_dependence
Definition: interface.h:1352
std::unique_ptr< Interface< dim > > create_material_model(const std::string &model_name)
void average(const AveragingOperation operation, const typename DoFHandler< dim >::active_cell_iterator &cell, const Quadrature< dim > &quadrature_formula, const Mapping< dim > &mapping, const MaterialProperties::Property &requested_properties, MaterialModelOutputs< dim > &values_out)
std::string get_valid_model_names_pattern()
std::vector< Tensor< 1, dim > > pressure_gradient
Definition: interface.h:382
std::vector< SymmetricTensor< 2, dim > > viscoelastic_strain_rate
Definition: interface.h:1197
MaterialProperties::Property requested_properties
Definition: interface.h:438
Definition: compat.h:42
void average(const MaterialAveraging::AveragingOperation, const FullMatrix< double > &, const FullMatrix< double > &) override
Definition: interface.h:907