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 SolverSSOLIAXI::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 SolverSSOLIAXI::fill_dirichlet_stack()
26 {
27  switch (type_of_bc) {
28  case Dirichlet:
29  Solver<2>::dirichlet_stack = { { bid_infty, &dirichlet_function },
30  { bid_axi, &dirichlet_function } };
31  break;
32  case Exact:
33  Solver<2>::dirichlet_stack = { { bid_infty, &exact_solution },
34  { bid_axi, &dirichlet_function } };
35  break;
36  }
37 }
38 
39 void
40 SolverSSOLIAXI::mark_materials()
41 {
42  Solver<2>::triangulation.reset_all_manifolds();
43 
44  for (auto cell : Solver<2>::triangulation.active_cell_iterators()) {
45  for (unsigned int f = 0; f < GeometryInfo<2>::faces_per_cell; f++) {
46  if (std::abs(cell->center().norm()) < a) {
47  if ((abs(cell->face(f)->vertex(0).norm() - a) < eps) &&
48  (abs(cell->face(f)->vertex(1).norm() - a) < eps)) {
49  cell->face(f)->set_user_index(1);
50  cell->set_user_index(1);
51  }
52  }
53 
54  if (cell->center().norm() > rd)
55  if (std::abs(cell->face(f)->vertex(1).norm() -
56  cell->face(f)->vertex(0).norm()) < eps)
57  cell->face(f)->set_all_manifold_ids(1);
58  }
59  }
60 
61  Solver<2>::triangulation.set_manifold(1, sphere);
62 }
63 
64 void
65 SolverSSOLIAXI::data_slice(std::string fname)
66 {
67  GridGenerator::hyper_cube(triangulation_slice, 0.0 + eps, b - eps);
68  triangulation_slice.refine_global(nr_slice_global_refs);
69 
70  dof_handler_slice.reinit(triangulation_slice);
71  dof_handler_slice.distribute_dofs(fe_slice);
72  solution_slice.reinit(dof_handler_slice.n_dofs());
73 
74  Functions::FEFieldFunction<2> potential(Solver<2>::dof_handler,
76 
77  VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
78 
79  DataOut<1, 2> data_out;
80 
81  data_out.attach_dof_handler(dof_handler_slice);
82  data_out.add_data_vector(solution_slice, "solution_slice");
83  data_out.build_patches();
84 
85  std::ofstream out(fname + "_slice" + ".gpi");
86 
87  data_out.write_gnuplot(out);
88  out.close();
89 }
90 
91 void
92 SolverSSOLIAXI::solve()
93 {
94  SolverControl control(Solver<2>::system_rhs.size(),
95  1e-8 * Solver<2>::system_rhs.l2_norm(),
96  false,
97  false);
98 
100  control.enable_history_data();
101 
102  GrowingVectorMemory<Vector<double>> memory;
103  SolverCG<Vector<double>> cg(control, memory);
104 
105  PreconditionJacobi<SparseMatrix<double>> preconditioner;
106  preconditioner.initialize(Solver<2>::system_matrix, 1.0);
107 
108  cg.solve(Solver<2>::system_matrix,
111  preconditioner);
112 
114 
115  if (log_cg_convergence) {
116  const std::vector<double> history_data = control.get_history_data();
117 
118  std::ofstream ofs(fname + "_cg_convergence.csv");
119 
120  unsigned int i = 1;
121  for (auto item : history_data) {
122  ofs << i << ", " << item << "\n";
123  i++;
124  }
125  ofs.close();
126  }
127 }
double b
The radius of the outer boundary of the problem domain.
Definition: settings.hpp:72
const types::boundary_id bid_infty
The ID of the curved section of the boundary of the problem domain. The boundary ID is set in the geo...
Definition: settings.hpp:79
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:86
const BoundaryConditionType type_of_bc
Switches between two boundary conditions options.
Definition: settings.hpp:91
double a
The radius of the coil.
Definition: settings.hpp:67
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:122
const double rd
The radius of the circle that encloses the rectangle in the middle of the mesh.
Definition: settings.hpp:62
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:97
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.