DOpE
pdeinterface.h
Go to the documentation of this file.
1 
24 #ifndef PDE_INTERFACE_H_
25 #define PDE_INTERFACE_H_
26 
27 #include <map>
28 #include <string>
29 
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>
35 
36 #include "fevalues_wrapper.h"
37 #include "elementdatacontainer.h"
38 #include "facedatacontainer.h"
41 
42 namespace DOpE
43 {
44 
62  template<
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>
67  {
68  public:
69  PDEInterface();
70  virtual
71  ~PDEInterface();
72 
73  /******************************************************/
74 
99  virtual void
100  ElementEquation(const EDC<DH, VECTOR, dealdim>& /*edc*/,
101  dealii::Vector<double> &/*local_vector*/,
102  double /*scale*/,
103  double /*scale_ico*/);
104 
105  /******************************************************/
106 
123  virtual void
124  StrongElementResidual(const EDC<DH, VECTOR, dealdim>& /*edc*/,
125  const EDC<DH, VECTOR, dealdim>& /*edc_weight*/,
126  double& /*ret*/,
127  double /*scale*/);
128 
129  /******************************************************/
130 
151  //Note that the _UU term is not needed, since we assume that ElementTimeEquation is linear!
152  virtual void
153  ElementTimeEquation(const EDC<DH, VECTOR, dealdim>& /*edc*/,
154  dealii::Vector<double> &/*local_vector*/,
155  double /*scale*/);
156 
157  /******************************************************/
177  virtual void
178  ElementTimeEquation_U(const EDC<DH, VECTOR, dealdim>& /*edc*/,
179  dealii::Vector<double> &/*local_vector*/,
180  double /*scale*/);
181 
182  /******************************************************/
199  virtual void
200  ElementTimeEquation_UT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
201  dealii::Vector<double> &/*local_vector*/,
202  double /*scale*/);
203 
204  /******************************************************/
224  virtual void
225  ElementTimeEquation_UTT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
226  dealii::Vector<double> &/*local_vector*/,
227  double /*scale*/);
228 
229  /******************************************************/
230 
248  virtual void
249  ElementTimeEquationExplicit(const EDC<DH, VECTOR, dealdim>& /*edc**/,
250  dealii::Vector<double> &/*local_vector*/,
251  double /*scale*/);
252  /******************************************************/
266  virtual void
267  ElementTimeEquationExplicit_U(const EDC<DH, VECTOR, dealdim>& /*edc*/,
268  dealii::Vector<double> &/*local_vector*/,
269  double /*scale*/);
270 
271  /******************************************************/
286  virtual void
287  ElementTimeEquationExplicit_UT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
288  dealii::Vector<double> &/*local_vector*/,
289  double /*scale*/);
290 
291  /******************************************************/
305  virtual void
306  ElementTimeEquationExplicit_UTT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
307  dealii::Vector<double> &/*local_vector*/,
308  double /*scale*/);
309 
310  /******************************************************/
326  virtual void
327  ElementTimeEquationExplicit_UU(const EDC<DH, VECTOR, dealdim>& /*edc*/,
328  dealii::Vector<double> &/*local_vector*/,
329  double /*scale*/);
330 
331  /******************************************************/
355  virtual void
356  ElementEquation_U(const EDC<DH, VECTOR, dealdim>& /*edc*/,
357  dealii::Vector<double> &/*local_vector*/,
358  double /*scale*/,
359  double /*scale_ico*/);
360 
361  /******************************************************/
378  virtual void
379  StrongElementResidual_U(const EDC<DH, VECTOR, dealdim>& /*edc*/,
380  const EDC<DH, VECTOR, dealdim>& /*edc_weight*/,
381  double& /*ret*/,
382  double /*scale*/);
383 
384  /******************************************************/
408  virtual void
409  ElementEquation_UT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
410  dealii::Vector<double> &/*local_vector*/,
411  double /*scale*/,
412  double /*scale_ico*/);
413 
414  /******************************************************/
440  virtual void
441  ElementEquation_UTT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
442  dealii::Vector<double> &/*local_vector*/,
443  double /*scale*/,
444  double /*scale_ico*/);
445 
446  /******************************************************/
447 
473  virtual void
474  ElementEquation_Q(const EDC<DH, VECTOR, dealdim>& /*edc*/,
475  dealii::Vector<double> &/*local_vector*/,
476  double /*scale*/,
477  double /*scale_ico*/);
478 
479  /******************************************************/
480 
505  virtual void
506  ElementEquation_QT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
507  dealii::Vector<double> &/*local_vector*/,
508  double /*scale*/,
509  double scale_ico);
510 
511  /******************************************************/
512 
533  virtual void
534  ElementEquation_QTT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
535  dealii::Vector<double> &/*local_vector*/,
536  double /*scale*/,
537  double /*scale_ico*/);
538 
539  /******************************************************/
540 
562  virtual void
563  ElementEquation_UU(const EDC<DH, VECTOR, dealdim>& /*edc*/,
564  dealii::Vector<double> &/*local_vector*/,
565  double /*scale*/,
566  double /*scale_ico*/);
567 
568  /******************************************************/
569 
593  virtual void
594  ElementEquation_QU(const EDC<DH, VECTOR, dealdim>& /*edc*/,
595  dealii::Vector<double> &/*local_vector*/,
596  double /*scale*/,
597  double /*scale_ico*/);
598 
599  /******************************************************/
621  virtual void
622  ElementEquation_UQ(const EDC<DH, VECTOR, dealdim>& /*edc*/,
623  dealii::Vector<double> &/*local_vector*/,
624  double /*scale*/,
625  double /*scale_ico*/);
626 
627  /******************************************************/
648  virtual void
649  ElementEquation_QQ(const EDC<DH, VECTOR, dealdim>& /*edc*/,
650  dealii::Vector<double> &/*local_vector*/,
651  double /*scale*/,
652  double /*scale_ico*/);
653 
654  /******************************************************/
668  virtual void
669  ElementRightHandSide(const EDC<DH, VECTOR, dealdim>& /*edc*/,
670  dealii::Vector<double> &/*local_vector*/,
671  double /*scale*/);
672 
673  /******************************************************/
674 
696  virtual void
697  ElementMatrix(const EDC<DH, VECTOR, dealdim>& /*edc*/,
698  dealii::FullMatrix<double> &/*local_entry_matrix*/,
699  double /*scale*/,
700  double /*scale_ico*/);
701 
702  /******************************************************/
716  virtual void
717  ElementTimeMatrix(const EDC<DH, VECTOR, dealdim>& /*edc*/,
718  dealii::FullMatrix<double> &/*local_entry_matrix*/);
719 
720  /******************************************************/
731  virtual void
732  ElementTimeMatrix_T(const EDC<DH, VECTOR, dealdim>& /*edc*/,
733  dealii::FullMatrix<double> &/*local_entry_matrix*/);
734 
735  /******************************************************/
749  virtual void
750  ElementTimeMatrixExplicit(const EDC<DH, VECTOR, dealdim>& /*edc*/,
751  dealii::FullMatrix<double> &/*local_entry_matrix*/);
752 
753  /******************************************************/
764  virtual void
765  ElementTimeMatrixExplicit_T(const EDC<DH, VECTOR, dealdim>& /*edc*/,
766  dealii::FullMatrix<double> &/*local_entry_matrix*/);
767 
768  /******************************************************/
793  virtual void
794  ElementMatrix_T(const EDC<DH, VECTOR, dealdim>& /*edc*/,
795  dealii::FullMatrix<double> &/*local_entry_matrix*/,
796  double /*scale*/,
797  double /*scale_ico*/);
798 
799  /******************************************************/
800 
816  virtual void
817  ControlElementEquation(const EDC<DH, VECTOR, dealdim>& /*edc*/,
818  dealii::Vector<double> &/*local_vector*/,
819  double /*scale*/);
820 
821  /******************************************************/
832  virtual void
833  ControlElementMatrix(const EDC<DH, VECTOR, dealdim>& /*edc*/,
834  dealii::FullMatrix<double> &/*local_entry_matrix*/,
835  double /*scale*/);
836  /******************************************************/
837 
853  virtual void
854  ControlBoundaryEquation(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
855  dealii::Vector<double> &/*local_vector*/,
856  double /*scale*/);
857 
858  /******************************************************/
869  virtual void
870  ControlBoundaryMatrix(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
871  dealii::FullMatrix<double> &/*local_entry_matrix*/,
872  double /*scale*/);
873  /******************************************************/
874 
890  virtual void
891  StrongElementResidual_Control(const EDC<DH, VECTOR, dealdim>& /*edc*/,
892  const EDC<DH, VECTOR, dealdim>& /*edc_weight*/,
893  double& /*ret*/,
894  double /*scale*/);
895  /******************************************************/
913  virtual void
914  StrongFaceResidual_Control(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
915  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
916  double& /*ret*/,
917  double /*scale*/);
918  /******************************************************/
936  virtual void
937  StrongBoundaryResidual_Control(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
938  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
939  double& /*ret*/,
940  double /*scale*/);
941 
942  /******************************************************/
943  // Functions for Face Integrals
944 
953  virtual void
954  FaceEquation(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
955  dealii::Vector<double> &/*local_vector*/,
956  double /*scale*/,
957  double /*scale_ico*/);
958  /******************************************************/
959 
960  virtual void
961  StrongFaceResidual(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
962  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
963  double& /*ret*/,
964  double /*scale*/);
965 
966  /******************************************************/
967 
968  virtual void
969  FaceEquation_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
970  dealii::Vector<double> &/*local_vector*/,
971  double /*scale*/,
972  double /*scale_ico*/);
973 
974  /******************************************************/
975 
976  virtual void
977  StrongFaceResidual_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
978  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
979  double& /*ret*/,
980  double /*scale*/);
981 
982  /******************************************************/
983 
984  virtual void
985  FaceEquation_UT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
986  dealii::Vector<double> &/*local_vector*/,
987  double /*scale*/,
988  double /*scale_ico*/);
989 
990  /******************************************************/
991 
992  virtual void
993  FaceEquation_UTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
994  dealii::Vector<double> &/*local_vector*/,
995  double /*scale*/,
996  double /*scale_ico*/);
997 
998  /******************************************************/
999 
1000  virtual void
1001  FaceEquation_Q(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1002  dealii::Vector<double> &/*local_vector*/,
1003  double /*scale*/,
1004  double /*scale_ico*/);
1005 
1006  /******************************************************/
1007 
1008  virtual void
1009  FaceEquation_QT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1010  dealii::Vector<double> &/*local_vector*/,
1011  double /*scale*/,
1012  double /*scale_ico*/);
1013 
1014  /******************************************************/
1015 
1016  virtual void
1017  FaceEquation_QTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1018  dealii::Vector<double> &/*local_vector*/,
1019  double /*scale*/,
1020  double /*scale_ico*/);
1021 
1022  /******************************************************/
1023 
1024  virtual void
1025  FaceEquation_UU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1026  dealii::Vector<double> &/*local_vector*/,
1027  double /*scale*/,
1028  double /*scale_ico*/);
1029 
1030  /******************************************************/
1031 
1032  virtual void
1033  FaceEquation_QU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1034  dealii::Vector<double> &/*local_vector*/,
1035  double /*scale*/,
1036  double /*scale_ico*/);
1037 
1038  /******************************************************/
1039 
1040  virtual void
1041  FaceEquation_UQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1042  dealii::Vector<double> &/*local_vector*/,
1043  double /*scale*/,
1044  double /*scale_ico*/);
1045 
1046  /******************************************************/
1047 
1048  virtual void
1049  FaceEquation_QQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1050  dealii::Vector<double> &/*local_vector*/,
1051  double /*scale*/,
1052  double /*scale_ico*/);
1053 
1054  /******************************************************/
1055 
1056  virtual void
1057  FaceRightHandSide(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1058  dealii::Vector<double> &/*local_vector*/,
1059  double /*scale*/);
1060 
1061  /******************************************************/
1062 
1066  virtual void
1067  FaceMatrix(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1068  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1069  double /*scale*/,
1070  double /*scale_ico*/);
1071 
1072  /******************************************************/
1073 
1074  virtual void
1075  FaceMatrix_T(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1076  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1077  double /*scale*/,
1078  double /*scale_ico*/);
1079 
1080  /******************************************************/
1081  //Functions for Interface Integrals
1082 
1083  virtual void
1084  InterfaceMatrix(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1085  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1086  double /*scale*/,
1087  double /*scale_ico*/);
1088 
1089  /******************************************************/
1090  //Functions for Interface Integrals
1091  virtual void
1092  InterfaceMatrix_T(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1093  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1094  double /*scale*/,
1095  double /*scale_ico*/);
1096 
1097  /******************************************************/
1098 
1099  virtual void
1100  InterfaceEquation(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1101  dealii::Vector<double> &/*local_vector*/,
1102  double /*scale*/,
1103  double /*scale_ico*/);
1104 
1105  /******************************************************/
1106 
1107  virtual void
1108  InterfaceEquation_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1109  dealii::Vector<double> &/*local_vector*/,
1110  double /*scale*/,
1111  double /*scale_ico*/);
1112 
1113  /******************************************************/
1114  // Functions for Boundary Integrals
1115 
1116  virtual void
1117  BoundaryEquation(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1118  dealii::Vector<double> &/*local_vector*/,
1119  double /*scale*/,
1120  double /*scale_ico*/);
1121 
1122  /******************************************************/
1123 
1124  virtual void
1125  StrongBoundaryResidual(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1126  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
1127  double& /*ret*/,
1128  double /*scale*/);
1129 
1130  /******************************************************/
1131 
1132  virtual void
1133  BoundaryEquation_U(const FDC<DH, VECTOR, dealdim>&/*fdc*/,
1134  dealii::Vector<double> &/*local_vector*/,
1135  double /*scale*/,
1136  double /*scale_ico*/);
1137 
1138  /******************************************************/
1139 
1140  virtual void
1141  StrongBoundaryResidual_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1142  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
1143  double& /*ret*/,
1144  double /*scale*/);
1145 
1146  /******************************************************/
1147 
1148  virtual void
1149  BoundaryEquation_UT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1150  dealii::Vector<double> &/*local_vector*/,
1151  double /*scale*/,
1152  double /*scale_ico*/);
1153 
1154  /******************************************************/
1155 
1156  virtual void
1157  BoundaryEquation_UTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1158  dealii::Vector<double> &/*local_vector*/,
1159  double /*scale*/,
1160  double /*scale_ico*/);
1161 
1162  /******************************************************/
1163 
1164  virtual void
1165  BoundaryEquation_Q(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1166  dealii::Vector<double> &/*local_vector*/,
1167  double /*scale*/,
1168  double /*scale_ico*/);
1169 
1170  /******************************************************/
1171 
1172  virtual void
1173  BoundaryEquation_QT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1174  dealii::Vector<double> &/*local_vector*/,
1175  double /*scale*/,
1176  double /*scale_ico*/);
1177 
1178  /******************************************************/
1179 
1180  virtual void
1181  BoundaryEquation_QTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1182  dealii::Vector<double> &/*local_vector*/,
1183  double /*scale*/,
1184  double /*scale_ico*/);
1185 
1186  /******************************************************/
1187 
1188  virtual void
1189  BoundaryEquation_UU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1190  dealii::Vector<double> &/*local_vector*/,
1191  double /*scale*/,
1192  double /*scale_ico*/);
1193 
1194  /******************************************************/
1195 
1196  virtual void
1197  BoundaryEquation_QU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1198  dealii::Vector<double> &/*local_vector*/,
1199  double /*scale*/,
1200  double /*scale_ico*/);
1201 
1202  /******************************************************/
1203 
1204  virtual void
1205  BoundaryEquation_UQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1206  dealii::Vector<double> &/*local_vector*/,
1207  double /*scale*/,
1208  double /*scale_ico*/);
1209 
1210  /******************************************************/
1211 
1212  virtual void
1213  BoundaryEquation_QQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1214  dealii::Vector<double> &/*local_vector*/,
1215  double /*scale*/,
1216  double /*scale_ico*/);
1217 
1218  /******************************************************/
1219 
1220  virtual void
1221  BoundaryRightHandSide(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1222  dealii::Vector<double> &/*local_vector*/,
1223  double /*scale*/);
1224 
1225  /******************************************************/
1226 
1227  virtual void
1228  BoundaryMatrix(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1229  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1230  double /*scale*/,
1231  double /*scale_ico*/);
1232 
1233  /******************************************************/
1234 
1235  virtual void
1236  BoundaryMatrix_T(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1237  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1238  double /*scale*/,
1239  double /*scale_ico*/);
1240 
1241  /******************************************************/
1242  /******************************************************/
1256  virtual void
1257  Init_ElementEquation(const EDC<DH, VECTOR, dealdim>& edc,
1258  dealii::Vector<double> &local_vector,
1259  double scale,
1260  double /*scale_ico*/)
1261  {
1262  const DOpEWrapper::FEValues<dealdim> & state_fe_values =
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);
1270 
1271  dealii::Vector<double> f_values(
1272  dealii::Vector<double>(this->GetStateNComponents()));
1273 
1274  for (unsigned int q_point = 0; q_point < n_q_points; q_point++)
1275  {
1276  for (unsigned int i = 0; i < n_dofs_per_element; i++)
1277  {
1278  for (unsigned int comp = 0; comp < this->GetStateNComponents();
1279  comp++)
1280  {
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);
1285  }
1286  }
1287  } //endfor q_point
1288  }
1289 
1290  virtual void
1291  Init_ElementRhs_Q(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1292  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1293  {
1294 
1295  }
1296  virtual void
1297  Init_ElementRhs_QT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1298  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1299  {
1300 
1301  }
1302  virtual void
1303  Init_ElementRhs_QTT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1304  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1305  {
1306 
1307  }
1308  virtual void
1309  Init_ElementRhs_QQ(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1310  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1311  {
1312 
1313  }
1314 
1315  virtual void
1316  Init_ElementRhs(const dealii::Function<dealdim>* init_values,
1317  const EDC<DH, VECTOR, dealdim>& edc,
1318  dealii::Vector<double> &local_vector,
1319  double scale)
1320  {
1321  const DOpEWrapper::FEValues<dealdim> & state_fe_values =
1322  edc.GetFEValuesState();
1323  unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1324  unsigned int n_q_points = edc.GetNQPoints();
1325 
1326  dealii::Vector<double> f_values(
1327  dealii::Vector<double>(this->GetStateNComponents()));
1328 
1329  for (unsigned int q_point = 0; q_point < n_q_points; q_point++)
1330  {
1331  init_values->vector_value(state_fe_values.quadrature_point(q_point),
1332  f_values);
1333 
1334  for (unsigned int i = 0; i < n_dofs_per_element; i++)
1335  {
1336  for (unsigned int comp = 0; comp < this->GetStateNComponents();
1337  comp++)
1338  {
1339  local_vector(i) += scale
1340  * (f_values(comp)
1341  * state_fe_values.shape_value_component(i, q_point,
1342  comp)) * state_fe_values.JxW(q_point);
1343  }
1344  }
1345  }
1346  }
1347 
1348  virtual void
1349  Init_ElementMatrix(const EDC<DH, VECTOR, dealdim>& edc,
1350  dealii::FullMatrix<double> &local_entry_matrix,
1351  double scale,
1352  double /*scale_ico*/)
1353  {
1354  const DOpEWrapper::FEValues<dealdim> & state_fe_values =
1355  edc.GetFEValuesState();
1356  unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1357  unsigned int n_q_points = edc.GetNQPoints();
1358 
1359  for (unsigned int q_point = 0; q_point < n_q_points; q_point++)
1360  {
1361  for (unsigned int i = 0; i < n_dofs_per_element; i++)
1362  {
1363  for (unsigned int j = 0; j < n_dofs_per_element; j++)
1364  {
1365  for (unsigned int comp = 0; comp < this->GetStateNComponents();
1366  comp++)
1367  {
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);
1372  }
1373  }
1374  }
1375  }
1376  }
1377  /*************************************************************/
1383  virtual dealii::UpdateFlags
1384  GetUpdateFlags() const;
1390  virtual dealii::UpdateFlags
1391  GetFaceUpdateFlags() const;
1392 
1399  virtual bool
1400  HasFaces() const;
1401  virtual bool
1402  HasInterfaces() const;
1403 
1404  /******************************************************/
1405 
1406  void
1407  SetProblemType(std::string type);
1408 
1409  /******************************************************/
1410 
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;
1423 
1424  /******************************************************/
1425 
1426  virtual void
1427  SetTime(double /*t*/) const
1428  {
1429  }
1430 
1431  /******************************************************/
1432 
1433  unsigned int
1434  GetStateNComponents() const;
1435 
1436  /******************************************************/
1443  boost::function1<void, double&> ResidualModifier;
1444  boost::function1<void, dealii::Vector<double>&> VectorResidualModifier;
1445 
1459  template<typename ELEMENTITERATOR>
1460  bool
1461  AtInterface(ELEMENTITERATOR& element, unsigned int face) const
1462  {
1463  if (element[0]->neighbor_index(face) != -1)
1464  if (element[0]->material_id()
1465  != element[0]->neighbor(face)->material_id())
1466  return true;
1467  return false;
1468  }
1469 
1470  protected:
1471  std::string problem_type_;
1472 
1473  private:
1474  };
1475 }
1476 
1477 #endif
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