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