ASPECT
interface.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 #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  {
59  namespace NonlinearDependence
60  {
83  {
85 
86  none = 1,
88  pressure = 4,
91 
93  };
94 
95 
100  const Dependence d2)
101  {
102  return Dependence(static_cast<int>(d1) | static_cast<int>(d2));
103  }
104 
106  const Dependence d2)
107  {
108  d1 = (d1 | d2);
109  return d1;
110  }
111 
117  {
123 
129 
135 
141 
147 
153  ModelDependence ();
154  };
155 
164  bool
165  identifies_single_variable(const Dependence dependence);
166  }
167 
172  namespace MaterialProperties
173  {
183  enum Property
184  {
186 
187  none = 1,
189  density = 4,
199 
202  specific_heat |
207  viscosity |
212  };
213 
218  inline Property operator | (const Property d1,
219  const Property d2)
220  {
221  return Property(static_cast<int>(d1) | static_cast<int>(d2));
222  }
223 
225  const Property d2)
226  {
227  return (d1 | d2);
228  }
229  }
230 
231 
232  // Forward declaration:
233  template <int dim>
235 
242  template <int dim>
244  {
245  public:
256  MaterialModelInputs(const unsigned int n_points,
257  const unsigned int n_comp);
258 
274  MaterialModelInputs(const DataPostprocessorInputs::Vector<dim> &input_data,
275  const Introspection<dim> &introspection,
276  const bool compute_strain_rate = true);
277 
278 
298  MaterialModelInputs(const FEValuesBase<dim,dim> &fe_values,
299  const typename DoFHandler<dim>::active_cell_iterator &cell,
300  const Introspection<dim> &introspection,
301  const LinearAlgebra::BlockVector &solution_vector,
302  const bool compute_strain_rate = true);
303 
318 
322  MaterialModelInputs (MaterialModelInputs &&) noexcept = default;
323 
328  MaterialModelInputs &operator= (const MaterialModelInputs &source) = delete;
329 
333  MaterialModelInputs &operator= (MaterialModelInputs &&) = default;
334 
340  void reinit(const FEValuesBase<dim,dim> &fe_values,
341  const typename DoFHandler<dim>::active_cell_iterator &cell,
342  const Introspection<dim> &introspection,
343  const LinearAlgebra::BlockVector &solution_vector,
344  const bool compute_strain_rate = true);
345 
350  unsigned int n_evaluation_points() const;
351 
358  bool requests_property(const MaterialProperties::Property &property) const;
359 
364  std::vector<Point<dim>> position;
365 
369  std::vector<double> temperature;
370 
374  std::vector<double> pressure;
375 
380  std::vector<Tensor<1,dim>> pressure_gradient;
381 
389  std::vector<Tensor<1,dim>> velocity;
390 
396  std::vector<std::vector<double>> composition;
397 
410  std::vector<SymmetricTensor<2,dim>> strain_rate;
411 
427  typename DoFHandler<dim>::active_cell_iterator current_cell;
428 
437 
447  template <class AdditionalInputType>
448  AdditionalInputType *get_additional_input();
449 
454  template <class AdditionalInputType>
455  const AdditionalInputType *get_additional_input() const;
456 
462  std::vector<std::unique_ptr<AdditionalMaterialInputs<dim>>> additional_inputs;
463  };
464 
465 
466  // Forward declaration:
467  template <int dim>
469 
470 
477  template <int dim>
479  {
480  public:
491  MaterialModelOutputs (const unsigned int n_points,
492  const unsigned int n_comp);
493 
508 
512  MaterialModelOutputs (MaterialModelOutputs &&) noexcept = default;
513 
518  MaterialModelOutputs &operator= (const MaterialModelOutputs &source) = delete;
519 
523  MaterialModelOutputs &operator= (MaterialModelOutputs &&) noexcept = default;
524 
529  unsigned int n_evaluation_points() const;
530 
534  std::vector<double> viscosities;
535 
539  std::vector<double> densities;
540 
545  std::vector<double> thermal_expansion_coefficients;
546 
550  std::vector<double> specific_heat;
551 
555  std::vector<double> thermal_conductivities;
556 
561  std::vector<double> compressibilities;
562 
568  std::vector<double> entropy_derivative_pressure;
569 
576  std::vector<double> entropy_derivative_temperature;
577 
614  std::vector<std::vector<double>> reaction_terms;
615 
621  std::vector<std::unique_ptr<AdditionalMaterialOutputs<dim>>> additional_outputs;
622 
630  template <class AdditionalOutputType>
631  AdditionalOutputType *get_additional_output();
632 
636  template <class AdditionalOutputType>
637  const AdditionalOutputType *get_additional_output() const;
638 
643  void move_additional_outputs_from(MaterialModelOutputs<dim> &other);
644  };
645 
646 
647 
662  namespace MaterialAveraging
663  {
703  {
715  };
716 
717 
724  std::string get_averaging_operation_names ();
725 
731  AveragingOperation parse_averaging_operation_name (const std::string &s);
732 
738  template <int dim>
739  void average (const AveragingOperation operation,
740  const typename DoFHandler<dim>::active_cell_iterator &cell,
741  const Quadrature<dim> &quadrature_formula,
742  const Mapping<dim> &mapping,
743  const MaterialProperties::Property &requested_properties,
744  MaterialModelOutputs<dim> &values_out);
745 
751  void average_property (const AveragingOperation operation,
752  const FullMatrix<double> &projection_matrix,
753  const FullMatrix<double> &expansion_matrix,
754  std::vector<double> &values_out);
755 
766  }
767 
786  template <int dim>
788  {
789  public:
794  virtual ~AdditionalMaterialInputs() = default;
795 
800  virtual void
801  fill (const LinearAlgebra::BlockVector &solution,
802  const FEValuesBase<dim> &fe_values,
803  const Introspection<dim> &introspection) = 0;
804  };
805 
806 
833  template <int dim>
835  {
836  public:
841  virtual ~AdditionalMaterialOutputs() = default;
842 
843  virtual void average (const MaterialAveraging::AveragingOperation /*operation*/,
844  const FullMatrix<double> &/*projection_matrix*/,
845  const FullMatrix<double> &/*expansion_matrix*/)
846  {}
847  };
848 
849 
871  template <int dim>
873  {
874  public:
883  NamedAdditionalMaterialOutputs(const std::vector<std::string> &output_names);
884 
895  NamedAdditionalMaterialOutputs(const std::vector<std::string> &output_names,
896  const unsigned int n_points);
897 
901  ~NamedAdditionalMaterialOutputs() override;
902 
907  const std::vector<std::string> &get_names() const;
908 
913  virtual std::vector<double> get_nth_output(const unsigned int idx) const;
914 
916  const FullMatrix<double> &/*projection_matrix*/,
917  const FullMatrix<double> &/*expansion_matrix*/) override
918  {}
919 
920 
925  std::vector<std::vector<double>> output_values;
926 
927  private:
928  const std::vector<std::string> names;
929  };
930 
931 
937  template <int dim>
939  {
940  public:
941  SeismicAdditionalOutputs(const unsigned int n_points);
942 
943  std::vector<double> get_nth_output(const unsigned int idx) const override;
944 
950  std::vector<double> vs;
951 
957  std::vector<double> vp;
958  };
959 
960 
981  template <int dim>
983  {
984  public:
985  ReactionRateOutputs (const unsigned int n_points,
986  const unsigned int n_comp);
987 
988  std::vector<double> get_nth_output(const unsigned int idx) const override;
989 
997  std::vector<std::vector<double>> reaction_rates;
998  };
999 
1000 
1001 
1007  template <int dim>
1009  {
1010  public:
1011  PhaseOutputs(const unsigned int n_points);
1012  };
1013 
1014 
1015 
1034  template <int dim>
1036  {
1037  public:
1038  PrescribedFieldOutputs (const unsigned int n_points,
1039  const unsigned int n_comp);
1040 
1041  std::vector<double> get_nth_output(const unsigned int idx) const override;
1042 
1050  std::vector<std::vector<double>> prescribed_field_outputs;
1051  };
1052 
1066  template <int dim>
1068  {
1069  public:
1070  PrescribedTemperatureOutputs (const unsigned int n_points);
1071 
1072  std::vector<double> get_nth_output(const unsigned int idx) const override;
1073 
1081  std::vector<double> prescribed_temperature_outputs;
1082  };
1083 
1090  template <int dim>
1092  {
1093  public:
1094  AdditionalMaterialOutputsStokesRHS(const unsigned int n_points)
1095  : rhs_u(n_points), rhs_p(n_points), rhs_melt_pc(n_points)
1096  {}
1097 
1099  = default;
1100 
1102  const FullMatrix<double> &/*projection_matrix*/,
1103  const FullMatrix<double> &/*expansion_matrix*/) override
1104  {
1105  // TODO: not implemented
1106  }
1107 
1113  std::vector<Tensor<1,dim>> rhs_u;
1114 
1119  std::vector<double> rhs_p;
1120 
1125  std::vector<double> rhs_melt_pc;
1126  };
1127 
1128 
1129 
1145  template <int dim>
1147  {
1148  public:
1152  explicit PrescribedPlasticDilation (const unsigned int n_points);
1153 
1157  std::vector<double> get_nth_output(const unsigned int idx) const override;
1158 
1163  std::vector<double> dilation;
1164  };
1165 
1166 
1167 
1174  template <int dim>
1176  {
1177  public:
1178  ElasticOutputs(const unsigned int n_points)
1179  : elastic_force(n_points, numbers::signaling_nan<SymmetricTensor<2,dim>>())
1180  , viscoelastic_strain_rate(n_points, numbers::signaling_nan<SymmetricTensor<2,dim>>())
1181  {}
1182 
1183  ~ElasticOutputs() override
1184  = default;
1185 
1187  const FullMatrix<double> &/*projection_matrix*/,
1188  const FullMatrix<double> &/*expansion_matrix*/) override
1189  {
1190  AssertThrow(operation == MaterialAveraging::AveragingOperation::none,ExcNotImplemented());
1191  return;
1192  }
1193 
1199  std::vector<SymmetricTensor<2,dim>> elastic_force;
1200 
1205  std::vector<SymmetricTensor<2,dim>> viscoelastic_strain_rate;
1206  };
1207 
1208 
1209 
1220  template <int dim>
1222  {
1223  public:
1224  EnthalpyOutputs(const unsigned int n_points)
1225  : enthalpies_of_fusion(n_points, numbers::signaling_nan<double>())
1226  {}
1227 
1228  virtual ~EnthalpyOutputs()
1229  = default;
1230 
1232  const FullMatrix<double> &/*projection_matrix*/,
1233  const FullMatrix<double> &/*expansion_matrix*/) override
1234  {
1235  AssertThrow(operation == MaterialAveraging::AveragingOperation::none,ExcNotImplemented());
1236  return;
1237  }
1238 
1244  std::vector<double> enthalpies_of_fusion;
1245  };
1246 
1247 
1248 
1264  template <int dim>
1265  class Interface : public Plugins::InterfaceBase
1266  {
1267  public:
1275 
1283 
1296  get_model_dependence () const;
1297 
1307  virtual bool is_compressible () const = 0;
1316  virtual
1317  void evaluate (const MaterialModel::MaterialModelInputs<dim> &in,
1319 
1325  virtual
1326  void
1327  create_additional_named_outputs (MaterialModelOutputs &outputs) const;
1328 
1329 
1337  virtual
1338  void
1339  fill_additional_material_model_inputs(MaterialModel::MaterialModelInputs<dim> &input,
1340  const LinearAlgebra::BlockVector &solution,
1341  const FEValuesBase<dim> &fe_values,
1342  const Introspection<dim> &introspection) const;
1343 
1344  protected:
1361  };
1362 
1378  template <int dim>
1379  void
1380  register_material_model (const std::string &name,
1381  const std::string &description,
1382  void (*declare_parameters_function) (ParameterHandler &),
1383  std::unique_ptr<Interface<dim>> (*factory_function) ());
1384 
1395  template <int dim>
1396  std::unique_ptr<Interface<dim>>
1397  create_material_model (const std::string &model_name);
1398 
1399 
1411  template <int dim>
1412  std::unique_ptr<Interface<dim>>
1413  create_material_model (ParameterHandler &prm);
1414 
1415 
1422  template <int dim>
1423  std::string
1425 
1426 
1432  template <int dim>
1433  void
1434  declare_parameters (ParameterHandler &prm);
1435 
1436 
1437 
1447  template <int dim>
1448  void
1449  write_plugin_graph (std::ostream &output_stream);
1450 
1451 
1452 
1453 // --------------------- template function definitions ----------------------------------
1454 
1455  template <int dim>
1456  template <class AdditionalInputType>
1458  {
1459  for (unsigned int i=0; i<additional_inputs.size(); ++i)
1460  {
1461  AdditionalInputType *result = dynamic_cast<AdditionalInputType *> (additional_inputs[i].get());
1462  if (result)
1463  return result;
1464  }
1465  return nullptr;
1466  }
1467 
1468 
1469  template <int dim>
1470  template <class AdditionalInputType>
1471  const AdditionalInputType *MaterialModelInputs<dim>::get_additional_input() const
1472  {
1473  for (unsigned int i=0; i<additional_inputs.size(); ++i)
1474  {
1475  const AdditionalInputType *result = dynamic_cast<const AdditionalInputType *> (additional_inputs[i].get());
1476  if (result)
1477  return result;
1478  }
1479  return nullptr;
1480  }
1481 
1482 
1483  template <int dim>
1484  template <class AdditionalOutputType>
1486  {
1487  for (unsigned int i=0; i<additional_outputs.size(); ++i)
1488  {
1489  AdditionalOutputType *result = dynamic_cast<AdditionalOutputType *> (additional_outputs[i].get());
1490  if (result)
1491  return result;
1492  }
1493  return nullptr;
1494  }
1495 
1496 
1497  template <int dim>
1498  template <class AdditionalOutputType>
1499  const AdditionalOutputType *MaterialModelOutputs<dim>::get_additional_output() const
1500  {
1501  for (unsigned int i=0; i<additional_outputs.size(); ++i)
1502  {
1503  const AdditionalOutputType *result = dynamic_cast<const AdditionalOutputType *> (additional_outputs[i].get());
1504  if (result)
1505  return result;
1506  }
1507  return nullptr;
1508  }
1509 
1510 
1511  template <int dim>
1513  {
1514  Assert(this->additional_outputs.empty(), ExcMessage("Destination of move needs to be empty!"));
1515  this->additional_outputs = std::move(other.additional_outputs);
1516  }
1517 
1518 
1526 #define ASPECT_REGISTER_MATERIAL_MODEL(classname,name,description) \
1527  template class classname<2>; \
1528  template class classname<3>; \
1529  namespace ASPECT_REGISTER_MATERIAL_MODEL_ ## classname \
1530  { \
1531  aspect::internal::Plugins::RegisterHelper<aspect::MaterialModel::Interface<2>,classname<2>> \
1532  dummy_ ## classname ## _2d (&aspect::MaterialModel::register_material_model<2>, \
1533  name, description); \
1534  aspect::internal::Plugins::RegisterHelper<aspect::MaterialModel::Interface<3>,classname<3>> \
1535  dummy_ ## classname ## _3d (&aspect::MaterialModel::register_material_model<3>, \
1536  name, description); \
1537  }
1538  }
1539 }
1540 
1541 
1542 #endif
void write_plugin_graph(std::ostream &output_stream)
AveragingOperation get_averaging_operation_for_viscosity(const AveragingOperation operation)
std::vector< double > entropy_derivative_pressure
Definition: interface.h:568
std::vector< double > entropy_derivative_temperature
Definition: interface.h:576
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:268
std::vector< std::vector< double > > composition
Definition: interface.h:396
std::vector< double > thermal_expansion_coefficients
Definition: interface.h:545
virtual void average(const MaterialAveraging::AveragingOperation, const FullMatrix< double > &, const FullMatrix< double > &)
Definition: interface.h:843
std::vector< Point< dim > > position
Definition: interface.h:364
std::vector< std::vector< double > > prescribed_field_outputs
Definition: interface.h:1050
AdditionalInputType * get_additional_input()
Definition: interface.h:1457
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:621
std::vector< double > thermal_conductivities
Definition: interface.h:555
ElasticOutputs(const unsigned int n_points)
Definition: interface.h:1178
std::vector< std::vector< double > > reaction_rates
Definition: interface.h:997
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:1101
std::vector< SymmetricTensor< 2, dim > > strain_rate
Definition: interface.h:410
std::vector< std::vector< double > > reaction_terms
Definition: interface.h:614
AdditionalOutputType * get_additional_output()
Definition: interface.h:1485
const std::vector< std::string > names
Definition: interface.h:928
void average(const MaterialAveraging::AveragingOperation operation, const FullMatrix< double > &, const FullMatrix< double > &) override
Definition: interface.h:1231
std::vector< Tensor< 1, dim > > velocity
Definition: interface.h:389
EnthalpyOutputs(const unsigned int n_points)
Definition: interface.h:1224
AveragingOperation parse_averaging_operation_name(const std::string &s)
AdditionalMaterialOutputsStokesRHS(const unsigned int n_points)
Definition: interface.h:1094
Dependence operator|=(Dependence &d1, const Dependence d2)
Definition: interface.h:105
void move_additional_outputs_from(MaterialModelOutputs< dim > &other)
Definition: interface.h:1512
std::vector< double > enthalpies_of_fusion
Definition: interface.h:1244
std::vector< SymmetricTensor< 2, dim > > elastic_force
Definition: interface.h:1199
std::vector< std::vector< double > > output_values
Definition: interface.h:925
std::vector< std::unique_ptr< AdditionalMaterialInputs< dim > > > additional_inputs
Definition: interface.h:462
void average(const MaterialAveraging::AveragingOperation operation, const FullMatrix< double > &, const FullMatrix< double > &) override
Definition: interface.h:1186
Definition: compat.h:59
DoFHandler< dim >::active_cell_iterator current_cell
Definition: interface.h:427
NonlinearDependence::ModelDependence model_dependence
Definition: interface.h:1360
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()
Dependence operator|(const Dependence d1, const Dependence d2)
Definition: interface.h:99
std::vector< Tensor< 1, dim > > pressure_gradient
Definition: interface.h:380
std::vector< SymmetricTensor< 2, dim > > viscoelastic_strain_rate
Definition: interface.h:1205
MaterialProperties::Property requested_properties
Definition: interface.h:436
void average(const MaterialAveraging::AveragingOperation, const FullMatrix< double > &, const FullMatrix< double > &) override
Definition: interface.h:915