24 #ifndef PDE_INTERFACE_H_
25 #define PDE_INTERFACE_H_
30 #include <deal.II/fe/fe_system.h>
31 #include <deal.II/fe/fe_values.h>
32 #include <deal.II/fe/mapping.h>
33 #include <deal.II/lac/full_matrix.h>
34 #include <deal.II/base/function.h>
63 template<
template<
int,
int>
class DH,
typename VECTOR,
int dealdim>
class EDC,
64 template<template<int, int> class DH, typename VECTOR, int dealdim> class FDC,
65 template<int, int> class DH, typename VECTOR, int dealdim>
100 ElementEquation(
const EDC<DH, VECTOR, dealdim>& ,
101 dealii::Vector<double> &,
124 StrongElementResidual(
const EDC<DH, VECTOR, dealdim>& ,
125 const EDC<DH, VECTOR, dealdim>& ,
153 ElementTimeEquation(
const EDC<DH, VECTOR, dealdim>& ,
154 dealii::Vector<double> &,
178 ElementTimeEquation_U(
const EDC<DH, VECTOR, dealdim>& ,
179 dealii::Vector<double> &,
200 ElementTimeEquation_UT(
const EDC<DH, VECTOR, dealdim>& ,
201 dealii::Vector<double> &,
225 ElementTimeEquation_UTT(
const EDC<DH, VECTOR, dealdim>& ,
226 dealii::Vector<double> &,
249 ElementTimeEquationExplicit(
const EDC<DH, VECTOR, dealdim>& ,
250 dealii::Vector<double> &,
267 ElementTimeEquationExplicit_U(
const EDC<DH, VECTOR, dealdim>& ,
268 dealii::Vector<double> &,
287 ElementTimeEquationExplicit_UT(
const EDC<DH, VECTOR, dealdim>& ,
288 dealii::Vector<double> &,
306 ElementTimeEquationExplicit_UTT(
const EDC<DH, VECTOR, dealdim>& ,
307 dealii::Vector<double> &,
327 ElementTimeEquationExplicit_UU(
const EDC<DH, VECTOR, dealdim>& ,
328 dealii::Vector<double> &,
356 ElementEquation_U(
const EDC<DH, VECTOR, dealdim>& ,
357 dealii::Vector<double> &,
379 StrongElementResidual_U(
const EDC<DH, VECTOR, dealdim>& ,
380 const EDC<DH, VECTOR, dealdim>& ,
409 ElementEquation_UT(
const EDC<DH, VECTOR, dealdim>& ,
410 dealii::Vector<double> &,
441 ElementEquation_UTT(
const EDC<DH, VECTOR, dealdim>& ,
442 dealii::Vector<double> &,
474 ElementEquation_Q(
const EDC<DH, VECTOR, dealdim>& ,
475 dealii::Vector<double> &,
506 ElementEquation_QT(
const EDC<DH, VECTOR, dealdim>& ,
507 dealii::Vector<double> &,
534 ElementEquation_QTT(
const EDC<DH, VECTOR, dealdim>& ,
535 dealii::Vector<double> &,
563 ElementEquation_UU(
const EDC<DH, VECTOR, dealdim>& ,
564 dealii::Vector<double> &,
594 ElementEquation_QU(
const EDC<DH, VECTOR, dealdim>& ,
595 dealii::Vector<double> &,
622 ElementEquation_UQ(
const EDC<DH, VECTOR, dealdim>& ,
623 dealii::Vector<double> &,
649 ElementEquation_QQ(
const EDC<DH, VECTOR, dealdim>& ,
650 dealii::Vector<double> &,
669 ElementRightHandSide(
const EDC<DH, VECTOR, dealdim>& ,
670 dealii::Vector<double> &,
697 ElementMatrix(
const EDC<DH, VECTOR, dealdim>& ,
698 dealii::FullMatrix<double> &,
717 ElementTimeMatrix(
const EDC<DH, VECTOR, dealdim>& ,
718 dealii::FullMatrix<double> &);
732 ElementTimeMatrix_T(
const EDC<DH, VECTOR, dealdim>& ,
733 dealii::FullMatrix<double> &);
750 ElementTimeMatrixExplicit(
const EDC<DH, VECTOR, dealdim>& ,
751 dealii::FullMatrix<double> &);
765 ElementTimeMatrixExplicit_T(
const EDC<DH, VECTOR, dealdim>& ,
766 dealii::FullMatrix<double> &);
794 ElementMatrix_T(
const EDC<DH, VECTOR, dealdim>& ,
795 dealii::FullMatrix<double> &,
817 ControlElementEquation(
const EDC<DH, VECTOR, dealdim>& ,
818 dealii::Vector<double> &,
833 ControlElementMatrix(
const EDC<DH, VECTOR, dealdim>& ,
834 dealii::FullMatrix<double> &,
854 ControlBoundaryEquation(
const FDC<DH, VECTOR, dealdim>& ,
855 dealii::Vector<double> &,
870 ControlBoundaryMatrix(
const FDC<DH, VECTOR, dealdim>& ,
871 dealii::FullMatrix<double> &,
891 StrongElementResidual_Control(
const EDC<DH, VECTOR, dealdim>& ,
892 const EDC<DH, VECTOR, dealdim>& ,
914 StrongFaceResidual_Control(
const FDC<DH, VECTOR, dealdim>& ,
915 const FDC<DH, VECTOR, dealdim>& ,
937 StrongBoundaryResidual_Control(
const FDC<DH, VECTOR, dealdim>& ,
938 const FDC<DH, VECTOR, dealdim>& ,
954 FaceEquation(
const FDC<DH, VECTOR, dealdim>& ,
955 dealii::Vector<double> &,
961 StrongFaceResidual(
const FDC<DH, VECTOR, dealdim>& ,
962 const FDC<DH, VECTOR, dealdim>& ,
969 FaceEquation_U(
const FDC<DH, VECTOR, dealdim>& ,
970 dealii::Vector<double> &,
977 StrongFaceResidual_U(
const FDC<DH, VECTOR, dealdim>& ,
978 const FDC<DH, VECTOR, dealdim>& ,
985 FaceEquation_UT(
const FDC<DH, VECTOR, dealdim>& ,
986 dealii::Vector<double> &,
993 FaceEquation_UTT(
const FDC<DH, VECTOR, dealdim>& ,
994 dealii::Vector<double> &,
1001 FaceEquation_Q(
const FDC<DH, VECTOR, dealdim>& ,
1002 dealii::Vector<double> &,
1009 FaceEquation_QT(
const FDC<DH, VECTOR, dealdim>& ,
1010 dealii::Vector<double> &,
1017 FaceEquation_QTT(
const FDC<DH, VECTOR, dealdim>& ,
1018 dealii::Vector<double> &,
1025 FaceEquation_UU(
const FDC<DH, VECTOR, dealdim>& ,
1026 dealii::Vector<double> &,
1033 FaceEquation_QU(
const FDC<DH, VECTOR, dealdim>& ,
1034 dealii::Vector<double> &,
1041 FaceEquation_UQ(
const FDC<DH, VECTOR, dealdim>& ,
1042 dealii::Vector<double> &,
1049 FaceEquation_QQ(
const FDC<DH, VECTOR, dealdim>& ,
1050 dealii::Vector<double> &,
1057 FaceRightHandSide(
const FDC<DH, VECTOR, dealdim>& ,
1058 dealii::Vector<double> &,
1067 FaceMatrix(
const FDC<DH, VECTOR, dealdim>& ,
1068 dealii::FullMatrix<double> &,
1075 FaceMatrix_T(
const FDC<DH, VECTOR, dealdim>& ,
1076 dealii::FullMatrix<double> &,
1084 InterfaceMatrix(
const FDC<DH, VECTOR, dealdim>& ,
1085 dealii::FullMatrix<double> &,
1092 InterfaceMatrix_T(
const FDC<DH, VECTOR, dealdim>& ,
1093 dealii::FullMatrix<double> &,
1100 InterfaceEquation(
const FDC<DH, VECTOR, dealdim>& ,
1101 dealii::Vector<double> &,
1108 InterfaceEquation_U(
const FDC<DH, VECTOR, dealdim>& ,
1109 dealii::Vector<double> &,
1117 BoundaryEquation(
const FDC<DH, VECTOR, dealdim>& ,
1118 dealii::Vector<double> &,
1125 StrongBoundaryResidual(
const FDC<DH, VECTOR, dealdim>& ,
1126 const FDC<DH, VECTOR, dealdim>& ,
1133 BoundaryEquation_U(
const FDC<DH, VECTOR, dealdim>&,
1134 dealii::Vector<double> &,
1141 StrongBoundaryResidual_U(
const FDC<DH, VECTOR, dealdim>& ,
1142 const FDC<DH, VECTOR, dealdim>& ,
1149 BoundaryEquation_UT(
const FDC<DH, VECTOR, dealdim>& ,
1150 dealii::Vector<double> &,
1157 BoundaryEquation_UTT(
const FDC<DH, VECTOR, dealdim>& ,
1158 dealii::Vector<double> &,
1165 BoundaryEquation_Q(
const FDC<DH, VECTOR, dealdim>& ,
1166 dealii::Vector<double> &,
1173 BoundaryEquation_QT(
const FDC<DH, VECTOR, dealdim>& ,
1174 dealii::Vector<double> &,
1181 BoundaryEquation_QTT(
const FDC<DH, VECTOR, dealdim>& ,
1182 dealii::Vector<double> &,
1189 BoundaryEquation_UU(
const FDC<DH, VECTOR, dealdim>& ,
1190 dealii::Vector<double> &,
1197 BoundaryEquation_QU(
const FDC<DH, VECTOR, dealdim>& ,
1198 dealii::Vector<double> &,
1205 BoundaryEquation_UQ(
const FDC<DH, VECTOR, dealdim>& ,
1206 dealii::Vector<double> &,
1213 BoundaryEquation_QQ(
const FDC<DH, VECTOR, dealdim>& ,
1214 dealii::Vector<double> &,
1221 BoundaryRightHandSide(
const FDC<DH, VECTOR, dealdim>& ,
1222 dealii::Vector<double> &,
1228 BoundaryMatrix(
const FDC<DH, VECTOR, dealdim>& ,
1229 dealii::FullMatrix<double> &,
1236 BoundaryMatrix_T(
const FDC<DH, VECTOR, dealdim>& ,
1237 dealii::FullMatrix<double> &,
1258 dealii::Vector<double> &local_vector,
1263 edc.GetFEValuesState();
1264 unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1265 unsigned int n_q_points = edc.GetNQPoints();
1266 std::vector<dealii::Vector<double> > uvalues;
1267 uvalues.resize(n_q_points,
1268 dealii::Vector<double>(this->GetStateNComponents()));
1269 edc.GetValuesState(
"last_newton_solution", uvalues);
1271 dealii::Vector<double> f_values(
1272 dealii::Vector<double>(this->GetStateNComponents()));
1274 for (
unsigned int q_point = 0; q_point < n_q_points; q_point++)
1276 for (
unsigned int i = 0; i < n_dofs_per_element; i++)
1278 for (
unsigned int comp = 0; comp < this->GetStateNComponents();
1281 local_vector(i) += scale
1282 * (state_fe_values.shape_value_component(i, q_point, comp)
1283 * uvalues[q_point](comp))
1284 * state_fe_values.JxW(q_point);
1292 dealii::Vector<double> &,
double )
1298 dealii::Vector<double> &,
double )
1304 dealii::Vector<double> &,
double )
1310 dealii::Vector<double> &,
double )
1317 const EDC<DH, VECTOR, dealdim>& edc,
1318 dealii::Vector<double> &local_vector,
1322 edc.GetFEValuesState();
1323 unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1324 unsigned int n_q_points = edc.GetNQPoints();
1326 dealii::Vector<double> f_values(
1327 dealii::Vector<double>(this->GetStateNComponents()));
1329 for (
unsigned int q_point = 0; q_point < n_q_points; q_point++)
1331 init_values->vector_value(state_fe_values.quadrature_point(q_point),
1334 for (
unsigned int i = 0; i < n_dofs_per_element; i++)
1336 for (
unsigned int comp = 0; comp < this->GetStateNComponents();
1339 local_vector(i) += scale
1341 * state_fe_values.shape_value_component(i, q_point,
1342 comp)) * state_fe_values.JxW(q_point);
1350 dealii::FullMatrix<double> &local_entry_matrix,
1355 edc.GetFEValuesState();
1356 unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1357 unsigned int n_q_points = edc.GetNQPoints();
1359 for (
unsigned int q_point = 0; q_point < n_q_points; q_point++)
1361 for (
unsigned int i = 0; i < n_dofs_per_element; i++)
1363 for (
unsigned int j = 0; j < n_dofs_per_element; j++)
1365 for (
unsigned int comp = 0; comp < this->GetStateNComponents();
1368 local_entry_matrix(i, j) += scale
1369 * state_fe_values.shape_value_component(i, q_point, comp)
1370 * state_fe_values.shape_value_component(j, q_point, comp)
1371 * state_fe_values.JxW(q_point);
1383 virtual dealii::UpdateFlags
1384 GetUpdateFlags()
const;
1390 virtual dealii::UpdateFlags
1391 GetFaceUpdateFlags()
const;
1402 HasInterfaces()
const;
1407 SetProblemType(std::string type);
1411 virtual unsigned int
1412 GetControlNBlocks()
const;
1413 virtual unsigned int
1414 GetStateNBlocks()
const;
1415 virtual std::vector<unsigned int>&
1416 GetControlBlockComponent();
1417 virtual const std::vector<unsigned int>&
1418 GetControlBlockComponent()
const;
1419 virtual std::vector<unsigned int>&
1420 GetStateBlockComponent();
1421 virtual const std::vector<unsigned int>&
1422 GetStateBlockComponent()
const;
1434 GetStateNComponents()
const;
1459 template<
typename ELEMENTITERATOR>
1463 if (element[0]->neighbor_index(face) != -1)
1464 if (element[0]->material_id()
1465 != element[0]->neighbor(face)->material_id())
virtual void Init_ElementRhs(const dealii::Function< dealdim > *init_values, const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
Definition: pdeinterface.h:1316
bool AtInterface(ELEMENTITERATOR &element, unsigned int face) const
Definition: pdeinterface.h:1461
boost::function1< void, dealii::Vector< double > & > VectorResidualModifier
Definition: pdeinterface.h:1444
virtual void Init_ElementRhs_QT(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1297
virtual void Init_ElementEquation(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale, double)
Definition: pdeinterface.h:1257
Definition: pdeinterface.h:66
boost::function1< void, double & > ResidualModifier
Definition: pdeinterface.h:1443
virtual void Init_ElementRhs_QTT(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1303
virtual void Init_ElementRhs_QQ(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1309
virtual void Init_ElementMatrix(const EDC< DH, VECTOR, dealdim > &edc, dealii::FullMatrix< double > &local_entry_matrix, double scale, double)
Definition: pdeinterface.h:1349
std::string problem_type_
Definition: pdeinterface.h:1471
Definition: fevalues_wrapper.h:43
virtual void SetTime(double) const
Definition: pdeinterface.h:1427
virtual void Init_ElementRhs_Q(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1291