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 <math.h>
17 
18 #pragma GCC diagnostic push
19 #pragma GCC diagnostic ignored "-Wunused-parameter"
20 #pragma GCC diagnostic ignored "-Wunused-but-set-variable"
21 
22 template<>
23 void
25  const std::vector<Point<3>>& r,
26  types::material_id mid,
27  unsigned int cuid,
28  std::vector<double>& values) const
29 {
30  Assert(r.size() == values.size(),
31  ExcDimensionMismatch(r.size(), values.size()));
32 
33  for (unsigned int i = 0; i < values.size(); i++)
34  values[i] = 1.0;
35 }
36 
37 template<>
38 void
40  const std::vector<Point<3>>& r,
41  types::material_id mid,
42  unsigned int cuid,
43  std::vector<Tensor<1, 3>>& values) const
44 {
45  Assert(r.size() == values.size(),
46  ExcDimensionMismatch(r.size(), values.size()));
47 
48  auto v = values.begin();
49  for (auto p : r) {
50  if (mid == mid_2) {
51  (*v)[0] = -p[1];
52  (*v)[1] = p[0];
53  (*v)[2] = 0.0;
54  } else {
55  (*v)[0] = 0.0;
56  (*v)[1] = 0.0;
57  (*v)[2] = 0.0;
58  }
59  v++;
60  }
61 }
62 
63 template<>
64 void
65 StaticVectorSolver::Gamma<3>::value_list(const std::vector<Point<3>>& r,
66  const std::vector<Tensor<1, 3>>& n,
67  types::boundary_id bid,
68  types::material_id mid,
69  unsigned int cuid,
70  unsigned int fuid,
71  std::vector<double>& values) const
72 {
73  Assert(r.size() == values.size(),
74  ExcDimensionMismatch(r.size(), values.size()));
75 
76  for (unsigned int i = 0; i < values.size(); i++)
77  values[i] = 0.0;
78 }
79 
80 template<>
81 void
83  const std::vector<Point<3>>& r,
84  const std::vector<Tensor<1, 3>>& n,
85  types::boundary_id bid,
86  types::material_id mid,
87  unsigned int cuid,
88  unsigned int fuid,
89  std::vector<Tensor<1, 3>>& values) const
90 {
91  Assert(r.size() == values.size(),
92  ExcDimensionMismatch(r.size(), values.size()));
93 
94  for (unsigned int i = 0; i < values.size(); i++) {
95  values[i][0] = 0.0;
96  values[i][1] = 0.0;
97  values[i][2] = 0.0;
98  }
99 }
100 
101 template<>
102 void
104  const std::vector<Point<3>>& r,
105  const std::vector<Tensor<1, 3>>& n,
106  types::material_id mid,
107  unsigned int cuid,
108  unsigned int fuid,
109  std::vector<Tensor<1, 3>>& values) const
110 {
111  Assert(r.size() == values.size(),
112  ExcDimensionMismatch(r.size(), values.size()));
113 
114  for (unsigned int i = 0; i < values.size(); i++) {
115  values[i][0] = 0.0;
116  values[i][1] = 0.0;
117  values[i][2] = 0.0;
118  }
119 }
120 
121 template<>
122 double
123 StaticVectorSolver::Weight<3>::value(const Point<3>& r,
124  const unsigned int component) const
125 {
126  return 1.0;
127 }
128 
129 #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...