Logbook  (07-04-2025)
Static problems
solver.hpp
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 #ifndef SolverSLDI_H__
13 #define SolverSLDI_H__
14 
15 #define BOOST_ALLOW_DEPRECATED_HEADERS
16 
17 #include <deal.II/base/function.h>
18 #include <deal.II/grid/grid_in.h>
19 #include <deal.II/grid/grid_out.h>
20 #include <deal.II/grid/grid_tools.h>
21 #include <deal.II/grid/manifold_lib.h>
22 
23 #include <deal.II/numerics/fe_field_function.h>
24 
25 #include "exact_solution.hpp"
26 #include "settings.hpp"
27 #include "static_scalar_solver.hpp"
28 
29 #define TMR(__name) TimerOutput::Scope timer_section(timer, __name)
30 
31 using namespace StaticScalarSolver;
32 
33 #pragma GCC diagnostic push
34 #pragma GCC diagnostic ignored "-Wunused-parameter"
41 template<int dim>
43  : public SettingsSLDI
44  , public Function<dim>
45 {
46 public:
47  DirichletSLDI(){};
48 
49  virtual double value(const Point<dim>& r,
50  const unsigned int component = 0) const override final;
51 };
52 
53 #pragma GCC diagnostic pop
54 
60 template<int dim>
62  : public SettingsSLDI
63  , public Solver<dim>
64 {
65 public:
66  SolverSLDI() = delete;
82  SolverSLDI(unsigned int p,
83  unsigned int mapping_degree,
84  unsigned int r,
85  std::string fname)
86  : Solver<dim>(p,
87  mapping_degree,
88  0,
89  fname,
90  &exact_solution,
91  false,
92  false,
93  print_time_tables,
94  project_exact_solution,
95  true)
96  , r(r)
97  , fname(fname)
98  {
100  }
101 
102  ~SolverSLDI() = default;
103 
104 private:
105  const unsigned int r;
106  const std::string fname;
107 
108  const ExactSolutionSLDI_PSI<dim> exact_solution;
109  const DirichletSLDI<dim> dirichlet;
110 
111  SphericalManifold<dim> sphere;
112 
113  virtual void make_mesh() override final;
114  virtual void fill_dirichlet_stack() override final;
115  virtual void solve() override final;
116 };
117 
118 template<int dim>
119 void
120 SolverSLDI<dim>::fill_dirichlet_stack()
121 {
122  switch (type_of_bc) {
123  case DirichletOnly:
124  Solver<dim>::dirichlet_stack = { { bid, &dirichlet } };
125  break;
126  case DirichletNeumann:
127  Solver<dim>::dirichlet_stack = { { bid, &dirichlet } };
128  break;
129  case Exact:
130  Solver<dim>::dirichlet_stack = { { bid, &exact_solution } };
131  break;
132  }
133 }
134 
135 template<int dim>
136 void
138 {
139  GridIn<dim> gridin;
140  gridin.attach_triangulation(Solver<dim>::triangulation);
141 
142  std::string fname_mesh =
143  (dim == 2) ? "../../gmsh/data/square_r" + std::to_string(r) + ".msh"
144  : "../../gmsh/data/cube_r" + std::to_string(r) + ".msh";
145 
146  std::ifstream ifs(fname_mesh);
147  gridin.read_msh(ifs);
148 
149  Solver<dim>::triangulation.reset_all_manifolds();
150 
151  for (auto cell : Solver<dim>::triangulation.active_cell_iterators()) {
152 
153  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; f++) {
154  if (((type_of_bc == DirichletOnly) || (type_of_bc == Exact)) &&
155  cell->face(f)->at_boundary())
156  cell->face(f)->set_boundary_id(bid);
157  }
158 
159  if (cell->center().norm() < a) {
160  // The cell is inside the shield.
161  cell->set_material_id(mid_1);
162  } else if (cell->center().norm() < b) {
163  // The cell is inside the wall of the shield.
164  cell->set_material_id(mid_2);
165  } else {
166  // The cell is outside the shield.
167  cell->set_material_id(mid_3);
168  }
169 
170  if ((cell->center().norm() > rd1) && (cell->center().norm() < b))
171  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; f++) {
172  double dif_norm = 0.0;
173  for (unsigned int v = 1; v < GeometryInfo<dim>::vertices_per_face; v++)
174  dif_norm += std::abs(cell->face(f)->vertex(0).norm() -
175  cell->face(f)->vertex(v).norm());
176 
177  if (dif_norm < eps)
178  cell->face(f)->set_all_manifold_ids(1);
179  }
180  }
181 
182  Solver<dim>::triangulation.set_manifold(1, sphere);
183 }
184 
185 template<int dim>
186 void
188 {
189  SolverControl control(Solver<dim>::system_rhs.size(),
190  1e-8 * Solver<dim>::system_rhs.l2_norm(),
191  false,
192  false);
193 
194  if (log_cg_convergence)
195  control.enable_history_data();
196 
197  GrowingVectorMemory<Vector<double>> memory;
198  SolverCG<Vector<double>> cg(control, memory);
199 
200  PreconditionJacobi<SparseMatrix<double>> preconditioner;
201  preconditioner.initialize(Solver<dim>::system_matrix, 1.0);
202 
206  preconditioner);
207 
209 
210  if (log_cg_convergence) {
211  const std::vector<double> history_data = control.get_history_data();
212 
213  std::ofstream ofs(fname + "_cg_convergence.csv");
214 
215  unsigned int i = 1;
216  for (auto item : history_data) {
217  ofs << i << ", " << item << "\n";
218  i++;
219  }
220  ofs.close();
221  }
222 }
223 
224 #endif
The Dirichlet boundary condition of the Magnetostatic shield - 1 (sld-i/) numerical experiment....
Definition: solver.hpp:45
Describes exact solution, , of the Magnetostatic shield - 1 (sld-i/) numerical experiment in two and ...
Global settings for the Magnetostatic shield - 1 (sld-i/) numerical experiment.
Definition: settings.hpp:33
Implements the solver that solves for in the Magnetostatic shield - 1 (sld-i/) numerical experiment.
Definition: solver.hpp:64
SolverSLDI(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Definition: solver.hpp:82
Solves static scalar boundary value problem.
void run()
Runs the simulation.