Logbook  (07-04-2025)
Static problems
static_vector_input.cpp
1 /******************************************************************************
2  * Copyright (C) Siarhei Uzunbajakau, 2023.
3  *
4  * This program is free software. You can use, modify, and redistribute it under
5  * the terms of the GNU Lesser General Public License as published by the Free
6  * Software Foundation, either version 3 or (at your option) any later version.
7  * This program is distributed without any warranty.
8  *
9  * Refer to COPYING.LESSER for more details.
10  ******************************************************************************/
11 
12 #include <deal.II/base/types.h>
13 #define BOOST_ALLOW_DEPRECATED_HEADERS
14 
15 #include "exact_solution.hpp"
16 #include "static_vector_input.hpp"
17 #include <math.h>
18 
19 #pragma GCC diagnostic push
20 #pragma GCC diagnostic ignored "-Wunused-parameter"
21 #pragma GCC diagnostic ignored "-Wunused-but-set-variable"
22 
23 //-----------------------------------------------------------------------------
24 //-------- Stage 0. Calculating current vector potential, T, given Jf ---------
25 //-----------------------------------------------------------------------------
26 
27 template<>
28 void
30  const std::vector<Point<3>>& r,
31  types::material_id mid,
32  unsigned int cuid,
33  std::vector<double>& values) const
34 {
35  Assert(r.size() == values.size(),
36  ExcDimensionMismatch(r.size(), values.size()));
37 
38  for (unsigned int i = 0; i < values.size(); i++)
39  values[i] = 1.0;
40 }
41 
42 template<>
43 void
45  const std::vector<Point<3>>& r,
46  types::material_id mid,
47  unsigned int cuid,
48  std::vector<Tensor<1, 3>>& values) const
49 {
50  Assert(r.size() == values.size(),
51  ExcDimensionMismatch(r.size(), values.size()));
52 
53  Tensor<1, 3> Jf;
54 
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);
57 
58  values[i][0] = Jf[0];
59  values[i][1] = Jf[1];
60  values[i][2] = Jf[2];
61  }
62 }
63 
64 template<>
65 void
66 StaticVectorSolver::Gamma<3, 0>::value_list(const std::vector<Point<3>>& r,
67  const std::vector<Tensor<1, 3>>& n,
68  types::boundary_id bid,
69  types::material_id mid,
70  unsigned int cuid,
71  unsigned int fuid,
72  std::vector<double>& values) const
73 {
74  Assert(r.size() == values.size(),
75  ExcDimensionMismatch(r.size(), values.size()));
76 
77  Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
78 
79  for (unsigned int i = 0; i < values.size(); i++)
80  values[i] = 0;
81 }
82 
83 template<>
84 void
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,
90  unsigned int cuid,
91  unsigned int fuid,
92  std::vector<Tensor<1, 3>>& values) const
93 {
94  Assert(r.size() == values.size(),
95  ExcDimensionMismatch(r.size(), values.size()));
96 
97  Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
98 
99  Tensor<1, 3> Jf;
100 
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);
103 
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]);
107  }
108 }
109 
110 template<>
111 void
113  const std::vector<Point<3>>& r,
114  const std::vector<Tensor<1, 3>>& n,
115  types::material_id mid,
116  unsigned int cuid,
117  unsigned int fuid,
118  std::vector<Tensor<1, 3>>& values) const
119 {
120  Assert(r.size() == values.size(),
121  ExcDimensionMismatch(r.size(), values.size()));
122 
123  Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
124 
125  for (unsigned int i = 0; i < values.size(); i++) {
126  values[i][0] = 0.0;
127  values[i][1] = 0.0;
128  values[i][2] = 0.0;
129  }
130 }
131 
132 template<>
133 double
135  const unsigned int component) const
136 {
137  return 1.0;
138 }
139 
140 //-----------------------------------------------------------------------------
141 // Stage 1. Calculating Jf given T (projection) to check if the Stage 0 was ok.
142 //-----------------------------------------------------------------------------
143 
144 template<>
145 double
147  const unsigned int component) const
148 {
149  return 1.0;
150 }
151 
152 //-----------------------------------------------------------------------------
153 // Stage 2. Calculating vector potential, A, given T --------------------------
154 //-----------------------------------------------------------------------------
155 
156 template<>
157 void
159  const std::vector<Point<3>>& r,
160  types::material_id mid,
161  unsigned int cuid,
162  std::vector<double>& values) const
163 {
164  Assert(r.size() == values.size(),
165  ExcDimensionMismatch(r.size(), values.size()));
166 
167  for (unsigned int i = 0; i < values.size(); i++)
168  values[i] = permeability(r[i][0], r[i][1], mu_0);
169 }
170 
171 template<>
172 void
174  const std::vector<Point<3>>& r,
175  types::material_id mid,
176  unsigned int cuid,
177  std::vector<Tensor<1, 3>>& values) const
178 {
179  Assert(r.size() == values.size(),
180  ExcDimensionMismatch(r.size(), values.size()));
181 
182  for (unsigned int i = 0; i < values.size(); i++) {
183  values[i][0] = 0.0;
184  values[i][1] = 0.0;
185  values[i][2] = 0.0;
186  }
187 }
188 
189 template<>
190 void
191 StaticVectorSolver::Gamma<3, 2>::value_list(const std::vector<Point<3>>& r,
192  const std::vector<Tensor<1, 3>>& n,
193  types::boundary_id bid,
194  types::material_id mid,
195  unsigned int cuid,
196  unsigned int fuid,
197  std::vector<double>& values) const
198 {
199  Assert(r.size() == values.size(),
200  ExcDimensionMismatch(r.size(), values.size()));
201 
202  Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
203 
204  for (unsigned int i = 0; i < values.size(); i++)
205  values[i] = robin_gamma(r[i][0], r[i][1], mu_0);
206 }
207 
208 template<>
209 void
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,
215  unsigned int cuid,
216  unsigned int fuid,
217  std::vector<Tensor<1, 3>>& values) const
218 {
219  Assert(r.size() == values.size(),
220  ExcDimensionMismatch(r.size(), values.size()));
221 
222  Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
223 
224  Tensor<1, 3> T, A;
225  double mu, gamma;
226 
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);
232 
233  values[i] = cross_product_3d(n[i], T) +
234  gamma * cross_product_3d(n[i], cross_product_3d(n[i], A));
235  }
236 }
237 
238 template<>
239 void
241  const std::vector<Point<3>>& r,
242  const std::vector<Tensor<1, 3>>& n,
243  types::material_id mid,
244  unsigned int cuid,
245  unsigned int fuid,
246  std::vector<Tensor<1, 3>>& values) const
247 {
248  Assert(r.size() == values.size(),
249  ExcDimensionMismatch(r.size(), values.size()));
250 
251  Assert(r.size() == n.size(), ExcDimensionMismatch(r.size(), n.size()));
252 
253  for (unsigned int i = 0; i < values.size(); i++) {
254  values[i][0] = 0.0;
255  values[i][1] = 0.0;
256  values[i][2] = 0.0;
257  }
258 }
259 
260 template<>
261 double
263  const unsigned int component) const
264 {
265  return 1.0;
266 }
267 
268 //-----------------------------------------------------------------------------
269 //-------- Stage 3. Calculating B given A -------------------------------------
270 //-----------------------------------------------------------------------------
271 
272 template<>
273 double
275  const unsigned int component) const
276 {
277  return 1.0;
278 }
279 
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...