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 SolverFLCAXI_H__
13 #define SolverFLCAXI_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/manifold_lib.h>
21 
22 #include <deal.II/numerics/fe_field_function.h>
23 
24 #include "exact_solution.hpp"
25 #include "settings.hpp"
26 #include "static_scalar_solver.hpp"
27 
28 #define TMR(__name) TimerOutput::Scope timer_section(timer, __name)
29 
30 using namespace StaticScalarSolver;
31 
37 template<bool is_cylinder>
39  : public SettingsFLCAXI
40  , public Solver<2>
41 {
42 public:
43  SolverFLCAXI() = delete;
44 
60  SolverFLCAXI(unsigned int p,
61  unsigned int mapping_degree,
62  unsigned int r,
63  std::string fname)
64  : Solver<2>(p,
65  mapping_degree,
66  0,
67  fname,
68  &exact_solution,
69  true,
70  false,
71  SettingsFLCAXI::print_time_tables,
72  SettingsFLCAXI::project_exact_solution,
73  true)
74  , r(r)
75  , fname(fname)
76  , dirichlet_function_in(1.0)
77  , fe_slice(1)
78  {
79  TimerOutput::OutputFrequency tf = (SettingsFLCAXI::print_time_tables)
80  ? TimerOutput::summary
81  : TimerOutput::never;
82 
83  TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
84 
85  {
86  TMR("Solver run");
88  }
89  {
90  TMR("Data slice");
91  data_slice(fname);
92  }
93  }
94 
95  ~SolverFLCAXI() = default;
96 
97 private:
98  const unsigned int r;
99  const std::string fname;
100 
101  const ExactSolutionFLCAXI_PHI<is_cylinder> exact_solution;
102  const dealii::Functions::ZeroFunction<2> dirichlet_function_out;
103  const dealii::Functions::ConstantFunction<2> dirichlet_function_in;
104 
105  // The amount of global mesh refinements that need to be done to the
106  // one- dimensional mesh used for the plot of potential vs. x coordinate.
107  const unsigned int nr_slice_global_refs = 3;
108 
109  // These four data members are needed for making the plot of potential
110  // vs. x coordinate.
111  Triangulation<1, 2> triangulation_slice;
112  FE_Q<1, 2> fe_slice;
113  DoFHandler<1, 2> dof_handler_slice;
114  Vector<double> solution_slice;
115 
116  SphericalManifold<2> sphere;
117 
118  virtual void make_mesh() override final;
119  virtual void fill_dirichlet_stack() override final;
120  virtual void solve() override final;
121 
122  // This function makes the plot of potential vs. x coordinate.
123  void data_slice(std::string fname);
124 };
125 
126 template<bool is_cylinder>
127 void
128 SolverFLCAXI<is_cylinder>::fill_dirichlet_stack()
129 {
130  Solver<2>::dirichlet_stack = { { bid_in, &dirichlet_function_in },
131  { bid_out, &dirichlet_function_out } };
132 }
133 
134 template<bool is_cylinder>
135 void
136 SolverFLCAXI<is_cylinder>::data_slice(std::string fname)
137 {
138  GridGenerator::hyper_cube(triangulation_slice, a + eps, b - eps);
139  triangulation_slice.refine_global(nr_slice_global_refs);
140 
141  dof_handler_slice.reinit(triangulation_slice);
142  dof_handler_slice.distribute_dofs(fe_slice);
143  solution_slice.reinit(dof_handler_slice.n_dofs());
144 
145  Functions::FEFieldFunction<2> potential(Solver<2>::dof_handler,
147 
148  VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
149 
150  // DataOut<1,DoFHandler<1,2>> data_out;
151  DataOut<1, 2> data_out;
152 
153  data_out.attach_dof_handler(dof_handler_slice);
154  data_out.add_data_vector(solution_slice, "solution_slice");
155  data_out.build_patches();
156 
157  std::ofstream out(fname + "_slice" + ".gpi");
158 
159  data_out.write_gnuplot(out);
160  out.close();
161 }
162 
163 template<bool is_cylinder>
164 void
166 {
167  SolverControl control(Solver<2>::system_rhs.size(),
168  1e-8 * Solver<2>::system_rhs.l2_norm(),
169  false,
170  false);
171 
172  if (log_cg_convergence)
173  control.enable_history_data();
174 
175  GrowingVectorMemory<Vector<double>> memory;
176  SolverCG<Vector<double>> cg(control, memory);
177 
178  PreconditionJacobi<SparseMatrix<double>> preconditioner;
179  preconditioner.initialize(Solver<2>::system_matrix, 1.0);
180 
181  cg.solve(Solver<2>::system_matrix,
184  preconditioner);
185 
187 
188  if (log_cg_convergence) {
189  const std::vector<double> history_data = control.get_history_data();
190 
191  std::ofstream ofs(fname + "_cg_convergence.csv");
192 
193  unsigned int i = 1;
194  for (auto item : history_data) {
195  ofs << i << ", " << item << "\n";
196  i++;
197  }
198  ofs.close();
199  }
200 }
201 
202 #endif
Describes the exact solution, , of the Axisymmetric - floating conductor (flc-axi/) numerical experim...
Global settings for the Axisymmetric - floating conductor (flc-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:113
Implements the solver of the Axisymmetric - floating conductor (flc-axi/) numerical experiment.
Definition: solver.hpp:41
SolverFLCAXI(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Definition: solver.hpp:60
Solves static scalar boundary value problem.
void run()
Runs the simulation.