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 "static_vector_input.hpp"
16 #include <iostream>
17 #include <math.h>
18 
19 using namespace StaticVectorSolver;
20 using namespace std;
21 
22 #pragma GCC diagnostic push
23 #pragma GCC diagnostic ignored "-Wunused-parameter"
24 #pragma GCC diagnostic ignored "-Wunused-but-set-variable"
25 
26 template<>
27 void
28 TheCoefficient<3>::value_list(const std::vector<Point<3>>& r,
29  types::material_id mid,
30  unsigned int cuid,
31  std::vector<double>& values) const
32 {
33  Assert(r.size() == values.size(),
34  ExcDimensionMismatch(r.size(), values.size()));
35 
36  auto v = values.begin();
37  for (auto p : r) {
38  *v = mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0);
39  v++;
40  }
41 }
42 
43 template<>
44 void
45 TheCoefficient<2>::value_list(const std::vector<Point<2>>& r,
46  types::material_id mid,
47  unsigned int cuid,
48  std::vector<double>& values) const
49 {
50  Assert(r.size() == values.size(),
51  ExcDimensionMismatch(r.size(), values.size()));
52 
53  auto v = values.begin();
54  for (auto p : r) {
55  *v = mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0);
56  v++;
57  }
58 }
59 
60 template<>
61 void
62 PdeRhs<3>::value_list(const std::vector<Point<3>>& r,
63  types::material_id mid,
64  unsigned int cuid,
65  std::vector<Tensor<1, 3>>& values) const
66 {
67  Assert(r.size() == values.size(),
68  ExcDimensionMismatch(r.size(), values.size()));
69 
70  auto v = values.begin();
71  for (auto p : r) {
72  (*v)[0] = 0.0;
73  (*v)[1] = 0.0;
74  (*v)[2] = (cos(k * p[0]) + cos(k * p[1])) /
75  ((mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0)));
76  v++;
77  }
78 }
79 
80 template<>
81 void
82 PdeRhs<2>::value_list(const std::vector<Point<2>>& r,
83  types::material_id mid,
84  unsigned int cuid,
85  std::vector<Tensor<1, 2>>& values) const
86 {
87  Assert(r.size() == values.size(),
88  ExcDimensionMismatch(r.size(), values.size()));
89 
90  auto v = values.begin();
91  for (auto p : r) {
92  (*v)[0] = (cos(k * p[0]) + cos(k * p[1])) /
93  ((mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0)));
94  (*v)[1] = 0.0;
95  v++;
96  }
97 }
98 
99 template<>
100 void
101 Gamma<3>::value_list(const std::vector<Point<3>>& r,
102  const std::vector<Tensor<1, 3>>& n,
103  types::boundary_id bid,
104  types::material_id mid,
105  unsigned int cuid,
106  unsigned int fuid,
107  std::vector<double>& values) const
108 {
109  Assert(r.size() == values.size(),
110  ExcDimensionMismatch(r.size(), values.size()));
111 
112  auto v = values.begin();
113  for (auto p : r) {
114  double mu = mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0);
115  *v = (1.0 / mu) * (sqrt(pow(p[0], 2) + pow(p[1], 2)) + 2.0);
116  v++;
117  }
118 }
119 
120 template<>
121 void
122 Gamma<2>::value_list(const std::vector<Point<2>>& r,
123  const std::vector<Tensor<1, 2>>& n,
124  types::boundary_id bid,
125  types::material_id mid,
126  unsigned int cuid,
127  unsigned int fuid,
128  std::vector<double>& values) const
129 {
130  Assert(r.size() == values.size(),
131  ExcDimensionMismatch(r.size(), values.size()));
132 
133  auto v = values.begin();
134  for (auto p : r) {
135  double mu = mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0);
136  *v = (1.0 / mu) * (sqrt(pow(p[0], 2) + pow(p[1], 2)) + 2.0);
137  v++;
138  }
139 }
140 
141 template<>
142 void
143 RobinRhs<3>::value_list(const std::vector<Point<3>>& r,
144  const std::vector<Tensor<1, 3>>& n,
145  types::boundary_id bid,
146  types::material_id mid,
147  unsigned int cuid,
148  unsigned int fuid,
149  std::vector<Tensor<1, 3>>& values) const
150 {
151  Assert(r.size() == values.size(),
152  ExcDimensionMismatch(r.size(), values.size()));
153 
154  double mu;
155  double gamma;
156  Tensor<1, 3> T;
157  Tensor<1, 3> A;
158  Tensor<1, 3> Q;
159 
160  auto v = values.begin();
161  auto nn = n.begin();
162  for (auto p : r) {
163  mu = mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0);
164  gamma = (1 / mu) * (sqrt(pow(p[0], 2) + pow(p[1], 2)) + 2.0);
165 
166  T[0] = 0.0;
167  T[1] = 0.0;
168  T[2] = (cos(k * p[0]) + cos(k * p[1])) / mu;
169 
170  A[0] = -sin(k * p[1]) / k;
171  A[1] = sin(k * p[0]) / k;
172  A[2] = 0.0;
173 
174  Q = cross_product_3d(*nn, T) +
175  gamma * cross_product_3d(*nn, cross_product_3d(*nn, A));
176 
177  (*v)[0] = Q[0];
178  (*v)[1] = Q[1];
179  (*v)[2] = Q[2];
180  v++;
181  nn++;
182  }
183 }
184 
185 template<>
186 void
187 RobinRhs<2>::value_list(const std::vector<Point<2>>& r,
188  const std::vector<Tensor<1, 2>>& n,
189  types::boundary_id bid,
190  types::material_id mid,
191  unsigned int cuid,
192  unsigned int fuid,
193  std::vector<Tensor<1, 2>>& values) const
194 {
195  Assert(r.size() == values.size(),
196  ExcDimensionMismatch(r.size(), values.size()));
197 
198  double mu;
199  double gamma;
200  double T;
201  Tensor<1, 2> A;
202  Tensor<1, 2> Q;
203 
204  auto v = values.begin();
205  auto nn = n.begin();
206  for (auto p : r) {
207  mu = mu_0 * (pow(p[0], 2) + pow(p[1], 2) + 1.0);
208  gamma = (1 / mu) * (sqrt(pow(p[0], 2) + pow(p[1], 2)) + 2.0);
209  T = (cos(k * p[0]) + cos(k * p[1])) / mu;
210 
211  A[0] = -sin(k * p[1]) / k;
212  A[1] = sin(k * p[0]) / k;
213 
214  Q[0] =
215  (*nn)[1] * T + gamma * (*nn)[1] * ((*nn)[0] * A[1] - (*nn)[1] * A[0]);
216  Q[1] =
217  -(*nn)[0] * T - gamma * (*nn)[0] * ((*nn)[0] * A[1] - (*nn)[1] * A[0]);
218 
219  (*v)[0] = Q[0];
220  (*v)[1] = Q[1];
221 
222  v++;
223  nn++;
224  }
225 }
226 
227 template<>
228 void
229 FreeSurfaceCurrent<3>::value_list(const std::vector<Point<3>>& r,
230  const std::vector<Tensor<1, 3>>& n,
231  types::material_id mid,
232  unsigned int cuid,
233  unsigned int fuid,
234  std::vector<Tensor<1, 3>>& values) const
235 {
236  Assert(r.size() == values.size(),
237  ExcDimensionMismatch(r.size(), values.size()));
238 
239  auto v = values.begin();
240  for (auto p : r) {
241  (*v)[0] = 0.0;
242  (*v)[1] = 0.0;
243  (*v)[2] = 0.0;
244  v++;
245  }
246 }
247 
248 template<>
249 void
250 FreeSurfaceCurrent<2>::value_list(const std::vector<Point<2>>& r,
251  const std::vector<Tensor<1, 2>>& n,
252  types::material_id mid,
253  unsigned int cuid,
254  unsigned int fuid,
255  std::vector<Tensor<1, 2>>& values) const
256 {
257  Assert(r.size() == values.size(),
258  ExcDimensionMismatch(r.size(), values.size()));
259 
260  auto v = values.begin();
261  for (auto p : r) {
262  (*v)[0] = 0.0;
263  (*v)[1] = 0.0;
264  v++;
265  }
266 }
267 
268 template<>
269 double
270 Weight<3>::value(const Point<3>& r, const unsigned int component) const
271 {
272  return 1.0;
273 }
274 
275 template<>
276 double
277 Weight<2>::value(const Point<2>& r, const unsigned int component) const
278 {
279  return 1.0;
280 }
281 
282 #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...