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 SolverMMSAXI_H__
13 #define SolverMMSAXI_H__
14 
15 #define BOOST_ALLOW_DEPRECATED_HEADERS
16 
17 #include <deal.II/base/function.h>
18 #include <deal.II/grid/grid_generator.h>
19 #include <deal.II/grid/grid_in.h>
20 #include <deal.II/grid/grid_out.h>
21 #include <deal.II/grid/grid_tools.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 
37 template<int dim>
39  : public SettingsMMSAXI
40  , public Solver<dim>
41 {
42 public:
43  SolverMMSAXI() = delete;
44 
58  SolverMMSAXI(unsigned int p, unsigned int r, std::string fname)
59  : Solver<dim>(p,
60  1,
61  1,
62  fname,
63  &exact_solution,
64 #if DIMENSION__ == 2
65  true,
66 #endif
67 #if DIMENSION__ == 3
68  false,
69 #endif
70  false,
71  print_time_tables,
72  project_exact_solution)
73  , r(r)
74  , fname(fname)
75  {
76 #if DIMENSION__ == 2
77  fname_mesh = "../../gmsh/data/cylinder_2d_r" + std::to_string(r) + ".msh";
78 #endif
79 #if DIMENSION__ == 3
80  fname_mesh = "../../gmsh/data/cylinder_3d_r" + std::to_string(r) + ".msh";
81 #endif
82 
83  TimerOutput::OutputFrequency tf =
84  (print_time_tables) ? TimerOutput::summary : TimerOutput::never;
85 
86  TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
87 
88  {
89  TMR("Solver run");
91  }
92  }
93 
94  ~SolverMMSAXI() = default;
95 
96 private:
97  const unsigned int r;
98  const std::string fname;
99  const ExactSolutionMMSAXI_PHI<dim> exact_solution;
100 
101  std::string fname_mesh;
102 
103  virtual void make_mesh() override final;
104  virtual void fill_dirichlet_stack() override final;
105  virtual void solve() override final;
106 };
107 
108 template<int dim>
109 void
110 SolverMMSAXI<dim>::fill_dirichlet_stack()
111 {
112  Solver<dim>::dirichlet_stack = { { bid_dirichlet, &exact_solution } };
113 }
114 
115 template<int dim>
116 void
118 {
119  GridIn<dim> gridin;
120  gridin.attach_triangulation(Solver<dim>::triangulation);
121 
122  std::ifstream ifs(fname_mesh);
123  gridin.read_msh(ifs);
124 }
125 
126 template<int dim>
127 void
129 {
130  SolverControl control(Solver<dim>::system_rhs.size(),
131  1e-8 * Solver<dim>::system_rhs.l2_norm(),
132  false,
133  false);
134 
135  if (log_cg_convergence)
136  control.enable_history_data();
137 
138  GrowingVectorMemory<Vector<double>> memory;
139  SolverCG<Vector<double>> cg(control, memory);
140 
141  PreconditionJacobi<SparseMatrix<double>> preconditioner;
142  preconditioner.initialize(Solver<dim>::system_matrix, 1.0);
143 
147  preconditioner);
148 
150 
151  if (log_cg_convergence) {
152  const std::vector<double> history_data = control.get_history_data();
153 
154  std::ofstream ofs(fname + "_cg_convergence.csv");
155 
156  unsigned int i = 1;
157  for (auto item : history_data) {
158  ofs << i << ", " << item << "\n";
159  i++;
160  }
161  ofs.close();
162  }
163 }
164 
165 #endif
Describes the exact solution, , of the Axisymmetric - method of manufactured solutions (mms-axi/) num...
Global settings for the Axisymmetric - method of manufactured solutions (mms-axi/) numerical experime...
Definition: settings.hpp:26
Implements the Axisymmetric - method of manufactured solutions (mms-axi/) numerical experiment.
Definition: solver.hpp:41
SolverMMSAXI(unsigned int p, unsigned int r, std::string fname)
Definition: solver.hpp:58
Solves static scalar boundary value problem.
void run()
Runs the simulation.