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 SolverINTAXI_H__
13 #define SolverINTAXI_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 #include <deal.II/grid/manifold_lib.h>
23 
24 #include <deal.II/numerics/fe_field_function.h>
25 
26 #include "exact_solution.hpp"
27 #include "settings.hpp"
28 #include "static_scalar_solver.hpp"
29 
30 #define TMR(__name) TimerOutput::Scope timer_section(timer, __name)
31 
32 using namespace StaticScalarSolver;
33 
38 template<bool is_cylinder>
40  : public SettingsINTAXI
41  , public Solver<2>
42 {
43 public:
44  SolverINTAXI() = delete;
45 
61  SolverINTAXI(unsigned int p,
62  unsigned int mapping_degree,
63  unsigned int r,
64  std::string fname)
65  : Solver<2>(p,
66  mapping_degree,
67  0,
68  fname,
69  &exact_solution,
70  true,
71  false,
72  SettingsINTAXI::print_time_tables,
73  SettingsINTAXI::project_exact_solution,
74  true)
75  , r(r)
76  , fname(fname)
77  , dirichlet_function_in(1.0)
78  , fe_slice(1)
79  {
80  if (CYLINDER__ == 1) {
81  fname_mesh = "../../gmsh/data/cylinder_r" + std::to_string(r) + ".msh";
82  } else {
83  fname_mesh = "../../gmsh/data/sphere_r" + std::to_string(r) + ".msh";
84  }
85 
86  TimerOutput::OutputFrequency tf = (SettingsINTAXI::print_time_tables)
87  ? TimerOutput::summary
88  : TimerOutput::never;
89 
90  TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
91 
92  {
93  TMR("Solver run");
95  }
96  {
97  TMR("Data slice");
98  data_slice(fname);
99  }
100  }
101 
102  ~SolverINTAXI() = default;
103 
104 private:
105  const unsigned int r;
106  const std::string fname;
107  std::string fname_mesh;
108 
109  const ExactSolutionINTAXI_PHI<is_cylinder> exact_solution;
110  const dealii::Functions::ZeroFunction<2> dirichlet_function_out;
111  const dealii::Functions::ConstantFunction<2> dirichlet_function_in;
112 
113  // The amount of global mesh refinements that need to be done to the
114  // one-dimensional mesh used for the plot of potential vs. x coordinate.
115  const unsigned int nr_slice_global_refs = 3;
116 
117  // These four data members are needed for making the plot of potential
118  // vs. x coordinate.
119  Triangulation<1, 2> triangulation_slice;
120  FE_Q<1, 2> fe_slice;
121  DoFHandler<1, 2> dof_handler_slice;
122  Vector<double> solution_slice;
123 
124  SphericalManifold<2> sphere;
125 
126  virtual void make_mesh() override final;
127  virtual void fill_dirichlet_stack() override final;
128  virtual void solve() override final;
129 
130  // This function makes the plot of potential vs. x coordinate.
131  void data_slice(std::string fname);
132 };
133 
134 template<bool is_cylinder>
135 void
136 SolverINTAXI<is_cylinder>::fill_dirichlet_stack()
137 {
138  Solver<2>::dirichlet_stack = { { bid_in, &dirichlet_function_in },
139  { bid_out, &dirichlet_function_out } };
140 }
141 
142 template<bool is_cylinder>
143 void
144 SolverINTAXI<is_cylinder>::data_slice(std::string fname)
145 {
146  GridGenerator::hyper_cube(triangulation_slice, a + eps, b - eps);
147  triangulation_slice.refine_global(nr_slice_global_refs);
148 
149  dof_handler_slice.reinit(triangulation_slice);
150  dof_handler_slice.distribute_dofs(fe_slice);
151  solution_slice.reinit(dof_handler_slice.n_dofs());
152 
153  Functions::FEFieldFunction<2> potential(Solver<2>::dof_handler,
155 
156  VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
157 
158  DataOut<1, 2> data_out;
159 
160  data_out.attach_dof_handler(dof_handler_slice);
161  data_out.add_data_vector(solution_slice, "solution_slice");
162  data_out.build_patches();
163 
164  std::ofstream out(fname + "_slice" + ".gpi");
165 
166  data_out.write_gnuplot(out);
167  out.close();
168 }
169 
170 template<bool is_cylinder>
171 void
173 {
174  SolverControl control(Solver<2>::system_rhs.size(),
175  1e-8 * Solver<2>::system_rhs.l2_norm(),
176  false,
177  false);
178 
179  if (log_cg_convergence)
180  control.enable_history_data();
181 
182  GrowingVectorMemory<Vector<double>> memory;
183  SolverCG<Vector<double>> cg(control, memory);
184 
185  PreconditionJacobi<SparseMatrix<double>> preconditioner;
186  preconditioner.initialize(Solver<2>::system_matrix, 1.0);
187 
188  cg.solve(Solver<2>::system_matrix,
191  preconditioner);
192 
194 
195  if (log_cg_convergence) {
196  const std::vector<double> history_data = control.get_history_data();
197 
198  std::ofstream ofs(fname + "_cg_convergence.csv");
199 
200  unsigned int i = 1;
201  for (auto item : history_data) {
202  ofs << i << ", " << item << "\n";
203  i++;
204  }
205  ofs.close();
206  }
207 }
208 
209 template<bool is_cylinder>
210 void
212 {
213  GridIn<2> gridin;
214  gridin.attach_triangulation(Solver<2>::triangulation);
215 
216  std::ifstream ifs(fname_mesh);
217  gridin.read_msh(ifs);
218 
219  if (is_cylinder) {
220  for (auto cell : Solver<2>::triangulation.active_cell_iterators()) {
221  if (cell->center()[0] < d) {
222  cell->set_material_id(mid_1);
223  } else {
224  cell->set_material_id(mid_2);
225  }
226  }
227  } else {
228  for (auto cell : Solver<2>::triangulation.active_cell_iterators()) {
229  if (cell->center().norm() < d) {
230  cell->set_material_id(mid_1);
231  } else {
232  cell->set_material_id(mid_2);
233  }
234  }
235 
236  Solver<2>::triangulation.set_all_manifold_ids(1);
237  Solver<2>::triangulation.set_manifold(1, sphere);
238  }
239 }
240 
241 #endif
Describes exact solution, , of the Axisymmetric - interface between dielectrics (int-axi/) numerical ...
Global settings for the Axisymmetric - interface between dielectrics (int-axi/) numerical experiment.
Definition: settings.hpp:26
const bool print_time_tables
If set to true, the program will print time tables on the screen.
Definition: settings.hpp:96
Implements the Axisymmetric - interface between dielectrics (int-axi/) numerical experiment.
Definition: solver.hpp:42
SolverINTAXI(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Definition: solver.hpp:61
Solves static scalar boundary value problem.
void run()
Runs the simulation.