Logbook  (07-04-2025)
Static problems
solver.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 "solver.hpp"
13 
14 void
15 SolverSSOLIIIAXI::make_mesh()
16 {
17  GridIn<2> gridin;
18  gridin.attach_triangulation(triangulation);
19  std::ifstream ifs("../../gmsh/data/circle_r" + std::to_string(r) + ".msh");
20  gridin.read_msh(ifs);
21  mark_materials();
22 }
23 
24 void
25 SolverSSOLIIIAXI::fill_dirichlet_stack()
26 {
27  Solver<2>::dirichlet_stack = { { bid_axi, &dirichlet_function } };
28 }
29 
30 void
31 SolverSSOLIIIAXI::mark_materials()
32 {
33  Solver<2>::triangulation.reset_all_manifolds();
34 
35  for (auto cell : Solver<2>::triangulation.active_cell_iterators()) {
36  cell->set_material_id(mid_1);
37 
38  if ((cell->center().norm() < b1) && (cell->center().norm() > a1))
39  cell->set_material_id(mid_2);
40 
41  if ((cell->center().norm() < b2) && (cell->center().norm() > a2))
42  cell->set_material_id(mid_3);
43 
44  if (cell->center().norm() > rd1)
45  for (unsigned int f = 0; f < GeometryInfo<2>::faces_per_cell; f++)
46  if (std::abs(cell->face(f)->vertex(1).norm() -
47  cell->face(f)->vertex(0).norm()) < eps)
48  cell->face(f)->set_all_manifold_ids(1);
49  }
50 
51  Solver<2>::triangulation.set_manifold(1, sphere);
52 }
53 
54 void
55 SolverSSOLIIIAXI::data_slice(std::string fname)
56 {
57  GridGenerator::hyper_cube(triangulation_slice, 0.0 + eps, d3 - eps);
58  triangulation_slice.refine_global(nr_slice_global_refs);
59 
60  dof_handler_slice.reinit(triangulation_slice);
61  dof_handler_slice.distribute_dofs(fe_slice);
62  solution_slice.reinit(dof_handler_slice.n_dofs());
63 
64  Functions::FEFieldFunction<2> potential(Solver<2>::dof_handler,
66 
67  VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
68 
69  DataOut<1, 2> data_out;
70 
71  data_out.attach_dof_handler(dof_handler_slice);
72  data_out.add_data_vector(solution_slice, "solution_slice");
73  data_out.build_patches();
74 
75  std::ofstream out(fname + "_slice" + ".gpi");
76 
77  data_out.write_gnuplot(out);
78  out.close();
79 }
80 
81 void
82 SolverSSOLIIIAXI::solve()
83 {
84  SolverControl control(Solver<2>::system_rhs.size(),
85  1e-8 * Solver<2>::system_rhs.l2_norm(),
86  false,
87  false);
88 
90  control.enable_history_data();
91 
92  GrowingVectorMemory<Vector<double>> memory;
93  SolverCG<Vector<double>> cg(control, memory);
94 
95  PreconditionJacobi<SparseMatrix<double>> preconditioner;
96  preconditioner.initialize(Solver<2>::system_matrix, 1.0);
97 
98  cg.solve(Solver<2>::system_matrix,
101  preconditioner);
102 
104 
105  if (log_cg_convergence) {
106  const std::vector<double> history_data = control.get_history_data();
107 
108  std::ofstream ofs(fname + "_cg_convergence.csv");
109 
110  unsigned int i = 1;
111  for (auto item : history_data) {
112  ofs << i << ", " << item << "\n";
113  i++;
114  }
115  ofs.close();
116  }
117 }
const double eps
Two values in double format are considered to be equal if the absolute value of their difference is l...
Definition: settings.hpp:135
const types::material_id mid_2
The ID of the material inside the magnetic core.
Definition: settings.hpp:112
const types::material_id mid_3
The ID of the material inside the coil, i.e., the windings.
Definition: settings.hpp:117
const double b2
The outer radius of the coil.
Definition: settings.hpp:86
const double b1
The outer radius of the magnetic core.
Definition: settings.hpp:76
const types::material_id mid_1
The ID of the material outside the coil and the magnetic core.
Definition: settings.hpp:107
const double a1
The inner radius of the magnetic core.
Definition: settings.hpp:71
const double rd1
The radius of the circle that encloses the rectangle in the middle of the mesh.
Definition: settings.hpp:66
const double d3
The radius of the problem domain.
Definition: settings.hpp:97
const double a2
The inner radius of the coil.
Definition: settings.hpp:81
const bool log_cg_convergence
If set to true, saves the residual at each iteration of the CG solver. The names of the files fit the...
Definition: settings.hpp:160
const types::boundary_id bid_axi
The ID of the straight segment of the boundary of the problem domain.
Definition: settings.hpp:129
std::map< types::boundary_id, const Function< dim > * > dirichlet_stack
A map that contains pairs of boundary IDs and the corresponding Dirichlet boundary conditions.
Triangulation< dim > triangulation
The mesh.
DoFHandler< dim > dof_handler
The degrees-of-freedom handler.
Vector< double > system_rhs
The system right-hand side vector.
Vector< double > solution
The solution vector, i.e., degrees of freedom yielded by the simulation.
AffineConstraints< double > constraints
The constraints associated with the Dirichlet boundary conditions.
SparseMatrix< double > system_matrix
The system matrix.