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 <fe/fe_system.h>
31 #include <fe/fe_values.h>
32 #include <fe/mapping.h>
33 #include <lac/full_matrix.h>
34 #include <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  /******************************************************/
836 
852  virtual void
853  StrongElementResidual_Control(const EDC<DH, VECTOR, dealdim>& /*edc*/,
854  const EDC<DH, VECTOR, dealdim>& /*edc_weight*/,
855  double& /*ret*/,
856  double /*scale*/);
857  /******************************************************/
875  virtual void
876  StrongFaceResidual_Control(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
877  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
878  double& /*ret*/,
879  double /*scale*/);
880  /******************************************************/
898  virtual void
899  StrongBoundaryResidual_Control(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
900  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
901  double& /*ret*/,
902  double /*scale*/);
903 
904  /******************************************************/
905  // Functions for Face Integrals
906 
915  virtual void
916  FaceEquation(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
917  dealii::Vector<double> &/*local_vector*/,
918  double /*scale*/,
919  double /*scale_ico*/);
920  /******************************************************/
921 
922  virtual void
923  StrongFaceResidual(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
924  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
925  double& /*ret*/,
926  double /*scale*/);
927 
928  /******************************************************/
929 
930  virtual void
931  FaceEquation_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
932  dealii::Vector<double> &/*local_vector*/,
933  double /*scale*/,
934  double /*scale_ico*/);
935 
936  /******************************************************/
937 
938  virtual void
939  StrongFaceResidual_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
940  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
941  double& /*ret*/,
942  double /*scale*/);
943 
944  /******************************************************/
945 
946  virtual void
947  FaceEquation_UT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
948  dealii::Vector<double> &/*local_vector*/,
949  double /*scale*/,
950  double /*scale_ico*/);
951 
952  /******************************************************/
953 
954  virtual void
955  FaceEquation_UTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
956  dealii::Vector<double> &/*local_vector*/,
957  double /*scale*/,
958  double /*scale_ico*/);
959 
960  /******************************************************/
961 
962  virtual void
963  FaceEquation_Q(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
964  dealii::Vector<double> &/*local_vector*/,
965  double /*scale*/,
966  double /*scale_ico*/);
967 
968  /******************************************************/
969 
970  virtual void
971  FaceEquation_QT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
972  dealii::Vector<double> &/*local_vector*/,
973  double /*scale*/,
974  double /*scale_ico*/);
975 
976  /******************************************************/
977 
978  virtual void
979  FaceEquation_QTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
980  dealii::Vector<double> &/*local_vector*/,
981  double /*scale*/,
982  double /*scale_ico*/);
983 
984  /******************************************************/
985 
986  virtual void
987  FaceEquation_UU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
988  dealii::Vector<double> &/*local_vector*/,
989  double /*scale*/,
990  double /*scale_ico*/);
991 
992  /******************************************************/
993 
994  virtual void
995  FaceEquation_QU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
996  dealii::Vector<double> &/*local_vector*/,
997  double /*scale*/,
998  double /*scale_ico*/);
999 
1000  /******************************************************/
1001 
1002  virtual void
1003  FaceEquation_UQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1004  dealii::Vector<double> &/*local_vector*/,
1005  double /*scale*/,
1006  double /*scale_ico*/);
1007 
1008  /******************************************************/
1009 
1010  virtual void
1011  FaceEquation_QQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1012  dealii::Vector<double> &/*local_vector*/,
1013  double /*scale*/,
1014  double /*scale_ico*/);
1015 
1016  /******************************************************/
1017 
1018  virtual void
1019  FaceRightHandSide(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1020  dealii::Vector<double> &/*local_vector*/,
1021  double /*scale*/);
1022 
1023  /******************************************************/
1024 
1028  virtual void
1029  FaceMatrix(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1030  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1031  double /*scale*/,
1032  double /*scale_ico*/);
1033 
1034  /******************************************************/
1035 
1036  virtual void
1037  FaceMatrix_T(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1038  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1039  double /*scale*/,
1040  double /*scale_ico*/);
1041 
1042  /******************************************************/
1043  //Functions for Interface Integrals
1044 
1045  virtual void
1046  InterfaceMatrix(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1047  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1048  double /*scale*/,
1049  double /*scale_ico*/);
1050 
1051  /******************************************************/
1052  //Functions for Interface Integrals
1053  virtual void
1054  InterfaceMatrix_T(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1055  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1056  double /*scale*/,
1057  double /*scale_ico*/);
1058 
1059  /******************************************************/
1060 
1061  virtual void
1062  InterfaceEquation(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1063  dealii::Vector<double> &/*local_vector*/,
1064  double /*scale*/,
1065  double /*scale_ico*/);
1066 
1067  /******************************************************/
1068 
1069  virtual void
1070  InterfaceEquation_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1071  dealii::Vector<double> &/*local_vector*/,
1072  double /*scale*/,
1073  double /*scale_ico*/);
1074 
1075  /******************************************************/
1076  // Functions for Boundary Integrals
1077 
1078  virtual void
1079  BoundaryEquation(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1080  dealii::Vector<double> &/*local_vector*/,
1081  double /*scale*/,
1082  double /*scale_ico*/);
1083 
1084  /******************************************************/
1085 
1086  virtual void
1087  StrongBoundaryResidual(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1088  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
1089  double& /*ret*/,
1090  double /*scale*/);
1091 
1092  /******************************************************/
1093 
1094  virtual void
1095  BoundaryEquation_U(const FDC<DH, VECTOR, dealdim>&/*fdc*/,
1096  dealii::Vector<double> &/*local_vector*/,
1097  double /*scale*/,
1098  double /*scale_ico*/);
1099 
1100  /******************************************************/
1101 
1102  virtual void
1103  StrongBoundaryResidual_U(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1104  const FDC<DH, VECTOR, dealdim>& /*fdc_weight*/,
1105  double& /*ret*/,
1106  double /*scale*/);
1107 
1108  /******************************************************/
1109 
1110  virtual void
1111  BoundaryEquation_UT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1112  dealii::Vector<double> &/*local_vector*/,
1113  double /*scale*/,
1114  double /*scale_ico*/);
1115 
1116  /******************************************************/
1117 
1118  virtual void
1119  BoundaryEquation_UTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1120  dealii::Vector<double> &/*local_vector*/,
1121  double /*scale*/,
1122  double /*scale_ico*/);
1123 
1124  /******************************************************/
1125 
1126  virtual void
1127  BoundaryEquation_Q(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1128  dealii::Vector<double> &/*local_vector*/,
1129  double /*scale*/,
1130  double /*scale_ico*/);
1131 
1132  /******************************************************/
1133 
1134  virtual void
1135  BoundaryEquation_QT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1136  dealii::Vector<double> &/*local_vector*/,
1137  double /*scale*/,
1138  double /*scale_ico*/);
1139 
1140  /******************************************************/
1141 
1142  virtual void
1143  BoundaryEquation_QTT(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1144  dealii::Vector<double> &/*local_vector*/,
1145  double /*scale*/,
1146  double /*scale_ico*/);
1147 
1148  /******************************************************/
1149 
1150  virtual void
1151  BoundaryEquation_UU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1152  dealii::Vector<double> &/*local_vector*/,
1153  double /*scale*/,
1154  double /*scale_ico*/);
1155 
1156  /******************************************************/
1157 
1158  virtual void
1159  BoundaryEquation_QU(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1160  dealii::Vector<double> &/*local_vector*/,
1161  double /*scale*/,
1162  double /*scale_ico*/);
1163 
1164  /******************************************************/
1165 
1166  virtual void
1167  BoundaryEquation_UQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1168  dealii::Vector<double> &/*local_vector*/,
1169  double /*scale*/,
1170  double /*scale_ico*/);
1171 
1172  /******************************************************/
1173 
1174  virtual void
1175  BoundaryEquation_QQ(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1176  dealii::Vector<double> &/*local_vector*/,
1177  double /*scale*/,
1178  double /*scale_ico*/);
1179 
1180  /******************************************************/
1181 
1182  virtual void
1183  BoundaryRightHandSide(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1184  dealii::Vector<double> &/*local_vector*/,
1185  double /*scale*/);
1186 
1187  /******************************************************/
1188 
1189  virtual void
1190  BoundaryMatrix(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1191  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1192  double /*scale*/,
1193  double /*scale_ico*/);
1194 
1195  /******************************************************/
1196 
1197  virtual void
1198  BoundaryMatrix_T(const FDC<DH, VECTOR, dealdim>& /*fdc*/,
1199  dealii::FullMatrix<double> &/*local_entry_matrix*/,
1200  double /*scale*/,
1201  double /*scale_ico*/);
1202 
1203  /******************************************************/
1204  /******************************************************/
1218  virtual void
1219  Init_ElementEquation(const EDC<DH, VECTOR, dealdim>& edc,
1220  dealii::Vector<double> &local_vector,
1221  double scale,
1222  double /*scale_ico*/)
1223  {
1224  const DOpEWrapper::FEValues<dealdim> & state_fe_values =
1225  edc.GetFEValuesState();
1226  unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1227  unsigned int n_q_points = edc.GetNQPoints();
1228  std::vector<dealii::Vector<double> > uvalues;
1229  uvalues.resize(n_q_points,
1230  dealii::Vector<double>(this->GetStateNComponents()));
1231  edc.GetValuesState("last_newton_solution", uvalues);
1232 
1233  dealii::Vector<double> f_values(
1234  dealii::Vector<double>(this->GetStateNComponents()));
1235 
1236  for (unsigned int q_point = 0; q_point < n_q_points; q_point++)
1237  {
1238  for (unsigned int i = 0; i < n_dofs_per_element; i++)
1239  {
1240  for (unsigned int comp = 0; comp < this->GetStateNComponents();
1241  comp++)
1242  {
1243  local_vector(i) += scale
1244  * (state_fe_values.shape_value_component(i, q_point, comp)
1245  * uvalues[q_point](comp))
1246  * state_fe_values.JxW(q_point);
1247  }
1248  }
1249  } //endfor q_point
1250  }
1251 
1252  virtual void
1253  Init_ElementRhs_Q(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1254  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1255  {
1256 
1257  }
1258  virtual void
1259  Init_ElementRhs_QT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1260  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1261  {
1262 
1263  }
1264  virtual void
1265  Init_ElementRhs_QTT(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1266  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1267  {
1268 
1269  }
1270  virtual void
1271  Init_ElementRhs_QQ(const EDC<DH, VECTOR, dealdim>& /*edc*/,
1272  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
1273  {
1274 
1275  }
1276 
1277  virtual void
1278  Init_ElementRhs(const dealii::Function<dealdim>* init_values,
1279  const EDC<DH, VECTOR, dealdim>& edc,
1280  dealii::Vector<double> &local_vector,
1281  double scale)
1282  {
1283  const DOpEWrapper::FEValues<dealdim> & state_fe_values =
1284  edc.GetFEValuesState();
1285  unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1286  unsigned int n_q_points = edc.GetNQPoints();
1287 
1288  dealii::Vector<double> f_values(
1289  dealii::Vector<double>(this->GetStateNComponents()));
1290 
1291  for (unsigned int q_point = 0; q_point < n_q_points; q_point++)
1292  {
1293  init_values->vector_value(state_fe_values.quadrature_point(q_point),
1294  f_values);
1295 
1296  for (unsigned int i = 0; i < n_dofs_per_element; i++)
1297  {
1298  for (unsigned int comp = 0; comp < this->GetStateNComponents();
1299  comp++)
1300  {
1301  local_vector(i) += scale
1302  * (f_values(comp)
1303  * state_fe_values.shape_value_component(i, q_point,
1304  comp)) * state_fe_values.JxW(q_point);
1305  }
1306  }
1307  }
1308  }
1309 
1310  virtual void
1311  Init_ElementMatrix(const EDC<DH, VECTOR, dealdim>& edc,
1312  dealii::FullMatrix<double> &local_entry_matrix,
1313  double scale,
1314  double /*scale_ico*/)
1315  {
1316  const DOpEWrapper::FEValues<dealdim> & state_fe_values =
1317  edc.GetFEValuesState();
1318  unsigned int n_dofs_per_element = edc.GetNDoFsPerElement();
1319  unsigned int n_q_points = edc.GetNQPoints();
1320 
1321  for (unsigned int q_point = 0; q_point < n_q_points; q_point++)
1322  {
1323  for (unsigned int i = 0; i < n_dofs_per_element; i++)
1324  {
1325  for (unsigned int j = 0; j < n_dofs_per_element; j++)
1326  {
1327  for (unsigned int comp = 0; comp < this->GetStateNComponents();
1328  comp++)
1329  {
1330  local_entry_matrix(i, j) += scale
1331  * state_fe_values.shape_value_component(i, q_point, comp)
1332  * state_fe_values.shape_value_component(j, q_point, comp)
1333  * state_fe_values.JxW(q_point);
1334  }
1335  }
1336  }
1337  }
1338  }
1339  /*************************************************************/
1345  virtual dealii::UpdateFlags
1346  GetUpdateFlags() const;
1352  virtual dealii::UpdateFlags
1353  GetFaceUpdateFlags() const;
1354 
1361  virtual bool
1362  HasFaces() const;
1363  virtual bool
1364  HasInterfaces() const;
1365 
1366  /******************************************************/
1367 
1368  void
1369  SetProblemType(std::string type);
1370 
1371  /******************************************************/
1372 
1373  virtual unsigned int
1374  GetControlNBlocks() const;
1375  virtual unsigned int
1376  GetStateNBlocks() const;
1377  virtual std::vector<unsigned int>&
1378  GetControlBlockComponent();
1379  virtual const std::vector<unsigned int>&
1380  GetControlBlockComponent() const;
1381  virtual std::vector<unsigned int>&
1382  GetStateBlockComponent();
1383  virtual const std::vector<unsigned int>&
1384  GetStateBlockComponent() const;
1385 
1386  /******************************************************/
1387 
1388  virtual void
1389  SetTime(double /*t*/) const
1390  {
1391  }
1392 
1393  /******************************************************/
1394 
1395  unsigned int
1396  GetStateNComponents() const;
1397 
1398  /******************************************************/
1405  boost::function1<void, double&> ResidualModifier;
1406  boost::function1<void, dealii::Vector<double>&> VectorResidualModifier;
1407 
1408  protected:
1409  std::string _problem_type;
1410 
1411  private:
1412  };
1413 }
1414 
1415 #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:1278
boost::function1< void, dealii::Vector< double > & > VectorResidualModifier
Definition: pdeinterface.h:1406
virtual void Init_ElementRhs_QT(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1259
virtual void Init_ElementEquation(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale, double)
Definition: pdeinterface.h:1219
Definition: pdeinterface.h:66
boost::function1< void, double & > ResidualModifier
Definition: pdeinterface.h:1405
virtual void Init_ElementRhs_QTT(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1265
virtual void Init_ElementRhs_QQ(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1271
virtual void Init_ElementMatrix(const EDC< DH, VECTOR, dealdim > &edc, dealii::FullMatrix< double > &local_entry_matrix, double scale, double)
Definition: pdeinterface.h:1311
Definition: fevalues_wrapper.h:43
virtual void SetTime(double) const
Definition: pdeinterface.h:1389
virtual void Init_ElementRhs_Q(const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
Definition: pdeinterface.h:1253
std::string _problem_type
Definition: pdeinterface.h:1409