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 SolverSCHAXI_H__
13 #define SolverSCHAXI_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 SettingsSCHAXI
40  , public Solver<2>
41 {
42 public:
43  SolverSCHAXI() = delete;
44 
60  SolverSCHAXI(unsigned int p,
61  unsigned int mapping_degree,
62  unsigned int r,
63  std::string fname)
64  : Solver<2>(p,
65  mapping_degree,
66  1,
67  fname,
68  &exact_solution,
69  true,
70  false,
71  SettingsSCHAXI::print_time_tables,
72  SettingsSCHAXI::project_exact_solution,
73  true)
74  , r(r)
75  , fname(fname)
76  , fe_slice(1)
77  {
78  TimerOutput::OutputFrequency tf = (SettingsSCHAXI::print_time_tables)
79  ? TimerOutput::summary
80  : TimerOutput::never;
81 
82  TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
83 
84  {
85  TMR("Solver run");
87  }
88  {
89  TMR("Data slice");
90  data_slice(fname);
91  }
92  }
93 
94  ~SolverSCHAXI() = default;
95 
96 private:
97  const unsigned int r;
98  const std::string fname;
99 
100  const ExactSolutionSCHAXI_PHI<static_cast<bool>(is_cylinder)> exact_solution;
101  const dealii::Functions::ZeroFunction<2> dirichlet_function;
102 
103  // The amount of global mesh refinements that need to be done to the
104  // one- dimensional mesh used for the plot of potential vs. x coordinate.
105  const unsigned int nr_slice_global_refs = 3;
106 
107  // These four data members are needed for making the plot of potential
108  // vs. x coordinate.
109  Triangulation<1, 2> triangulation_slice;
110  FE_Q<1, 2> fe_slice;
111  DoFHandler<1, 2> dof_handler_slice;
112  Vector<double> solution_slice;
113 
114  SphericalManifold<2> sphere;
115 
116  virtual void make_mesh() override final;
117  virtual void fill_dirichlet_stack() override final;
118  virtual void solve() override final;
119 
120  // This function makes the plot of potential vs. \f$x\f$ coordinate.
121  void data_slice(std::string fname);
122 };
123 
124 template<bool is_cylinder>
125 void
126 SolverSCHAXI<is_cylinder>::data_slice(std::string fname)
127 {
128  triangulation_slice.refine_global(nr_slice_global_refs);
129 
130  dof_handler_slice.reinit(triangulation_slice);
131  dof_handler_slice.distribute_dofs(fe_slice);
132  solution_slice.reinit(dof_handler_slice.n_dofs());
133 
134  Functions::FEFieldFunction<2> potential(Solver<2>::dof_handler,
136 
137  VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
138 
139  // DataOut<1,DoFHandler<1, 2>> data_out;
140  DataOut<1, 2> data_out;
141 
142  data_out.attach_dof_handler(dof_handler_slice);
143  data_out.add_data_vector(solution_slice, "solution_slice");
144  data_out.build_patches();
145 
146  std::ofstream out(fname + "_slice" + ".gpi");
147 
148  data_out.write_gnuplot(out);
149  out.close();
150 }
151 
152 template<bool is_cylinder>
153 void
155 {
156  Solver<2>::dirichlet_stack = { { bid, &dirichlet_function } };
157 }
158 
159 template<bool is_cylinder>
160 void
162 {
163  SolverControl control(Solver<2>::system_rhs.size(),
164  1e-8 * Solver<2>::system_rhs.l2_norm(),
165  false,
166  false);
167 
168  if (log_cg_convergence)
169  control.enable_history_data();
170 
171  GrowingVectorMemory<Vector<double>> memory;
172  SolverCG<Vector<double>> cg(control, memory);
173 
174  PreconditionJacobi<SparseMatrix<double>> preconditioner;
175  preconditioner.initialize(Solver<2>::system_matrix, 1.0);
176 
177  cg.solve(Solver<2>::system_matrix,
180  preconditioner);
181 
183 
184  if (log_cg_convergence) {
185  const std::vector<double> history_data = control.get_history_data();
186 
187  std::ofstream ofs(fname + "_cg_convergence.csv");
188 
189  unsigned int i = 1;
190  for (auto item : history_data) {
191  ofs << i << ", " << item << "\n";
192  i++;
193  }
194  ofs.close();
195  }
196 }
197 
198 #endif
Describes exact solution, , of the Axisymmetric - surface charge (sch-axi/) numerical experiment.
Global settings for the Axisymmetric - surface charge (sch-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:92
Implements the Axisymmetric - surface charge (sch-axi/) numerical experiment.
Definition: solver.hpp:41
SolverSCHAXI(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.