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 SolverRHO_H__
13 #define SolverRHO_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 
36 template<int dim>
37 class SolverRHO
38  : public SettingsRHO
39  , public Solver<dim>
40 {
41 public:
42  SolverRHO() = delete;
43 
59  SolverRHO(unsigned int p,
60  unsigned int mapping_degree,
61  unsigned int r,
62  std::string fname)
63  : Solver<dim>(p,
64  mapping_degree,
65  1,
66  fname,
67  &exact_solution,
68  false,
69  false,
70  print_time_tables,
71  project_exact_solution,
72  true)
73  , r(r)
74  , fname(fname)
75  , fe_slice(1)
76  {
77  TimerOutput::OutputFrequency tf =
78  (print_time_tables) ? TimerOutput::summary : TimerOutput::never;
79 
80  TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
81 
82  {
83  TMR("Solver run");
85  }
86  {
87  TMR("Data slice");
88  data_slice(fname);
89  }
90  }
91 
92  ~SolverRHO() = default;
93 
94 private:
95  const unsigned int r;
96  const std::string fname;
97 
98  const ExactSolutionRHO_PHI<dim> exact_solution;
99  const dealii::Functions::ZeroFunction<dim> dirichlet_function;
100 
101  // The amount of global mesh refinements that need to be done to the
102  // one-dimensional mesh used for the plot of potential vs. x coordinate.
103  const unsigned int nr_slice_global_refs = 10;
104 
105  // These four data members are needed for making the plot of potential
106  // vs. \f$x\f$ coordinate.
107  Triangulation<1, dim> triangulation_slice;
108  FE_Q<1, dim> fe_slice;
109  DoFHandler<1, dim> dof_handler_slice;
110  Vector<double> solution_slice;
111 
112  SphericalManifold<dim> sphere;
113 
114  virtual void make_mesh() override final;
115  virtual void fill_dirichlet_stack() override final;
116  virtual void solve() override final;
117 
118  void mark_materials();
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<int dim>
125 void
126 SolverRHO<dim>::fill_dirichlet_stack()
127 {
128  Solver<dim>::dirichlet_stack = { { bid, &dirichlet_function } };
129 }
130 
131 template<int dim>
132 void
134 {
135  Solver<dim>::triangulation.reset_all_manifolds();
136 
137  for (auto cell : Solver<dim>::triangulation.active_cell_iterators()) {
138  if (cell->center().norm() < a) {
139  cell->set_material_id(mid_2);
140  } else {
141  cell->set_material_id(mid_1);
142  }
143 
144  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; f++) {
145  double dif_norm = 0.0;
146  for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; v++)
147  dif_norm += std::abs(cell->face(f)->vertex(0).norm() -
148  cell->face(f)->vertex(v).norm());
149 
150  if ((dif_norm < eps) && (cell->center().norm() > rd1))
151  cell->face(f)->set_all_manifold_ids(1);
152  }
153  }
154 
155  Solver<dim>::triangulation.set_manifold(1, sphere);
156 }
157 
158 template<int dim>
159 void
160 SolverRHO<dim>::data_slice(std::string fname)
161 {
162  GridGenerator::hyper_cube(triangulation_slice, 0.0 + eps, b - eps);
163  triangulation_slice.refine_global(nr_slice_global_refs);
164 
165  dof_handler_slice.reinit(triangulation_slice);
166  dof_handler_slice.distribute_dofs(fe_slice);
167  solution_slice.reinit(dof_handler_slice.n_dofs());
168 
169  Functions::FEFieldFunction<dim> potential(Solver<dim>::dof_handler,
171 
172  VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
173 
174  // DataOut<1,DoFHandler<1,dim>> data_out;
175  DataOut<1, dim> data_out;
176 
177  data_out.attach_dof_handler(dof_handler_slice);
178  data_out.add_data_vector(solution_slice, "solution_slice");
179  data_out.build_patches();
180 
181  std::ofstream out(fname + "_slice" + ".gpi");
182 
183  data_out.write_gnuplot(out);
184  out.close();
185 }
186 
187 template<int dim>
188 void
190 {
191  SolverControl control(Solver<dim>::system_rhs.size(),
192  1e-8 * Solver<dim>::system_rhs.l2_norm(),
193  false,
194  false);
195 
196  if (log_cg_convergence)
197  control.enable_history_data();
198 
199  GrowingVectorMemory<Vector<double>> memory;
200  SolverCG<Vector<double>> cg(control, memory);
201 
202  PreconditionJacobi<SparseMatrix<double>> preconditioner;
203  preconditioner.initialize(Solver<dim>::system_matrix, 1.0);
204 
208  preconditioner);
209 
211 
212  if (log_cg_convergence) {
213  const std::vector<double> history_data = control.get_history_data();
214 
215  std::ofstream ofs(fname + "_cg_convergence.csv");
216 
217  unsigned int i = 1;
218  for (auto item : history_data) {
219  ofs << i << ", " << item << "\n";
220  i++;
221  }
222  ofs.close();
223  }
224 }
225 
226 #endif
Describes exact solution, , of the Volume charge (rho/) numerical experiment in two and three dimensi...
Global settings for the Volume charge (rho/) numerical experiment.
Definition: settings.hpp:25
Implements the Volume charge (rho/) numerical experiment.
Definition: solver.hpp:40
SolverRHO(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Definition: solver.hpp:59
Solves static scalar boundary value problem.
void run()
Runs the simulation.