1 #pragma once // Ensure that file is included only once in a single compilation.
19 template <
unsigned int space_dimT,
typename param_
float_t =
double>
20 struct DiffusionParametersDefault
29 static constexpr std::array<unsigned int, 0U>
neumann_nodes{};
34 const param_float_t = 0.)
42 const param_float_t = 0.)
50 const param_float_t = 0.)
58 const param_float_t = 0.)
66 const param_float_t = 0.)
103 template <
unsigned int hyEdge_dimT,
104 unsigned int poly_deg,
105 unsigned int quad_deg,
106 template <
unsigned int,
typename>
typename parametersT = DiffusionParametersDefault,
107 typename lSol_float_t =
double>
140 std::array<lSol_float_t, 1U> summed_error;
141 summed_error.fill(0.);
149 for (
unsigned int k = 0; k < summed_error.size(); ++k)
150 summed_error[k] += new_error[k];
158 for (
unsigned int k = 0; k < summed_error.size(); ++k)
159 summed_error[k] = std::sqrt(summed_error[k]);
171 static constexpr
unsigned int hyEdge_dim() {
return hyEdge_dimT; }
225 template <
typename parameters>
228 return std::find(parameters::dirichlet_nodes.begin(), parameters::dirichlet_nodes.end(),
229 node_type) != parameters::dirichlet_nodes.end();
262 template <
typename hyEdgeT>
264 assemble_loc_matrix(
const lSol_float_t tau, hyEdgeT& hyper_edge,
const lSol_float_t time)
const;
285 template <
typename hyEdgeT,
typename SmallMatT>
287 const SmallMatT& lambda_values,
288 hyEdgeT& hyper_edge)
const;
307 template <
typename hyEdgeT>
310 const lSol_float_t time)
const;
330 template <
typename hyEdgeT,
typename SmallVecT>
332 hyEdgeT& hyper_edge)
const;
344 template <
typename hyEdgeT,
typename SmallMatT>
346 const unsigned int solution_type,
348 const lSol_float_t time)
const
353 if (solution_type == 0)
355 else if (solution_type == 1)
359 hy_assert(0 == 1,
"This has not been implemented!");
365 <<
"This can happen if quadrature is too inaccurate!");
381 template <
typename hyEdgeT,
typename SmallMatT>
384 const lSol_float_t time)
const
394 <<
"This can happen if quadrature is too inaccurate!");
403 template <
typename hyEdgeT,
typename SmallMatT,
typename SmallVecT>
405 const SmallMatT& lambda_values,
406 const SmallVecT& coeffs,
408 const lSol_float_t
UNUSED(delta_time),
409 const lSol_float_t time)
const
422 <<
"This can happen if quadrature is too inaccurate!");
437 template <
typename hyEdgeT>
440 hyEdgeT& hyper_edge)
const;
452 template <
typename hyEdgeT>
455 hyEdgeT& hyper_edge)
const;
489 template <
typename hyEdgeT,
typename SmallMatInT,
typename SmallMatOutT>
491 SmallMatOutT& lambda_values_out,
493 const lSol_float_t time = 0.)
const
495 hy_assert(lambda_values_in.size() == lambda_values_out.size() &&
496 lambda_values_in.size() == 2 * hyEdge_dimT,
497 "Both matrices must be of same size which corresponds to the number of faces!");
498 for (
unsigned int i = 0; i < lambda_values_in.size(); ++i)
500 lambda_values_in[i].size() == lambda_values_out[i].size() &&
502 "Both matrices must be of same size which corresponds to the number of dofs per face!");
512 for (
unsigned int i = 0; i < lambda_values_out.size(); ++i)
513 if (is_dirichlet<parameters>(hyper_edge.node_descriptor[i]))
514 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
515 lambda_values_out[i][j] = 0.;
517 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
518 lambda_values_out[i][j] =
519 duals(i, j) +
tau_ * primals(i, j) -
520 tau_ * lambda_values_in[i][j] * hyper_edge.geometry.face_area(i);
522 return lambda_values_out;
527 template <
typename hyEdgeT>
528 std::array<unsigned int, 2 * hyEdge_dimT>
node_types(hyEdgeT& hyper_edge)
const
532 std::array<unsigned int, 2 * hyEdge_dimT> result;
535 for (
unsigned int i = 0; i < 2 * hyEdge_dimT; ++i)
536 if (is_dirichlet<parameters>(hyper_edge.node_descriptor[i]))
553 template <
typename hyEdgeT,
typename SmallMatInT,
typename SmallMatOutT>
555 SmallMatOutT& lambda_values_out,
557 const lSol_float_t time = 0.)
const
559 hy_assert(lambda_values_in.size() == lambda_values_out.size() &&
560 lambda_values_in.size() == 2 * hyEdge_dimT,
561 "Both matrices must be of same size which corresponds to the number of faces!");
562 for (
unsigned int i = 0; i < lambda_values_in.size(); ++i)
564 lambda_values_in[i].size() == lambda_values_out[i].size() &&
566 "Both matrices must be of same size which corresponds to the number of dofs per face!");
576 for (
unsigned int i = 0; i < lambda_values_out.size(); ++i)
578 if (is_dirichlet<parameters>(hyper_edge.node_descriptor[i]))
579 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
580 lambda_values_out[i][j] = 0.;
582 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
583 lambda_values_out[i][j] =
584 duals(i, j) +
tau_ * primals(i, j) -
585 tau_ * lambda_values_in[i][j] * hyper_edge.geometry.face_area(i);
588 return lambda_values_out;
602 template <
typename hyEdgeT,
typename SmallMatT>
605 const lSol_float_t time = 0.)
const
607 hy_assert(lambda_values.size() == 2 * hyEdge_dimT,
"Matrix must have appropriate size!");
608 for (
unsigned int i = 0; i < lambda_values.size(); ++i)
610 "Matrix must have appropriate size!");
613 for (
unsigned int i = 0; i < lambda_values.size(); ++i)
615 if (is_dirichlet<parameters>(hyper_edge.node_descriptor[i]))
616 for (
unsigned int j = 0; j < lambda_values[i].size(); ++j)
617 lambda_values[i][j] = 0.;
619 for (
unsigned int j = 0; j < lambda_values[i].size(); ++j)
620 lambda_values[i][j] =
622 parameters::initial>(
623 j, i, hyper_edge.geometry, time);
626 return lambda_values;
718 template <
typename hyEdgeT,
typename SmallMatInT,
typename SmallMatOutT>
720 SmallMatOutT& lambda_values_out,
722 const lSol_float_t time = 0.)
const
724 hy_assert(lambda_values_in.size() == lambda_values_out.size() &&
725 lambda_values_in.size() == 2 * hyEdge_dimT,
726 "Both matrices must be of same size which corresponds to the number of faces!");
727 for (
unsigned int i = 0; i < lambda_values_in.size(); ++i)
729 lambda_values_in[i].size() == lambda_values_out[i].size() &&
731 "Both matrices must be of same size which corresponds to the number of dofs per face!");
742 for (
unsigned int i = 0; i < lambda_values_out.size(); ++i)
744 if (is_dirichlet<parameters>(hyper_edge.node_descriptor[i]))
745 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
746 lambda_values_out[i][j] = 0.;
748 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
749 lambda_values_out[i][j] = duals(i, j) +
tau_ * primals(i, j);
752 return lambda_values_out;
768 template <
typename hyEdgeT,
typename SmallMatInT,
typename SmallMatOutT>
770 SmallMatOutT& lambda_values_out,
772 const lSol_float_t time = 0.)
const
774 hy_assert(lambda_values_in.size() == lambda_values_out.size() &&
775 lambda_values_in.size() == 2 * hyEdge_dimT,
776 "Both matrices must be of same size which corresponds to the number of faces!");
777 for (
unsigned int i = 0; i < lambda_values_in.size(); ++i)
779 lambda_values_in[i].size() == lambda_values_out[i].size() &&
781 "Both matrices must be of same size which corresponds to the number of dofs per face!");
792 for (
unsigned int i = 0; i < lambda_values_out.size(); ++i)
794 if (is_dirichlet<parameters>(hyper_edge.node_descriptor[i]))
795 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
796 lambda_values_out[i][j] = 0.;
798 for (
unsigned int j = 0; j < lambda_values_out[i].size(); ++j)
799 lambda_values_out[i][j] = duals(i, j) +
tau_ * primals(i, j);
802 return lambda_values_out;
814 template <
typename hyEdgeT,
typename SmallMatT>
815 std::array<lSol_float_t, 1U>
errors(
const SmallMatT& lambda_values,
817 const lSol_float_t time = 0.)
const
819 hy_assert(lambda_values.size() == 2 * hyEdge_dimT,
"Matrix must have appropriate size!");
820 for (
unsigned int i = 0; i < lambda_values.size(); ++i)
822 "Matrix must have appropriate size!");
827 std::array<lSol_float_t, n_shape_fct_> coeffs;
828 for (
unsigned int i = 0; i < coeffs.size(); ++i)
829 coeffs[i] = coefficients[i + hyEdge_dimT *
n_shape_fct_];
831 return std::array<lSol_float_t, 1U>({integrator::template integrate_vol_diffsquare_discana<
850 template <
typename hyEdgeT,
typename SmallMatT>
852 const SmallMatT& lambda_values_old,
854 const lSol_float_t delta_t,
855 const lSol_float_t time)
const
857 hy_assert(lambda_values_new.size() == lambda_values_old.size() &&
858 lambda_values_new.size() == 2 * hyEdge_dimT,
859 "Both matrices must be of same size which corresponds to the number of faces!");
860 for (
unsigned int i = 0; i < lambda_values_new.size(); ++i)
862 lambda_values_new[i].size() == lambda_values_old[i].size() &&
864 "Both matrices must be of same size which corresponds to the number of dofs per face!");
872 for (
unsigned int i = 0; i < coeffs_old.
size(); ++i)
873 coeffs_old[i] = (coeffs_old[i] - coeffs_new[i]) / delta_t;
877 std::array<lSol_float_t, n_shape_fct_> coeffs;
878 for (
unsigned int i = 0; i < coeffs.size(); ++i)
879 coeffs[i] = coefficients[i + hyEdge_dimT *
n_shape_fct_];
880 return integrator::template integrate_vol_diffsquare_discana<decltype(
hyEdgeT::geometry),
881 parameters::analytic_result>(
882 coeffs, hy_edge.geometry, time);
897 template <
typename abscissa_float_t,
898 std::size_t abscissas_sizeT,
901 std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(abscissas_sizeT)>,
system_dim>
902 bulk_values(
const std::array<abscissa_float_t, abscissas_sizeT>& abscissas,
903 const input_array_t& lambda_values,
905 const lSol_float_t time = 0.)
const;
924 template <
unsigned int hyEdge_dimT,
925 unsigned int poly_deg,
926 unsigned int quad_deg,
927 template <
unsigned int,
typename>
928 typename parametersT,
929 typename lSol_float_t>
930 template <
typename hyEdgeT>
935 const lSol_float_t tau,
937 const lSol_float_t time)
const
941 lSol_float_t vol_integral, face_integral, helper;
944 for (
unsigned int i = 0; i < n_shape_fct_; ++i)
946 for (
unsigned int j = 0; j < n_shape_fct_; ++j)
949 vol_integral = integrator::template integrate_vol_phiphifunc<
951 parameters::inverse_diffusion_coeff,
Point<hyEdge_dimT, lSol_float_t> >(
956 integrator::template integrate_vol_nablaphiphi<SmallVec<hyEdge_dimT, lSol_float_t>,
958 i, j, hyper_edge.geometry);
962 for (
unsigned int face = 0; face < 2 * hyEdge_dimT; ++face)
964 helper = integrator::template integrate_bdr_phiphi<decltype(hyEdgeT::geometry)>(
965 i, j, face, hyper_edge.geometry);
966 face_integral += helper;
967 normal_int_vec += helper * hyper_edge.geometry.local_normal(face);
970 local_mat(hyEdge_dimT * n_shape_fct_ + i, hyEdge_dimT * n_shape_fct_ + j) +=
972 for (
unsigned int dim = 0; dim < hyEdge_dimT; ++dim)
974 local_mat(dim * n_shape_fct_ + i, dim * n_shape_fct_ + j) += vol_integral;
975 local_mat(hyEdge_dimT * n_shape_fct_ + i, dim * n_shape_fct_ + j) -= grad_int_vec[dim];
976 local_mat(dim * n_shape_fct_ + i, hyEdge_dimT * n_shape_fct_ + j) -= grad_int_vec[dim];
977 local_mat(hyEdge_dimT * n_shape_fct_ + i, dim * n_shape_fct_ + j) += normal_int_vec[dim];
989 template <
unsigned int hyEdge_dimT,
990 unsigned int poly_deg,
991 unsigned int quad_deg,
992 template <
unsigned int,
typename>
993 typename parametersT,
994 typename lSol_float_t>
995 template <
typename hyEdgeT,
typename SmallMatT>
999 const SmallMatT& lambda_values,
1000 hyEdgeT& hyper_edge)
const
1002 static_assert(std::is_same<typename SmallMatT::value_type::value_type, lSol_float_t>::value,
1003 "Lambda values should have same floating point arithmetics as local solver!");
1004 hy_assert(lambda_values.size() == 2 * hyEdge_dimT,
1005 "The size of the lambda values should be twice the dimension of a hyperedge.");
1006 for (
unsigned int i = 0; i < 2 * hyEdge_dimT; ++i)
1007 hy_assert(lambda_values[i].size() == n_shape_bdr_,
1008 "The size of lambda should be the amount of ansatz functions at boundary.");
1011 lSol_float_t integral;
1013 for (
unsigned int i = 0; i < n_shape_fct_; ++i)
1014 for (
unsigned int j = 0; j < n_shape_bdr_; ++j)
1015 for (
unsigned int face = 0; face < 2 * hyEdge_dimT; ++face)
1017 integral = integrator::template integrate_bdr_phipsi<decltype(hyEdgeT::geometry)>(
1018 i, j, face, hyper_edge.geometry);
1019 right_hand_side[hyEdge_dimT * n_shape_fct_ + i] += tau_ * lambda_values[face][j] * integral;
1020 for (
unsigned int dim = 0; dim < hyEdge_dimT; ++dim)
1021 right_hand_side[dim * n_shape_fct_ + i] -=
1022 hyper_edge.geometry.local_normal(face).operator[](dim) * lambda_values[face][j] *
1026 return right_hand_side;
1033 template <
unsigned int hyEdge_dimT,
1034 unsigned int poly_deg,
1035 unsigned int quad_deg,
1036 template <
unsigned int,
typename>
1037 typename parametersT,
1038 typename lSol_float_t>
1039 template <
typename hyEdgeT>
1043 hyEdgeT& hyper_edge,
1044 const lSol_float_t time)
const
1048 lSol_float_t integral;
1049 for (
unsigned int i = 0; i < n_shape_fct_; ++i)
1051 right_hand_side[hyEdge_dimT * n_shape_fct_ + i] = integrator::template integrate_vol_phifunc<
1054 for (
unsigned int face = 0; face < 2 * hyEdge_dimT; ++face)
1056 if (!is_dirichlet<parameters>(hyper_edge.node_descriptor[face]))
1058 integral = integrator::template integrate_bdr_phifunc<
1060 parameters::dirichlet_value,
Point<hyEdge_dimT, lSol_float_t> >(i, face,
1062 right_hand_side[hyEdge_dimT * n_shape_fct_ + i] += tau_ * integral;
1063 for (
unsigned int dim = 0; dim < hyEdge_dimT; ++dim)
1064 right_hand_side[dim * n_shape_fct_ + i] -=
1065 hyper_edge.geometry.local_normal(face).operator[](dim) * integral;
1069 return right_hand_side;
1076 template <
unsigned int hyEdge_dimT,
1077 unsigned int poly_deg,
1078 unsigned int quad_deg,
1079 template <
unsigned int,
typename>
1080 typename parametersT,
1081 typename lSol_float_t>
1082 template <
typename hyEdgeT,
typename SmallVecT>
1086 const SmallVecT& coeffs,
1087 hyEdgeT& hyper_edge)
const
1090 for (
unsigned int i = 0; i < n_shape_fct_; ++i)
1091 for (
unsigned int j = 0; j < n_shape_fct_; ++j)
1092 right_hand_side[hyEdge_dimT * n_shape_fct_ + i] +=
1093 coeffs[hyEdge_dimT * n_shape_fct_ + j] *
1094 integrator::template integrate_vol_phiphi(i, j, hyper_edge.geometry);
1096 return right_hand_side;
1103 template <
unsigned int hyEdge_dimT,
1104 unsigned int poly_deg,
1105 unsigned int quad_deg,
1106 template <
unsigned int,
typename>
1107 typename parametersT,
1108 typename lSol_float_t>
1109 template <
typename hyEdgeT>
1115 hyEdgeT& hyper_edge)
const
1119 for (
unsigned int i = 0; i < n_shape_fct_; ++i)
1120 for (
unsigned int j = 0; j < n_shape_bdr_; ++j)
1121 for (
unsigned int face = 0; face < 2 * hyEdge_dimT; ++face)
1122 bdr_values(face, j) +=
1123 coeffs[hyEdge_dimT * n_shape_fct_ + i] *
1124 integrator::template integrate_bdr_phipsi<decltype(hyEdgeT::geometry)>(
1125 i, j, face, hyper_edge.geometry);
1134 template <
unsigned int hyEdge_dimT,
1135 unsigned int poly_deg,
1136 unsigned int quad_deg,
1137 template <
unsigned int,
typename>
1138 typename parametersT,
1139 typename lSol_float_t>
1140 template <
typename hyEdgeT>
1146 hyEdgeT& hyper_edge)
const
1149 lSol_float_t integral;
1151 for (
unsigned int i = 0; i < n_shape_fct_; ++i)
1152 for (
unsigned int j = 0; j < n_shape_bdr_; ++j)
1153 for (
unsigned int face = 0; face < 2 * hyEdge_dimT; ++face)
1155 integral = integrator::template integrate_bdr_phipsi<decltype(hyEdgeT::geometry)>(
1156 i, j, face, hyper_edge.geometry);
1157 for (
unsigned int dim = 0; dim < hyEdge_dimT; ++dim)
1158 bdr_values(face, j) += hyper_edge.geometry.local_normal(face).operator[](dim) * integral *
1159 coeffs[dim * n_shape_fct_ + i];
1169 template <
unsigned int hyEdge_dimT,
1170 unsigned int poly_deg,
1171 unsigned int quad_deg,
1172 template <
unsigned int,
typename>
1173 typename parametersT,
1174 typename lSol_float_t>
1175 template <
typename abscissa_float_t,
1176 std::size_t abscissas_sizeT,
1177 class input_array_t,
1179 std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(abscissas_sizeT)>,
1182 const std::array<abscissa_float_t, abscissas_sizeT>& abscissas,
1183 const input_array_t& lambda_values,
1184 hyEdgeT& hyper_edge,
1185 const lSol_float_t time)
const
1188 solve_local_problem(lambda_values, 1U, hyper_edge, time);
1192 std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(abscissas_sizeT)>,
1196 for (
unsigned int d = 0; d < system_dim; ++d)
1198 for (
unsigned int i = 0; i < coeffs.
size(); ++i)
1199 coeffs[i] = coefficients[d * n_shape_fct_ + i];
1200 for (
unsigned int pt = 0; pt < Hypercube<hyEdge_dimT>::pow(abscissas_sizeT); ++pt)
1201 point_vals[d][pt] = integrator::shape_fun_t::template lin_comb_fct_val<float>(