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