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 #define BOOST_ALLOW_DEPRECATED_HEADERS
13 
14 #include "solver.hpp"
15 
16 void
17 SolverMMSVTII_T::make_mesh()
18 {
19  GridIn<2> gridin;
20 
21  gridin.attach_triangulation(Solver<2, 0>::triangulation);
22 #if DOMAIN__ == 0
23  std::ifstream ifs("../../gmsh/data/square_r" + std::to_string(r) + ".msh");
24 #elif DOMAIN__ == 1
25  std::ifstream ifs("../../gmsh/data/circle_r" + std::to_string(r) + ".msh");
26 #endif
27  gridin.read_msh(ifs);
28 
29 #if DOMAIN__ == 1
30  Solver<2, 0>::triangulation.reset_all_manifolds();
31 
32  double dif_norm;
33 
34  for (auto cell : Solver<2, 0>::triangulation.active_cell_iterators()) {
35  for (unsigned int f = 0; f < GeometryInfo<2>::faces_per_cell; f++) {
36  dif_norm = std::abs(cell->face(f)->vertex(0).norm() -
37  cell->face(f)->vertex(1).norm());
38 
39  if ((dif_norm < eps) && (cell->center().norm() > rd1))
40  cell->face(f)->set_all_manifold_ids(1);
41  }
42  }
43 
44  Solver<2, 0>::triangulation.set_manifold(1, sphere);
45 #endif
46 }
47 
48 void
49 SolverMMSVTII_T::fill_dirichlet_stack()
50 {
51  Solver<2, 0>::dirichlet_stack = { { bid_dirichlet, &dirichlet_bc } };
52 }
53 
54 void
55 SolverMMSVTII_T::solve()
56 {
57  SolverControl control(1000 * Solver<2, 0>::system_rhs.size(),
58  1e-12 * Solver<2, 0>::system_rhs.l2_norm(),
59  false,
60  false);
61 
63  control.enable_history_data();
64 
65  GrowingVectorMemory<Vector<double>> memory;
66  SolverCG<Vector<double>> cg(control, memory);
67 
68  PreconditionSSOR<SparseMatrix<double>> preconditioner;
69  preconditioner.initialize(Solver<2, 0>::system_matrix, 1.2);
70 
74  preconditioner);
75 
77 
79  const std::vector<double> history_data = control.get_history_data();
80 
81  std::ofstream ofs(fname + "_cg_convergence.csv");
82 
83  unsigned int i = 1;
84  for (auto item : history_data) {
85  ofs << i << ", " << item << "\n";
86  i++;
87  }
88  ofs.close();
89  }
90 }
91 
92 //-----------------------------------------------------------------------------
93 //-----------------------------------------------------------------------------
94 //-----------------------------------------------------------------------------
95 
96 void
97 SolverMMSVTII_A::fill_dirichlet_stack()
98 {
99  Solver2<2, 2>::dirichlet_stack = { { bid_dirichlet, &dirichlet_bc } };
100 }
101 
102 void
103 SolverMMSVTII_A::solve()
104 {
105  SolverControl control(1000 * Solver2<2, 2>::system_rhs.size(),
106  1e-12 * Solver2<2, 2>::system_rhs.l2_norm(),
107  false,
108  false);
109 
110  if (log_cg_convergence)
111  control.enable_history_data();
112 
113  GrowingVectorMemory<Vector<double>> memory;
114  SolverCG<Vector<double>> cg(control, memory);
115 
116  PreconditionSSOR<SparseMatrix<double>> preconditioner;
117  preconditioner.initialize(Solver2<2, 2>::system_matrix, 1.2);
118 
122  preconditioner);
123 
125 
127  const std::vector<double> history_data = control.get_history_data();
128 
129  std::ofstream ofs(fname + "_cg_convergence.csv");
130 
131  unsigned int i = 1;
132  for (auto item : history_data) {
133  ofs << i << ", " << item << "\n";
134  i++;
135  }
136  ofs.close();
137  }
138 }
const types::boundary_id bid_dirichlet
The Dirichlet boundary condition will be applied to the boundaries with ID = 1.
Definition: settings.hpp:61
const double rd1
The radius of the circle (sphere) that encloses the square (cube) in the middle of the mesh.
Definition: settings.hpp:55
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:73
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:98
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.
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.
std::map< types::boundary_id, const Function< dim > * > dirichlet_stack
A map that contains pairs of boundary IDs and the corresponding Dirichlet boundary conditions....
SparseMatrix< double > system_matrix
The system matrix.
Vector< double > solution
The solution vector, that is, degrees of freedom yielded by the simulation.
Vector< double > system_rhs
The system right-hand side vector.
AffineConstraints< double > constraints
The constraints associated with the Dirichlet boundary conditions.