HyperHDG
bernoulli_beams.hxx
Go to the documentation of this file.
1 #pragma once // Ensure that file is included only once in a single compilation.
2 
4 #include <HyperHDG/dense_la.hxx>
5 #include <HyperHDG/hypercube.hxx>
6 //#include <tpp/quadrature/tensorial.hxx>
8 
9 #include <tuple>
10 
11 namespace LocalSolver
12 {
13 /*!*************************************************************************************************
14  * \brief Local solver for the equation that governs the lengthening of an elastic beam.
15  *
16  * \authors Guido Kanschat, Heidelberg University, 2019--2020.
17  * \authors Andreas Rupp, Heidelberg University, 2019--2020.
18  **************************************************************************************************/
19 template <unsigned int hyEdge_dimT,
20  unsigned int space_dim,
21  unsigned int poly_deg,
22  unsigned int quad_deg,
23  typename lSol_float_t = double,
24  typename diffusion_sol_t =
25  DiffusionUniform<hyEdge_dimT, poly_deg, quad_deg, lSol_float_t> >
27 {
28  /*!***********************************************************************************************
29  * \brief Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
30  ************************************************************************************************/
31  HAS_MEMBER_FUNCTION(trace_to_flux, has_trace_to_flux);
32  /*!***********************************************************************************************
33  * \brief Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
34  ************************************************************************************************/
35  HAS_MEMBER_FUNCTION(residual_flux, has_residual_flux);
36  /*!***********************************************************************************************
37  * \brief Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
38  ************************************************************************************************/
39  HAS_MEMBER_FUNCTION(errors, has_errors);
40 
41  public:
42  /*!***********************************************************************************************
43  * \brief Define type of (hyperedge related) data that is stored in HyDataContainer.
44  ************************************************************************************************/
45  struct data_type
46  {
47  };
48  /*!***********************************************************************************************
49  * \brief Define type of node elements, especially with respect to nodal shape functions.
50  ************************************************************************************************/
51  struct node_element
52  {
53  typedef std::tuple<TPP::ShapeFunction<
56  };
57  /*!***********************************************************************************************
58  * \brief Define how errors are evaluated.
59  ************************************************************************************************/
60  struct error_def
61  {
62  /*!*********************************************************************************************
63  * \brief Define the typename returned by function errors.
64  **********************************************************************************************/
65  typedef std::array<lSol_float_t, 1U> error_t;
66  /*!*********************************************************************************************
67  * \brief Define how initial error is generated.
68  **********************************************************************************************/
70  {
71  std::array<lSol_float_t, 1U> summed_error;
72  summed_error.fill(0.);
73  return summed_error;
74  }
75  /*!*********************************************************************************************
76  * \brief Define how local errors should be accumulated.
77  **********************************************************************************************/
78  static error_t sum_error(error_t& summed_error, const error_t& new_error)
79  {
80  for (unsigned int k = 0; k < summed_error.size(); ++k)
81  summed_error[k] += new_error[k];
82  return summed_error;
83  }
84  /*!*********************************************************************************************
85  * \brief Define how global errors should be postprocessed.
86  **********************************************************************************************/
87  static error_t postprocess_error(error_t& summed_error)
88  {
89  for (unsigned int k = 0; k < summed_error.size(); ++k)
90  summed_error[k] = std::sqrt(summed_error[k]);
91  return summed_error;
92  }
93  };
94  /*!***********************************************************************************************
95  * \brief Return template parameter \c hyEdge_dimT.
96  *
97  * \retval hyEdge_dimT Dimension of hypergraph's hyperedges.
98  ************************************************************************************************/
99  static constexpr unsigned int hyEdge_dim() { return hyEdge_dimT; }
100  /*!***********************************************************************************************
101  * \brief Evaluate amount of global degrees of freedom per hypernode.
102  *
103  * This number must be equal to HyperNodeFactory::n_glob_dofs_per_node()() of the HyperNodeFactory
104  * cooperating with this object.
105  *
106  * \retval n_dofs Number of global degrees of freedom per hypernode.
107  ************************************************************************************************/
108  static constexpr unsigned int n_glob_dofs_per_node()
109  {
110  return 2 * space_dim * Hypercube<hyEdge_dimT - 1>::pow(poly_deg + 1);
111  }
112  /*!***********************************************************************************************
113  * \brief Dimension of of the solution evaluated with respect to a hyperedge.
114  ************************************************************************************************/
115  static constexpr unsigned int system_dimension() { return space_dim; }
116  /*!***********************************************************************************************
117  * \brief Dimension of of the solution evaluated with respect to a hypernode.
118  ************************************************************************************************/
119  static constexpr unsigned int node_system_dimension() { return space_dim * 2; }
120 
121  private:
122  /*!***********************************************************************************************
123  * \brief The diffusion solver that solves the locally defined PDE.
124  ************************************************************************************************/
125  const diffusion_sol_t diffusion;
126  /*!***********************************************************************************************
127  * \brief Do the pretprocessing to transfer global to local dofs.
128  ************************************************************************************************/
129  template <class hyEdgeT, typename SmallMatT>
130  inline std::array<std::array<double, diffusion_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT>
131  node_dof_to_edge_dof(const SmallMatT& lambda, hyEdgeT& hyper_edge) const
132  {
133  std::array<std::array<double, diffusion_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT> result;
134  hy_assert(result.size() == 2, "Only implemented in one dimension!");
135  for (unsigned int i = 0; i < result.size(); ++i)
136  {
137  hy_assert(result[i].size() == 1, "Only implemented in one dimension!");
138  result[i].fill(0.);
139  }
140 
141  Point<space_dim, lSol_float_t> normal_vector =
142  (Point<space_dim, lSol_float_t>)hyper_edge.geometry.inner_normal(1);
143 
144  for (unsigned int i = 0; i < 2 * hyEdge_dimT; ++i)
145  for (unsigned int dim = 0; dim < space_dim; ++dim)
146  result[i][0] += normal_vector[dim] * lambda[i][dim];
147 
148  return result;
149  }
150  /*!***********************************************************************************************
151  * \brief Do the postprocessing to transfer local to global dofs.
152  ************************************************************************************************/
153  template <class hyEdgeT, typename SmallMatInT, typename SmallMatOutT>
154  inline SmallMatOutT& edge_dof_to_node_dof(const SmallMatInT& lambda,
155  SmallMatOutT& lambda_values_out,
156  hyEdgeT& hyper_edge) const
157  {
158  hy_assert(diffusion_sol_t::n_glob_dofs_per_node() == 1, "This should be 1!");
159  Point<space_dim, lSol_float_t> normal_vector =
160  (Point<space_dim, lSol_float_t>)hyper_edge.geometry.inner_normal(1);
161 
162  for (unsigned int i = 0; i < 2 * hyEdge_dimT; ++i)
163  for (unsigned int dim = 0; dim < space_dim; ++dim)
164  lambda_values_out[i][dim] += normal_vector[dim] * lambda[i][0];
165 
166  return lambda_values_out;
167  }
168 
169  public:
170  /*!***********************************************************************************************
171  * \brief Class is constructed using a single double indicating the penalty parameter.
172  ************************************************************************************************/
173  typedef lSol_float_t constructor_value_type;
174  /*!***********************************************************************************************
175  * \brief Constructor for local solver.
176  *
177  * \param tau Penalty parameter of HDG scheme.
178  ************************************************************************************************/
180  /*!***********************************************************************************************
181  * \brief Evaluate local contribution to matrix--vector multiplication.
182  *
183  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
184  * \tparam SmallMatInT Data type of \c lambda_values_in.
185  * \tparam SmallMatOutT Data type of \c lambda_values_out.
186  * \param lambda_values_in Local part of vector x.
187  * \param lambda_values_out Local part that will be added to A * x.
188  * \param hyper_edge HyperEdge that is considered.
189  * \param time Time at which analytic functions are evaluated.
190  * \retval vecAx Local part of vector A * x.
191  ************************************************************************************************/
192  template <class hyEdgeT, typename SmallMatInT, typename SmallMatOutT>
193  SmallMatOutT& trace_to_flux(const SmallMatInT& lambda_values_in,
194  SmallMatOutT& lambda_values_out,
195  hyEdgeT& hyper_edge,
196  const lSol_float_t time = 0.) const
197  {
198  static_assert(hyEdge_dimT == 1, "Elastic graphs must be graphs, not hypergraphs!");
199  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>, 2 *hyEdge_dimT>
200  lambda_old = node_dof_to_edge_dof(lambda_values_in, hyper_edge),
201  lambda_new;
202  for (unsigned int i = 0; i < lambda_new.size(); ++i)
203  lambda_new[i].fill(0.);
204 
205  if constexpr (has_trace_to_flux<
206  diffusion_sol_t,
207  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
208  2 * hyEdge_dimT>&(
209  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
210  2 * hyEdge_dimT>&,
211  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
212  2 * hyEdge_dimT>&)>::value)
213  diffusion.trace_to_flux(lambda_old, lambda_new, time);
214  else
215  diffusion.trace_to_flux(lambda_old, lambda_new, hyper_edge, time);
216 
217  return lambda_values_out = edge_dof_to_node_dof(lambda_new, lambda_values_out, hyper_edge);
218  }
219  /*!***********************************************************************************************
220  * \brief Evaluate local contribution to residual.
221  ************************************************************************************************/
222  template <class hyEdgeT, typename SmallMatInT, typename SmallMatOutT>
223  SmallMatOutT& residual_flux(const SmallMatInT& lambda_values_in,
224  SmallMatOutT& lambda_values_out,
225  hyEdgeT& hyper_edge,
226  const lSol_float_t time = 0.) const
227  {
228  static_assert(hyEdge_dimT == 1, "Elastic graphs must be graphs, not hypergraphs!");
229  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>, 2 *hyEdge_dimT>
230  lambda_old = node_dof_to_edge_dof(lambda_values_in, hyper_edge),
231  lambda_new;
232  for (unsigned int i = 0; i < lambda_new.size(); ++i)
233  lambda_new[i].fill(0.);
234 
235  if constexpr (has_residual_flux<
236  diffusion_sol_t,
237  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
238  2 * hyEdge_dimT>&(
239  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
240  2 * hyEdge_dimT>&,
241  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
242  2 * hyEdge_dimT>&)>::value)
243  diffusion.residual_flux(lambda_old, lambda_new, time);
244  else
245  diffusion.residual_flux(lambda_old, lambda_new, hyper_edge, time);
246 
247  return lambda_values_out = edge_dof_to_node_dof(lambda_new, lambda_values_out, hyper_edge);
248  }
249  /*!***********************************************************************************************
250  * \brief Local squared contribution to the L2 error.
251  *
252  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
253  * \tparam SmallMatT The data type of \c lambda_values.
254  * \param lambda_values The values of the skeletal variable's coefficients.
255  * \param hyper_edge The geometry of the considered hyperedge (of typename GeomT).
256  * \param time Time at which analytic functions are evaluated.
257  * \retval vec_b Local part of vector b.
258  ************************************************************************************************/
259  template <class hyEdgeT, typename SmallMatT>
260  std::array<lSol_float_t, 1U> errors(const SmallMatT& lambda_values,
261  hyEdgeT& hyper_edge,
262  const lSol_float_t time = 0.) const
263  {
264  lSol_float_t error = 0.;
265  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT>
266  lambda = node_dof_to_edge_dof(lambda_values, hyper_edge);
267 
268  if constexpr (has_errors<
269  diffusion_sol_t,
270  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
271  2 * hyEdge_dimT>&(
272  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
273  2 * hyEdge_dimT>&,
274  std::array<std::array<lSol_float_t, diffusion_sol_t::n_glob_dofs_per_node()>,
275  2 * hyEdge_dimT>&)>::value)
276  error = diffusion.errors(lambda, time);
277  else
278  error = diffusion.errors(lambda, hyper_edge, time);
279 
280  return std::array<lSol_float_t, 1U>({error});
281  }
282  /*!***********************************************************************************************
283  * \brief Evaluate local local reconstruction at tensorial products of abscissas.
284  *
285  * \tparam absc_float_t Floating type for the abscissa values.
286  * \tparam abscissas_sizeT Size of the array of array of abscissas.
287  * \tparam input_array_t Type of input array.
288  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
289  * \param abscissas Abscissas of the supporting points.
290  * \param lambda_values The values of the skeletal variable's coefficients.
291  * \param hyper_edge The geometry of the considered hyperedge (of typename GeomT).
292  * \param time Time.
293  * \retval func_values Function values at tensorial points.
294  ************************************************************************************************/
295  template <typename abscissa_float_t, std::size_t sizeT, class input_array_t, class hyEdgeT>
296  std::array<
297  std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(sizeT)>,
299  bulk_values(const std::array<abscissa_float_t, sizeT>& abscissas,
300  const input_array_t& lambda_values,
301  hyEdgeT& hyper_edge,
302  const lSol_float_t time = 0.) const
303  {
304  std::array<
305  std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(sizeT)>,
307  result;
308 
309  auto bulk =
310  diffusion.bulk_values(abscissas, node_dof_to_edge_dof(lambda_values, hyper_edge), time);
311  Point<space_dim, lSol_float_t> normal_vector =
312  (Point<space_dim, lSol_float_t>)hyper_edge.geometry.inner_normal(1);
313 
314  for (unsigned int dim = 0; dim < result.size(); ++dim)
315  for (unsigned int q = 0; q < result[dim].size(); ++q)
316  result[dim][q] = bulk[1][q] * normal_vector[dim];
317 
318  return result;
319  }
320 }; // end of class LengtheningBeam
321 
322 /*!*************************************************************************************************
323  * \brief Local solver for the equation that governs the bending of an elastic Bernoulli beam.
324  *
325  * \authors Guido Kanschat, Heidelberg University, 2019--2020.
326  * \authors Andreas Rupp, Heidelberg University, 2019--2020.
327  **************************************************************************************************/
328 template <unsigned int hyEdge_dimT,
329  unsigned int space_dim,
330  unsigned int poly_deg,
331  unsigned int quad_deg,
332  typename lSol_float_t = double,
333  typename bilaplacian_sol_t =
334  BilaplacianUniform<hyEdge_dimT, poly_deg, quad_deg, lSol_float_t> >
336 {
337  /*!***********************************************************************************************
338  * \brief Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
339  ************************************************************************************************/
340  HAS_MEMBER_FUNCTION(trace_to_flux, has_trace_to_flux);
341  /*!***********************************************************************************************
342  * \brief Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
343  ************************************************************************************************/
344  HAS_MEMBER_FUNCTION(residual_flux, has_residual_flux);
345  /*!***********************************************************************************************
346  * \brief Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
347  ************************************************************************************************/
348  HAS_MEMBER_FUNCTION(errors, has_errors);
349 
350  public:
351  /*!***********************************************************************************************
352  * \brief Define type of (hyperedge related) data that is stored in HyDataContainer.
353  ************************************************************************************************/
354  struct data_type
355  {
356  };
357  /*!***********************************************************************************************
358  * \brief Define type of node elements, especially with respect to nodal shape functions.
359  ************************************************************************************************/
361  {
362  typedef std::tuple<TPP::ShapeFunction<
365  };
366  /*!***********************************************************************************************
367  * \brief Define how errors are evaluated.
368  ************************************************************************************************/
369  struct error_def
370  {
371  /*!*********************************************************************************************
372  * \brief Define the typename returned by function errors.
373  **********************************************************************************************/
374  typedef std::array<lSol_float_t, 1U> error_t;
375  /*!*********************************************************************************************
376  * \brief Define how initial error is generated.
377  **********************************************************************************************/
379  {
380  std::array<lSol_float_t, 1U> summed_error;
381  summed_error.fill(0.);
382  return summed_error;
383  }
384  /*!*********************************************************************************************
385  * \brief Define how local errors should be accumulated.
386  **********************************************************************************************/
387  static error_t sum_error(error_t& summed_error, const error_t& new_error)
388  {
389  for (unsigned int k = 0; k < summed_error.size(); ++k)
390  summed_error[k] += new_error[k];
391  return summed_error;
392  }
393  /*!*********************************************************************************************
394  * \brief Define how global errors should be postprocessed.
395  **********************************************************************************************/
396  static error_t postprocess_error(error_t& summed_error)
397  {
398  for (unsigned int k = 0; k < summed_error.size(); ++k)
399  summed_error[k] = std::sqrt(summed_error[k]);
400  return summed_error;
401  }
402  };
403  /*!***********************************************************************************************
404  * \brief Return template parameter \c hyEdge_dimT.
405  *
406  * \retval hyEdge_dimT Dimension of hypergraph's hyperedges.
407  ************************************************************************************************/
408  static constexpr unsigned int hyEdge_dim() { return hyEdge_dimT; }
409  /*!***********************************************************************************************
410  * \brief Evaluate amount of global degrees of freedom per hypernode.
411  *
412  * This number must be equal to HyperNodeFactory::n_glob_dofs_per_node()() of the HyperNodeFactory
413  * cooperating with this object.
414  *
415  * \retval n_dofs Number of global degrees of freedom per hypernode.
416  ************************************************************************************************/
417  static constexpr unsigned int n_glob_dofs_per_node()
418  {
419  return 2 * space_dim * Hypercube<hyEdge_dimT - 1>::pow(poly_deg + 1);
420  }
421  /*!***********************************************************************************************
422  * \brief Dimension of of the solution evaluated with respect to a hyperedge.
423  ************************************************************************************************/
424  static constexpr unsigned int system_dimension() { return space_dim; }
425  /*!***********************************************************************************************
426  * \brief Dimension of of the solution evaluated with respect to a hypernode.
427  ************************************************************************************************/
428  static constexpr unsigned int node_system_dimension() { return space_dim * 2; }
429 
430  private:
431  /*!***********************************************************************************************
432  * \brief The bilaplacian solver that solves the locally defined PDE.
433  ************************************************************************************************/
434  const bilaplacian_sol_t bilaplacian_solver;
435  /*!***********************************************************************************************
436  * \brief Do the preprocessing to transfer global to local dofs.
437  ************************************************************************************************/
438  template <class hyEdgeT, typename SmallMatT>
439  inline std::array<std::array<double, bilaplacian_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT>
440  node_dof_to_edge_dof(const SmallMatT& lambda,
441  hyEdgeT& hyper_edge,
442  const unsigned int outer_index) const
443  {
444  std::array<std::array<double, bilaplacian_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT>
445  result;
446  hy_assert(result.size() == 2, "Only implemented in one dimension!");
447  for (unsigned int i = 0; i < result.size(); ++i)
448  {
449  hy_assert(result[i].size() == 2, "Only implemented in one dimension!");
450  result[i].fill(0.);
451  }
452 
453  Point<space_dim, lSol_float_t> normal_vector =
454  (Point<space_dim, lSol_float_t>)hyper_edge.geometry.outer_normal(outer_index);
455 
456  for (unsigned int i = 0; i < 2 * hyEdge_dimT; ++i)
457  for (unsigned int dim = 0; dim < space_dim; ++dim)
458  {
459  result[i][0] += normal_vector[dim] * lambda[i][dim];
460  result[i][1] += normal_vector[dim] * lambda[i][space_dim + dim];
461  }
462 
463  return result;
464  }
465  /*!***********************************************************************************************
466  * \brief Do the postprocessing to transfer local to global dofs.
467  ************************************************************************************************/
468  template <class hyEdgeT, typename SmallMatT>
469  inline SmallMatT& edge_dof_to_node_dof(
470  const std::array<std::array<double, bilaplacian_sol_t::n_glob_dofs_per_node()>,
471  2 * hyEdge_dimT>& lambda,
472  SmallMatT& lambda_values_out,
473  hyEdgeT& hyper_edge,
474  const unsigned int outer_index) const
475  {
476  hy_assert(bilaplacian_sol_t::n_glob_dofs_per_node() == 2, "This should be 1*2!");
477  Point<space_dim, lSol_float_t> normal_vector =
478  (Point<space_dim, lSol_float_t>)hyper_edge.geometry.outer_normal(outer_index);
479 
480  for (unsigned int i = 0; i < 2 * hyEdge_dimT; ++i)
481  for (unsigned int dim = 0; dim < space_dim; ++dim)
482  {
483  lambda_values_out[i][dim] += normal_vector[dim] * lambda[i][0];
484  lambda_values_out[i][space_dim + dim] += normal_vector[dim] * lambda[i][1];
485  }
486 
487  return lambda_values_out;
488  }
489 
490  public:
491  /*!***********************************************************************************************
492  * \brief Class is constructed using a single double indicating the penalty parameter.
493  ************************************************************************************************/
494  typedef lSol_float_t constructor_value_type;
495  /*!***********************************************************************************************
496  * \brief Constructor for local solver.
497  *
498  * \param tau Penalty parameter of HDG scheme.
499  ************************************************************************************************/
501  /*!***********************************************************************************************
502  * \brief Evaluate local contribution to matrix--vector multiplication.
503  *
504  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
505  * \tparam SmallMatInT Data type of \c lambda_values_in.
506  * \tparam SmallMatOutT Data type of \c lambda_values_out.
507  * \param lambda_values_in Local part of vector x.
508  * \param lambda_values_out Local part that will be added to A * x.
509  * \param hyper_edge HyperEdge that is considered.
510  * \param time Time at which analytic functions are evaluated.
511  * \retval vecAx Local part of vector A * x.
512  ************************************************************************************************/
513  template <class hyEdgeT, typename SmallMatInT, typename SmallMatOutT>
514  SmallMatOutT& trace_to_flux(const SmallMatInT& lambda_values_in,
515  SmallMatOutT& lambda_values_out,
516  hyEdgeT& hyper_edge,
517  const lSol_float_t time = 0.) const
518  {
519  static_assert(hyEdge_dimT == 1, "The beam must be one-dimensional!");
520  std::array<std::array<lSol_float_t, bilaplacian_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT>
521  lambda_old, lambda_new;
522 
523  for (unsigned int dim = 0; dim < space_dim - hyEdge_dimT; ++dim)
524  {
525  lambda_old = node_dof_to_edge_dof(lambda_values_in, hyper_edge, dim);
526  for (unsigned int i = 0; i < lambda_new.size(); ++i)
527  lambda_new[i].fill(0.);
528 
529  if constexpr (
530  has_trace_to_flux<
531  bilaplacian_sol_t,
532  std::array<std::array<lSol_float_t, bilaplacian_sol_t::n_glob_dofs_per_node()>,
533  2 * hyEdge_dimT>&(
534  std::array<std::array<lSol_float_t, bilaplacian_sol_t::n_glob_dofs_per_node()>,
535  2 * hyEdge_dimT>&,
536  std::array<std::array<lSol_float_t, bilaplacian_sol_t::n_glob_dofs_per_node()>,
537  2 * hyEdge_dimT>&)>::value)
538  bilaplacian_solver.trace_to_flux(lambda_old, lambda_new, time);
539  else
540  bilaplacian_solver.trace_to_flux(lambda_old, lambda_new, hyper_edge, time);
541 
542  edge_dof_to_node_dof(lambda_new, lambda_values_out, hyper_edge, dim);
543  }
544 
545  return lambda_values_out;
546  }
547  /*!***********************************************************************************************
548  * \brief Evaluate local contribution to residual.
549  ************************************************************************************************/
550  template <class hyEdgeT, typename SmallMatInT, typename SmallMatOutT>
551  SmallMatOutT& residual_flux(const SmallMatInT& lambda_values_in,
552  SmallMatOutT& lambda_values_out,
553  hyEdgeT& hyper_edge,
554  const lSol_float_t time = 0.) const
555  {
556  static_assert(hyEdge_dimT == 1, "The beam must be one-dimensional!");
557  std::array<std::array<lSol_float_t, bilaplacian_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT>
558  lambda_old, lambda_new;
559 
560  for (unsigned int dim = 0; dim < space_dim - hyEdge_dimT; ++dim)
561  {
562  lambda_old = node_dof_to_edge_dof(lambda_values_in, hyper_edge, dim);
563  for (unsigned int i = 0; i < lambda_new.size(); ++i)
564  lambda_new[i].fill(0.);
565 
566  if constexpr (
567  has_residual_flux<
568  bilaplacian_sol_t,
569  std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&(
570  std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&,
571  std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&)>::value)
572  bilaplacian_solver.residual_flux(lambda_old, lambda_new, time);
573  else
574  bilaplacian_solver.residual_flux(lambda_old, lambda_new, hyper_edge, time);
575 
576  edge_dof_to_node_dof(lambda_new, lambda_values_out, hyper_edge, dim);
577  }
578 
579  return lambda_values_out;
580  }
581  /*!***********************************************************************************************
582  * \brief Local squared contribution to the L2 error.
583  *
584  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
585  * \param lambda_values The values of the skeletal variable's coefficients.
586  * \param hyper_edge The geometry of the considered hyperedge (of typename GeomT).
587  * \param time Time at which analytic functions are evaluated.
588  * \retval vec_b Local part of vector b.
589  ************************************************************************************************/
590  template <class hyEdgeT>
591  std::array<lSol_float_t, 1U> errors(
592  const std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&
593  lambda_values,
594  hyEdgeT& hyper_edge,
595  const lSol_float_t time = 0.) const
596  {
597  lSol_float_t error = 0.;
598  std::array<std::array<lSol_float_t, bilaplacian_sol_t::n_glob_dofs_per_node()>, 2 * hyEdge_dimT>
599  lambda;
600  for (unsigned int dim = 0; dim < space_dim - hyEdge_dimT; ++dim)
601  {
602  lambda = node_dof_to_edge_dof(lambda_values, hyper_edge, dim);
603 
604  if constexpr (
605  has_errors<
606  bilaplacian_sol_t,
607  std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&(
608  std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&,
609  std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&)>::value)
610  error += bilaplacian_solver.errors(lambda, time);
611  else
612  error += bilaplacian_solver.errors(lambda, hyper_edge, time);
613  }
614 
615  return std::array<lSol_float_t, 1U>({error});
616  }
617  /*!***********************************************************************************************
618  * \brief Evaluate local local reconstruction at tensorial products of abscissas.
619  *
620  * \tparam absc_float_t Floating type for the abscissa values.
621  * \tparam abscissas_sizeT Size of the array of array of abscissas.
622  * \tparam input_array_t Type of input array.
623  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
624  * \param abscissas Abscissas of the supporting points.
625  * \param lambda_values The values of the skeletal variable's coefficients.
626  * \param hyper_edge The geometry of the considered hyperedge (of typename GeomT).
627  * \param time Time.
628  * \retval func_values Function values at tensorial points.
629  ************************************************************************************************/
630  template <typename abscissa_float_t, std::size_t sizeT, class input_array_t, class hyEdgeT>
631  std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(sizeT)>,
634  bulk_values(const std::array<abscissa_float_t, sizeT>& abscissas,
635  const input_array_t& lambda_values,
636  hyEdgeT& hyper_edge,
637  const lSol_float_t time = 0.) const
638  {
639  std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(sizeT)>, system_dimension()>
640  values;
641  for (unsigned int i = 0; i < values.size(); ++i)
642  values[i].fill(0.);
643 
644  for (unsigned int dim_on = 0; dim_on < space_dim - hyEdge_dimT; ++dim_on)
645  {
646  auto bulk = bilaplacian_solver.bulk_values(
647  abscissas, node_dof_to_edge_dof(lambda_values, hyper_edge, dim_on), time);
648  Point<space_dim, lSol_float_t> normal_vector =
649  (Point<space_dim, lSol_float_t>)hyper_edge.geometry.outer_normal(dim_on);
650 
651  for (unsigned int dim = 0; dim < values.size(); ++dim)
652  for (unsigned int q = 0; q < values[dim].size(); ++q)
653  values[dim][q] += bulk[1][q] * normal_vector[dim];
654  }
655 
656  return values;
657  }
658 }; // end of class BernoulliBendingBeam
659 
660 /*!*************************************************************************************************
661  * \brief Local solver for the equation that governs the bending and change of length of an
662  * elastic Bernoulli beam.
663  *
664  * \authors Guido Kanschat, Heidelberg University, 2019--2020.
665  * \authors Andreas Rupp, Heidelberg University, 2019--2020.
666  **************************************************************************************************/
667 template <unsigned int hyEdge_dimT,
668  unsigned int space_dim,
669  unsigned int poly_deg,
670  unsigned int quad_deg,
671  typename lSol_float_t = double>
673 {
674  public:
675  /*!***********************************************************************************************
676  * \brief Define type of (hyperedge related) data that is stored in HyDataContainer.
677  ************************************************************************************************/
678  struct data_type
679  {
680  };
681  /*!***********************************************************************************************
682  * \brief Define type of node elements, especially with respect to nodal shape functions.
683  ************************************************************************************************/
685  {
686  typedef std::tuple<TPP::ShapeFunction<
689  };
690  /*!***********************************************************************************************
691  * \brief Define how errors are evaluated.
692  ************************************************************************************************/
693  struct error_def
694  {
695  /*!*********************************************************************************************
696  * \brief Define the typename returned by function errors.
697  **********************************************************************************************/
698  typedef std::array<lSol_float_t, 1U> error_t;
699  /*!*********************************************************************************************
700  * \brief Define how initial error is generated.
701  **********************************************************************************************/
703  {
704  std::array<lSol_float_t, 1U> summed_error;
705  summed_error.fill(0.);
706  return summed_error;
707  }
708  /*!*********************************************************************************************
709  * \brief Define how local errors should be accumulated.
710  **********************************************************************************************/
711  static error_t sum_error(error_t& summed_error, const error_t& new_error)
712  {
713  for (unsigned int k = 0; k < summed_error.size(); ++k)
714  summed_error[k] += new_error[k];
715  return summed_error;
716  }
717  /*!*********************************************************************************************
718  * \brief Define how global errors should be postprocessed.
719  **********************************************************************************************/
720  static error_t postprocess_error(error_t& summed_error)
721  {
722  for (unsigned int k = 0; k < summed_error.size(); ++k)
723  summed_error[k] = std::sqrt(summed_error[k]);
724  return summed_error;
725  }
726  };
727  /*!***********************************************************************************************
728  * \brief Return template parameter \c hyEdge_dimT.
729  *
730  * \retval hyEdge_dimT Dimension of hypergraph's hyperedges.
731  ************************************************************************************************/
732  static constexpr unsigned int hyEdge_dim() { return hyEdge_dimT; }
733  /*!***********************************************************************************************
734  * \brief Evaluate amount of global degrees of freedom per hypernode.
735  *
736  * This number must be equal to HyperNodeFactory::n_glob_dofs_per_node()() of the HyperNodeFactory
737  * cooperating with this object.
738  *
739  * \retval n_dofs Number of global degrees of freedom per hypernode.
740  ************************************************************************************************/
741  static constexpr unsigned int n_glob_dofs_per_node()
742  {
743  return 2 * space_dim * Hypercube<hyEdge_dimT - 1>::pow(poly_deg + 1);
744  }
745  /*!***********************************************************************************************
746  * \brief Dimension of of the solution evaluated with respect to a hyperedge.
747  ************************************************************************************************/
748  static constexpr unsigned int system_dimension() { return space_dim; }
749  /*!***********************************************************************************************
750  * \brief Dimension of of the solution evaluated with respect to a hypernode.
751  ************************************************************************************************/
752  static constexpr unsigned int node_system_dimension() { return space_dim * 2; }
753 
754  private:
755  /*!***********************************************************************************************
756  * \brief The lengthening beam solver that does the lengthening of the beam.
757  ************************************************************************************************/
759  /*!***********************************************************************************************
760  * \brief The bending beam solver that does the bending of the beam.
761  ************************************************************************************************/
763 
764  public:
765  /*!***********************************************************************************************
766  * \brief Class is constructed using a single double indicating the penalty parameter.
767  ************************************************************************************************/
768  typedef lSol_float_t constructor_value_type;
769  /*!***********************************************************************************************
770  * \brief Constructor for local solver.
771  *
772  * \param tau Penalty parameter of HDG scheme.
773  ************************************************************************************************/
775  : len_beam(tau), ben_beam(tau)
776  {
777  }
778  /*!***********************************************************************************************
779  * \brief Evaluate local contribution to matrix--vector multiplication.
780  ************************************************************************************************/
781  template <class hyEdgeT, typename SmallMatInT, typename SmallMatOutT>
782  SmallMatOutT& trace_to_flux(const SmallMatInT& lambda_values_in,
783  SmallMatOutT& lambda_values_out,
784  hyEdgeT& hyper_edge,
785  const lSol_float_t time = 0.) const
786  {
787  static_assert(hyEdge_dimT == 1, "A beam must be one-dimensional!");
788 
789  len_beam.trace_to_flux(lambda_values_in, lambda_values_out, hyper_edge, time);
790  ben_beam.trace_to_flux(lambda_values_in, lambda_values_out, hyper_edge, time);
791 
792  return lambda_values_out;
793  }
794  /*!***********************************************************************************************
795  * \brief Evaluate local contribution to residual.
796  ************************************************************************************************/
797  template <class hyEdgeT, typename SmallMatInT, typename SmallMatOutT>
798  SmallMatOutT& residual_flux(const SmallMatInT& lambda_values_in,
799  SmallMatOutT& lambda_values_out,
800  hyEdgeT& hyper_edge,
801  const lSol_float_t time = 0.) const
802  {
803  static_assert(hyEdge_dimT == 1, "A beam must be one-dimensional!");
804 
805  len_beam.residual_flux(lambda_values_in, lambda_values_out, hyper_edge, time);
806  ben_beam.residual_flux(lambda_values_in, lambda_values_out, hyper_edge, time);
807 
808  return lambda_values_out;
809  }
810  /*!***********************************************************************************************
811  * \brief Local squared contribution to the L2 error.
812  *
813  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
814  * \param lambda_values The values of the skeletal variable's coefficients.
815  * \param hyper_edge The geometry of the considered hyperedge (of typename GeomT).
816  * \param time Time at which analytic functions are evaluated.
817  * \retval vec_b Local part of vector b.
818  ************************************************************************************************/
819  template <class hyEdgeT>
820  std::array<lSol_float_t, 1U> errors(
821  const std::array<std::array<lSol_float_t, n_glob_dofs_per_node()>, 2 * hyEdge_dimT>&
822  lambda_values,
823  hyEdgeT& hyper_edge,
824  const lSol_float_t time = 0.) const
825  {
826  lSol_float_t error = 0.;
827  error += len_beam.errors(lambda_values, hyper_edge, time);
828  error += ben_beam.errors(lambda_values, hyper_edge, time);
829  return std::array<lSol_float_t, 1U>({error});
830  }
831  /*!***********************************************************************************************
832  * \brief Evaluate local local reconstruction at tensorial products of abscissas.
833  *
834  * \tparam absc_float_t Floating type for the abscissa values.
835  * \tparam abscissas_sizeT Size of the array of array of abscissas.
836  * \tparam input_array_t Type of input array.
837  * \tparam hyEdgeT The geometry type / typename of the considered hyEdge's geometry.
838  * \param abscissas Abscissas of the supporting points.
839  * \param lambda_values The values of the skeletal variable's coefficients.
840  * \param hyper_edge The geometry of the considered hyperedge (of typename GeomT).
841  * \param time Time.
842  * \retval func_values Function values at tensorial points.
843  ************************************************************************************************/
844  template <typename abscissa_float_t, std::size_t sizeT, class input_array_t, class hyEdgeT>
845  std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(sizeT)>, system_dimension()>
846  bulk_values(const std::array<abscissa_float_t, sizeT>& abscissas,
847  const input_array_t& lambda_values,
848  hyEdgeT& hyper_edge,
849  const lSol_float_t time = 0.) const
850  {
851  std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(sizeT)>, system_dimension()>
852  result, auxiliary;
853  result = len_beam.bulk_values(abscissas, lambda_values, hyper_edge, time);
854  auxiliary = ben_beam.bulk_values(abscissas, lambda_values, hyper_edge, time);
855 
856  for (unsigned int i = 0; i < system_dimension(); ++i)
857  for (unsigned int j = 0; j < Hypercube<hyEdge_dimT>::pow(sizeT); ++j)
858  result[i][j] += auxiliary[i][j];
859 
860  return result;
861  }
862 }; // end of class LengtheningBernoulliBendingBeam
863 
864 } // namespace LocalSolver
LocalSolver::BernoulliBendingBeam::residual_flux
SmallMatOutT & residual_flux(const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local contribution to residual.
Definition: bernoulli_beams.hxx:551
LocalSolver::BernoulliBendingBeam::error_def::initial_error
static error_t initial_error()
Define how initial error is generated.
Definition: bernoulli_beams.hxx:378
LocalSolver::BernoulliBendingBeam::node_dof_to_edge_dof
std::array< std::array< double, bilaplacian_sol_t::n_glob_dofs_per_node()>, 2 *hyEdge_dimT > node_dof_to_edge_dof(const SmallMatT &lambda, hyEdgeT &hyper_edge, const unsigned int outer_index) const
Do the preprocessing to transfer global to local dofs.
Definition: bernoulli_beams.hxx:440
compile_time_tricks.hxx
LocalSolver::BernoulliBendingBeam::edge_dof_to_node_dof
SmallMatT & edge_dof_to_node_dof(const std::array< std::array< double, bilaplacian_sol_t::n_glob_dofs_per_node()>, 2 *hyEdge_dimT > &lambda, SmallMatT &lambda_values_out, hyEdgeT &hyper_edge, const unsigned int outer_index) const
Do the postprocessing to transfer local to global dofs.
Definition: bernoulli_beams.hxx:469
LocalSolver::LengtheningBeam::bulk_values
std::array< std::array< lSol_float_t, Hypercube< hyEdge_dimT >::pow(sizeT)>, LengtheningBeam< hyEdge_dimT, space_dim, poly_deg, quad_deg, lSol_float_t >::system_dimension()> bulk_values(const std::array< abscissa_float_t, sizeT > &abscissas, const input_array_t &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local local reconstruction at tensorial products of abscissas.
Definition: bernoulli_beams.hxx:299
shape_function.hxx
LocalSolver::LengtheningBernoulliBendingBeam::bulk_values
std::array< std::array< lSol_float_t, Hypercube< hyEdge_dimT >::pow(sizeT)>, system_dimension()> bulk_values(const std::array< abscissa_float_t, sizeT > &abscissas, const input_array_t &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local local reconstruction at tensorial products of abscissas.
Definition: bernoulli_beams.hxx:846
LocalSolver::LengtheningBernoulliBendingBeam::LengtheningBernoulliBendingBeam
LengtheningBernoulliBendingBeam(const constructor_value_type &tau=1.)
Constructor for local solver.
Definition: bernoulli_beams.hxx:774
LocalSolver::LengtheningBernoulliBendingBeam::error_def
Define how errors are evaluated.
Definition: bernoulli_beams.hxx:693
LocalSolver::LengtheningBeam::error_def::error_t
std::array< lSol_float_t, 1U > error_t
Define the typename returned by function errors.
Definition: bernoulli_beams.hxx:65
LocalSolver::LengtheningBernoulliBendingBeam::node_element
Define type of node elements, especially with respect to nodal shape functions.
Definition: bernoulli_beams.hxx:684
LocalSolver
A namespace for local solvers.
Definition: advection_parab_ldgh.hxx:11
LocalSolver::BernoulliBendingBeam::error_def::error_t
std::array< lSol_float_t, 1U > error_t
Define the typename returned by function errors.
Definition: bernoulli_beams.hxx:374
LocalSolver::LengtheningBernoulliBendingBeam::data_type
Define type of (hyperedge related) data that is stored in HyDataContainer.
Definition: bernoulli_beams.hxx:678
LocalSolver::LengtheningBeam::errors
std::array< lSol_float_t, 1U > errors(const SmallMatT &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Local squared contribution to the L2 error.
Definition: bernoulli_beams.hxx:260
LocalSolver::LengtheningBernoulliBendingBeam
Local solver for the equation that governs the bending and change of length of an elastic Bernoulli b...
Definition: bernoulli_beams.hxx:672
LocalSolver::LengtheningBernoulliBendingBeam::error_def::error_t
std::array< lSol_float_t, 1U > error_t
Define the typename returned by function errors.
Definition: bernoulli_beams.hxx:698
LocalSolver::LengtheningBernoulliBendingBeam::n_glob_dofs_per_node
static constexpr unsigned int n_glob_dofs_per_node()
Evaluate amount of global degrees of freedom per hypernode.
Definition: bernoulli_beams.hxx:741
LocalSolver::LengtheningBeam::trace_to_flux
SmallMatOutT & trace_to_flux(const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local contribution to matrix–vector multiplication.
Definition: bernoulli_beams.hxx:193
LocalSolver::LengtheningBeam::error_def::postprocess_error
static error_t postprocess_error(error_t &summed_error)
Define how global errors should be postprocessed.
Definition: bernoulli_beams.hxx:87
LocalSolver::LengtheningBernoulliBendingBeam::error_def::sum_error
static error_t sum_error(error_t &summed_error, const error_t &new_error)
Define how local errors should be accumulated.
Definition: bernoulli_beams.hxx:711
LocalSolver::LengtheningBeam::error_def::sum_error
static error_t sum_error(error_t &summed_error, const error_t &new_error)
Define how local errors should be accumulated.
Definition: bernoulli_beams.hxx:78
LocalSolver::LengtheningBernoulliBendingBeam::error_def::postprocess_error
static error_t postprocess_error(error_t &summed_error)
Define how global errors should be postprocessed.
Definition: bernoulli_beams.hxx:720
LocalSolver::BernoulliBendingBeam::data_type
Define type of (hyperedge related) data that is stored in HyDataContainer.
Definition: bernoulli_beams.hxx:354
LocalSolver::BernoulliBendingBeam
Local solver for the equation that governs the bending of an elastic Bernoulli beam.
Definition: bernoulli_beams.hxx:335
LocalSolver::LengtheningBeam::data_type
Define type of (hyperedge related) data that is stored in HyDataContainer.
Definition: bernoulli_beams.hxx:45
LocalSolver::LengtheningBernoulliBendingBeam::ben_beam
const BernoulliBendingBeam< hyEdge_dimT, space_dim, poly_deg, quad_deg, lSol_float_t > ben_beam
The bending beam solver that does the bending of the beam.
Definition: bernoulli_beams.hxx:762
TPP::ShapeType::Tensorial
Struct that handles the evaluation of tensorial shape functions.
Definition: tensorial.hxx:22
LocalSolver::LengtheningBernoulliBendingBeam::node_element::functions
std::tuple< TPP::ShapeFunction< TPP::ShapeType::Tensorial< TPP::ShapeType::Legendre< poly_deg >, hyEdge_dimT - 1 > > > functions
Definition: bernoulli_beams.hxx:688
LocalSolver::BernoulliBendingBeam::node_system_dimension
static constexpr unsigned int node_system_dimension()
Dimension of of the solution evaluated with respect to a hypernode.
Definition: bernoulli_beams.hxx:428
LocalSolver::BernoulliBendingBeam::error_def::sum_error
static error_t sum_error(error_t &summed_error, const error_t &new_error)
Define how local errors should be accumulated.
Definition: bernoulli_beams.hxx:387
LocalSolver::BernoulliBendingBeam::node_element::functions
std::tuple< TPP::ShapeFunction< TPP::ShapeType::Tensorial< TPP::ShapeType::Legendre< poly_deg >, hyEdge_dimT - 1 > > > functions
Definition: bernoulli_beams.hxx:364
LocalSolver::LengtheningBernoulliBendingBeam::error_def::initial_error
static error_t initial_error()
Define how initial error is generated.
Definition: bernoulli_beams.hxx:702
LocalSolver::LengtheningBeam::system_dimension
static constexpr unsigned int system_dimension()
Dimension of of the solution evaluated with respect to a hyperedge.
Definition: bernoulli_beams.hxx:115
LocalSolver::LengtheningBernoulliBendingBeam::residual_flux
SmallMatOutT & residual_flux(const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local contribution to residual.
Definition: bernoulli_beams.hxx:798
TPP::ShapeFunction
Struct that handles different types of evaluation of shape functions.
Definition: shape_function.hxx:15
LocalSolver::LengtheningBeam
Local solver for the equation that governs the lengthening of an elastic beam.
Definition: bernoulli_beams.hxx:26
LocalSolver::BernoulliBendingBeam::bilaplacian_solver
const bilaplacian_sol_t bilaplacian_solver
The bilaplacian solver that solves the locally defined PDE.
Definition: bernoulli_beams.hxx:434
LocalSolver::BernoulliBendingBeam::n_glob_dofs_per_node
static constexpr unsigned int n_glob_dofs_per_node()
Evaluate amount of global degrees of freedom per hypernode.
Definition: bernoulli_beams.hxx:417
LocalSolver::LengtheningBeam::node_element::functions
std::tuple< TPP::ShapeFunction< TPP::ShapeType::Tensorial< TPP::ShapeType::Legendre< poly_deg >, hyEdge_dimT - 1 > > > functions
Definition: bernoulli_beams.hxx:55
LocalSolver::LengtheningBeam::diffusion
const diffusion_sol_t diffusion
The diffusion solver that solves the locally defined PDE.
Definition: bernoulli_beams.hxx:125
LocalSolver::LengtheningBernoulliBendingBeam::system_dimension
static constexpr unsigned int system_dimension()
Dimension of of the solution evaluated with respect to a hyperedge.
Definition: bernoulli_beams.hxx:748
LocalSolver::BernoulliBendingBeam::error_def::postprocess_error
static error_t postprocess_error(error_t &summed_error)
Define how global errors should be postprocessed.
Definition: bernoulli_beams.hxx:396
LocalSolver::BernoulliBendingBeam::bulk_values
std::array< std::array< lSol_float_t, Hypercube< hyEdge_dimT >::pow(sizeT)>, BernoulliBendingBeam< hyEdge_dimT, space_dim, poly_deg, quad_deg, lSol_float_t >::system_dimension()> bulk_values(const std::array< abscissa_float_t, sizeT > &abscissas, const input_array_t &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local local reconstruction at tensorial products of abscissas.
Definition: bernoulli_beams.hxx:634
LocalSolver::LengtheningBeam::node_element
Define type of node elements, especially with respect to nodal shape functions.
Definition: bernoulli_beams.hxx:51
LocalSolver::LengtheningBeam::LengtheningBeam
LengtheningBeam(const constructor_value_type &tau=1.)
Constructor for local solver.
Definition: bernoulli_beams.hxx:179
LocalSolver::LengtheningBeam::hyEdge_dim
static constexpr unsigned int hyEdge_dim()
Return template parameter hyEdge_dimT.
Definition: bernoulli_beams.hxx:99
LocalSolver::BernoulliBendingBeam::trace_to_flux
SmallMatOutT & trace_to_flux(const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local contribution to matrix–vector multiplication.
Definition: bernoulli_beams.hxx:514
LocalSolver::LengtheningBeam::n_glob_dofs_per_node
static constexpr unsigned int n_glob_dofs_per_node()
Evaluate amount of global degrees of freedom per hypernode.
Definition: bernoulli_beams.hxx:108
LocalSolver::BernoulliBendingBeam::BernoulliBendingBeam
BernoulliBendingBeam(const constructor_value_type &tau=1.)
Constructor for local solver.
Definition: bernoulli_beams.hxx:500
hypercube.hxx
LocalSolver::LengtheningBernoulliBendingBeam::trace_to_flux
SmallMatOutT & trace_to_flux(const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local contribution to matrix–vector multiplication.
Definition: bernoulli_beams.hxx:782
LocalSolver::LengtheningBernoulliBendingBeam::hyEdge_dim
static constexpr unsigned int hyEdge_dim()
Return template parameter hyEdge_dimT.
Definition: bernoulli_beams.hxx:732
LocalSolver::LengtheningBernoulliBendingBeam::constructor_value_type
lSol_float_t constructor_value_type
Class is constructed using a single double indicating the penalty parameter.
Definition: bernoulli_beams.hxx:768
LocalSolver::LengtheningBernoulliBendingBeam::len_beam
const LengtheningBeam< hyEdge_dimT, space_dim, poly_deg, quad_deg, lSol_float_t > len_beam
The lengthening beam solver that does the lengthening of the beam.
Definition: bernoulli_beams.hxx:758
Hypercube::pow
static constexpr unsigned int pow(unsigned int base)
Return n to the power dim.
Definition: hypercube.hxx:25
LocalSolver::LengtheningBeam::HAS_MEMBER_FUNCTION
HAS_MEMBER_FUNCTION(trace_to_flux, has_trace_to_flux)
Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
LocalSolver::BernoulliBendingBeam::errors
std::array< lSol_float_t, 1U > errors(const std::array< std::array< lSol_float_t, n_glob_dofs_per_node()>, 2 *hyEdge_dimT > &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Local squared contribution to the L2 error.
Definition: bernoulli_beams.hxx:591
LocalSolver::BernoulliBendingBeam::hyEdge_dim
static constexpr unsigned int hyEdge_dim()
Return template parameter hyEdge_dimT.
Definition: bernoulli_beams.hxx:408
LocalSolver::LengtheningBeam::node_dof_to_edge_dof
std::array< std::array< double, diffusion_sol_t::n_glob_dofs_per_node()>, 2 *hyEdge_dimT > node_dof_to_edge_dof(const SmallMatT &lambda, hyEdgeT &hyper_edge) const
Do the pretprocessing to transfer global to local dofs.
Definition: bernoulli_beams.hxx:131
LocalSolver::LengtheningBeam::error_def
Define how errors are evaluated.
Definition: bernoulli_beams.hxx:60
hy_assert
#define hy_assert(Expr, Msg)
The assertion to be used within HyperHDG — deactivate using -DNDEBUG compile flag.
Definition: hy_assert.hxx:38
LocalSolver::BernoulliBendingBeam::HAS_MEMBER_FUNCTION
HAS_MEMBER_FUNCTION(trace_to_flux, has_trace_to_flux)
Prepare struct to check for function to exist (cf. compile_time_tricks.hxx).
dense_la.hxx
LocalSolver::LengtheningBeam::node_system_dimension
static constexpr unsigned int node_system_dimension()
Dimension of of the solution evaluated with respect to a hypernode.
Definition: bernoulli_beams.hxx:119
LocalSolver::BernoulliBendingBeam::system_dimension
static constexpr unsigned int system_dimension()
Dimension of of the solution evaluated with respect to a hyperedge.
Definition: bernoulli_beams.hxx:424
LocalSolver::LengtheningBeam::residual_flux
SmallMatOutT & residual_flux(const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Evaluate local contribution to residual.
Definition: bernoulli_beams.hxx:223
LocalSolver::LengtheningBeam::edge_dof_to_node_dof
SmallMatOutT & edge_dof_to_node_dof(const SmallMatInT &lambda, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge) const
Do the postprocessing to transfer local to global dofs.
Definition: bernoulli_beams.hxx:154
LocalSolver::BernoulliBendingBeam::error_def
Define how errors are evaluated.
Definition: bernoulli_beams.hxx:369
LocalSolver::LengtheningBeam::constructor_value_type
lSol_float_t constructor_value_type
Class is constructed using a single double indicating the penalty parameter.
Definition: bernoulli_beams.hxx:173
LocalSolver::BernoulliBendingBeam::constructor_value_type
lSol_float_t constructor_value_type
Class is constructed using a single double indicating the penalty parameter.
Definition: bernoulli_beams.hxx:494
LocalSolver::BernoulliBendingBeam::node_element
Define type of node elements, especially with respect to nodal shape functions.
Definition: bernoulli_beams.hxx:360
LocalSolver::LengtheningBernoulliBendingBeam::errors
std::array< lSol_float_t, 1U > errors(const std::array< std::array< lSol_float_t, n_glob_dofs_per_node()>, 2 *hyEdge_dimT > &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
Local squared contribution to the L2 error.
Definition: bernoulli_beams.hxx:820
SmallMat
This class implements a small/dense matrix.
Definition: dense_la.hxx:34
LocalSolver::LengtheningBernoulliBendingBeam::node_system_dimension
static constexpr unsigned int node_system_dimension()
Dimension of of the solution evaluated with respect to a hypernode.
Definition: bernoulli_beams.hxx:752
LocalSolver::LengtheningBeam::error_def::initial_error
static error_t initial_error()
Define how initial error is generated.
Definition: bernoulli_beams.hxx:69