12 #include <deal.II/base/types.h>
13 #define BOOST_ALLOW_DEPRECATED_HEADERS
15 #include "exact_solution.hpp"
16 #include "static_vector_input.hpp"
19 #pragma GCC diagnostic push
20 #pragma GCC diagnostic ignored "-Wunused-parameter"
21 #pragma GCC diagnostic ignored "-Wunused-but-set-variable"
30 const std::vector<Point<3>>& r,
31 types::material_id mid,
33 std::vector<double>& values)
const
35 Assert(r.size() == values.size(),
36 ExcDimensionMismatch(r.size(), values.size()));
38 for (
unsigned int i = 0; i < values.size(); i++)
45 const std::vector<Point<3>>& r,
46 types::material_id mid,
48 std::vector<Tensor<1, 3>>& values)
const
50 Assert(r.size() == values.size(),
51 ExcDimensionMismatch(r.size(), values.size()));
55 for (
unsigned int i = 0; i < values.size(); i++) {
56 Jf = volume_free_current_density(r[i][0], r[i][1], mu_0, k);
67 const std::vector<Tensor<1, 3>>& n,
68 types::boundary_id bid,
69 types::material_id mid,
72 std::vector<double>& values)
const
74 Assert(r.size() == values.size(),
75 ExcDimensionMismatch(r.size(), values.size()));
77 Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
79 for (
unsigned int i = 0; i < values.size(); i++)
86 const std::vector<Point<3>>& r,
87 const std::vector<Tensor<1, 3>>& n,
88 types::boundary_id bid,
89 types::material_id mid,
92 std::vector<Tensor<1, 3>>& values)
const
94 Assert(r.size() == values.size(),
95 ExcDimensionMismatch(r.size(), values.size()));
97 Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
101 for (
unsigned int i = 0; i < r.size(); i++) {
102 Jf = volume_free_current_density(r[i][0], r[i][1], mu_0, k);
104 values[i][0] = (n[i][1] * Jf[2] - n[i][2] * Jf[1]);
105 values[i][1] = -(n[i][0] * Jf[2] - n[i][2] * Jf[0]);
106 values[i][2] = (n[i][0] * Jf[1] - n[i][1] * Jf[0]);
113 const std::vector<Point<3>>& r,
114 const std::vector<Tensor<1, 3>>& n,
115 types::material_id mid,
118 std::vector<Tensor<1, 3>>& values)
const
120 Assert(r.size() == values.size(),
121 ExcDimensionMismatch(r.size(), values.size()));
123 Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
125 for (
unsigned int i = 0; i < values.size(); i++) {
135 const unsigned int component)
const
147 const unsigned int component)
const
159 const std::vector<Point<3>>& r,
160 types::material_id mid,
162 std::vector<double>& values)
const
164 Assert(r.size() == values.size(),
165 ExcDimensionMismatch(r.size(), values.size()));
167 for (
unsigned int i = 0; i < values.size(); i++)
168 values[i] = permeability(r[i][0], r[i][1], mu_0);
174 const std::vector<Point<3>>& r,
175 types::material_id mid,
177 std::vector<Tensor<1, 3>>& values)
const
179 Assert(r.size() == values.size(),
180 ExcDimensionMismatch(r.size(), values.size()));
182 for (
unsigned int i = 0; i < values.size(); i++) {
192 const std::vector<Tensor<1, 3>>& n,
193 types::boundary_id bid,
194 types::material_id mid,
197 std::vector<double>& values)
const
199 Assert(r.size() == values.size(),
200 ExcDimensionMismatch(r.size(), values.size()));
202 Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
204 for (
unsigned int i = 0; i < values.size(); i++)
205 values[i] = robin_gamma(r[i][0], r[i][1], mu_0);
211 const std::vector<Point<3>>& r,
212 const std::vector<Tensor<1, 3>>& n,
213 types::boundary_id bid,
214 types::material_id mid,
217 std::vector<Tensor<1, 3>>& values)
const
219 Assert(r.size() == values.size(),
220 ExcDimensionMismatch(r.size(), values.size()));
222 Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
227 for (
unsigned int i = 0; i < r.size(); i++) {
228 mu = permeability(r[i][0], r[i][1], mu_0);
229 gamma = robin_gamma(r[i][0], r[i][1], mu_0);
230 T = current_vector_potential(r[i][0], r[i][1], mu_0, k);
231 A = magnetic_vector_potential(r[i][0], r[i][1], k);
233 values[i] = cross_product_3d(n[i], T) +
234 gamma * cross_product_3d(n[i], cross_product_3d(n[i], A));
241 const std::vector<Point<3>>& r,
242 const std::vector<Tensor<1, 3>>& n,
243 types::material_id mid,
246 std::vector<Tensor<1, 3>>& values)
const
248 Assert(r.size() == values.size(),
249 ExcDimensionMismatch(r.size(), values.size()));
251 Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
253 for (
unsigned int i = 0; i < values.size(); i++) {
263 const unsigned int component)
const
275 const unsigned int component)
const
280 #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< Tensor< 1, dim >> &values) const
Computes values of the surface free-current density, , on the right-hand side of the continuity condi...
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 values of the parameter of the Robin boundary condition 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 vector field on the right-hand side of the partial differential equation.
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< Tensor< 1, dim >> &values) const
Computes values of the vector field on the right-hand side of the Robin boundary condition at quadra...
void value_list(const std::vector< Point< dim >> &r, types::material_id mid, unsigned int cuid, std::vector< double > &values) const
Computes values of permeability, , 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...