12 #include <deal.II/base/types.h>
13 #define BOOST_ALLOW_DEPRECATED_HEADERS
15 #include "static_scalar_input.hpp"
18 using namespace StaticScalarSolver;
21 #pragma GCC diagnostic push
22 #pragma GCC diagnostic ignored "-Wunused-parameter"
27 types::material_id mid,
29 std::vector<double>& values)
const
31 Assert(r.size() == values.size(),
32 ExcDimensionMismatch(r.size(), values.size()));
34 for (
unsigned int i = 0; i < values.size(); i++)
35 values[i] = 1.0 / (mu_0 * r.at(i)[0]);
41 types::material_id mid,
43 std::vector<double>& values)
const
45 Assert(r.size() == values.size(),
46 ExcDimensionMismatch(r.size(), values.size()));
48 for (
unsigned int i = 0; i < values.size(); i++)
55 types::material_id mid,
57 std::vector<Tensor<1, 2>>& values)
const
59 Assert(r.size() == values.size(),
60 ExcDimensionMismatch(r.size(), values.size()));
62 for (
unsigned int i = 0; i < values.size(); i++) {
63 values.at(i)[0] = 0.0;
64 values.at(i)[1] = 0.0;
71 const std::vector<Tensor<1, 2>>& n,
72 types::boundary_id bid,
73 types::material_id mid,
76 std::vector<double>& values)
const
78 Assert(r.size() == values.size(),
79 ExcDimensionMismatch(r.size(), values.size()));
81 for (
unsigned int i = 0; i < values.size(); i++)
88 const std::vector<Tensor<1, 2>>& n,
89 types::boundary_id bid,
90 types::material_id mid,
93 std::vector<double>& values)
const
96 Assert(r.size() == values.size(),
97 ExcDimensionMismatch(r.size(), values.size()));
99 for (
unsigned int i = 0; i < values.size(); i++)
106 const std::vector<Tensor<1, 2>>& n,
107 types::material_id mid,
110 std::vector<double>& values)
const
112 Assert(r.size() == values.size(),
113 ExcDimensionMismatch(r.size(), values.size()));
115 if ((cuid == 1) && (fuid == 1)) {
116 for (
unsigned int i = 0; i < values.size(); i++) {
117 values[i] = K_0 * r.at(i)[0];
120 for (
unsigned int i = 0; i < values.size(); i++)
132 #pragma GCC diagnostic pop
void value_list(const std::vector< Point< dim >> &r, const std::vector< Tensor< 1, dim >> &n, types::material_id mid, unsigned int cuid, unsigned int fuid, std::vector< double > &values) const
Computes the right-hand side of the second continuity condition ( , , , or ).
void value_list(const std::vector< Point< dim >> &r, const std::vector< Tensor< 1, dim >> &n, types::boundary_id bid, types::material_id mid, unsigned int cuid, unsigned int fuid, std::vector< double > &values) const
Computes the coefficient at quadrature points.
void value_list(const std::vector< Point< dim >> &r, types::material_id mid, unsigned int cuid, std::vector< Tensor< 1, dim >> &values) const
Computes the two-dimensional free-current density on the right-hand side of the partial differential...
void value_list(const std::vector< Point< dim >> &r, types::material_id mid, unsigned int cuid, std::vector< double > &values) const
Computes the right-hand side of the div-grad partial differential equation at quadrature points.
void value_list(const std::vector< Point< dim >> &r, const std::vector< Tensor< 1, dim >> &n, types::boundary_id bid, types::material_id mid, unsigned int cuid, unsigned int fuid, std::vector< double > &values) const
Computes the right-hand side of the Robin boundary condition ( or ).
void value_list(const std::vector< Point< dim >> &r, types::material_id mid, unsigned int cuid, std::vector< double > &values) const
Computes the values of the coefficient at quadrature points.
virtual double value(const Point< dim > &r, const unsigned int component=0) const override final
Returns the value of weight at point r. All error norms, , , and , at point r will be multiplied by t...