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 SolverSSOLIIAXI::make_mesh()
16 {
17  GridIn<2> gridin;
18  gridin.attach_triangulation(triangulation);
19 
20  std::ifstream ifs("../../gmsh/data/circle_r" + std::to_string(r) + ".msh");
21  gridin.read_msh(ifs);
22 
23  mark_materials();
24 }
25 
26 void
27 SolverSSOLIIAXI::fill_dirichlet_stack()
28 {
29  Solver<2>::dirichlet_stack = { { bid_axi, &dirichlet_function } };
30 }
31 
32 void
33 SolverSSOLIIAXI::mark_materials()
34 {
35  Solver<2>::triangulation.reset_all_manifolds();
36 
37  for (auto cell : Solver<2>::triangulation.active_cell_iterators()) {
38  cell->set_material_id(mid_1);
39 
40  if ((cell->center().norm() < b) && (cell->center().norm() > a))
41  cell->set_material_id(mid_2);
42 
43  if (cell->center().norm() > rd1)
44  for (unsigned int f = 0; f < GeometryInfo<2>::faces_per_cell; f++)
45  if (std::abs(cell->face(f)->vertex(1).norm() -
46  cell->face(f)->vertex(0).norm()) < eps)
47  cell->face(f)->set_all_manifold_ids(1);
48  }
49 
50  Solver<2>::triangulation.set_manifold(1, sphere);
51 }
52 
53 void
54 SolverSSOLIIAXI::data_slice(std::string fname)
55 {
56  GridGenerator::hyper_cube(triangulation_slice, 0.0 + eps, d3 - eps);
57  triangulation_slice.refine_global(nr_slice_global_refs);
58 
59  dof_handler_slice.reinit(triangulation_slice);
60  dof_handler_slice.distribute_dofs(fe_slice);
61  solution_slice.reinit(dof_handler_slice.n_dofs());
62 
63  Functions::FEFieldFunction<2> potential(Solver<2>::dof_handler,
65 
66  VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
67 
68  DataOut<1, 2> data_out;
69 
70  data_out.attach_dof_handler(dof_handler_slice);
71  data_out.add_data_vector(solution_slice, "solution_slice");
72  data_out.build_patches();
73 
74  std::ofstream out(fname + "_slice" + ".gpi");
75 
76  data_out.write_gnuplot(out);
77  out.close();
78 }
79 
80 void
81 SolverSSOLIIAXI::solve()
82 {
83  SolverControl control(Solver<2>::system_rhs.size(),
84  1e-8 * Solver<2>::system_rhs.l2_norm(),
85  false,
86  false);
87 
89  control.enable_history_data();
90 
91  GrowingVectorMemory<Vector<double>> memory;
92  SolverCG<Vector<double>> cg(control, memory);
93 
94  PreconditionJacobi<SparseMatrix<double>> preconditioner;
95  preconditioner.initialize(Solver<2>::system_matrix, 1.0);
96 
97  cg.solve(Solver<2>::system_matrix,
100  preconditioner);
101 
103 
104  if (log_cg_convergence) {
105  const std::vector<double> history_data = control.get_history_data();
106 
107  std::ofstream ofs(fname + "_cg_convergence.csv");
108 
109  unsigned int i = 1;
110  for (auto item : history_data) {
111  ofs << i << ", " << item << "\n";
112  i++;
113  }
114  ofs.close();
115  }
116 }
double a
The inner radius of the coil.
Definition: settings.hpp:61
const double rd1
The radius of the circle that encloses the rectangle in the middle of the mesh.
Definition: settings.hpp:56
double b
The outer radius of the coil.
Definition: settings.hpp:66
double d3
The radius of the problem domain.
Definition: settings.hpp:76
const types::material_id mid_1
The ID of the material outside the coil, J_f is zero in this region.
Definition: settings.hpp:81
const types::boundary_id bid_axi
The ID of the straight section of the boundary of the problem domain. The boundary ID is set in the g...
Definition: settings.hpp:101
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:107
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:132
const types::material_id mid_2
The ID of the material inside the coil, J_f is nonzero in this region.
Definition: settings.hpp:87
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.