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 SolverABC_H__
13 #define SolverABC_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<int dim>
38 class SolverABC
39  : public SettingsABC
40  , public Solver<dim>
41 {
42 public:
43  SolverABC() = delete;
44 
62  SolverABC(unsigned int m,
63  unsigned int p,
64  unsigned int mapping_degree,
65  unsigned int r,
66  std::string fname)
67  : Solver<dim>(p,
68  mapping_degree,
69  0,
70  fname,
71  &exact_solution,
72  false,
73  false,
74  print_time_tables,
75  project_exact_solution,
76  true)
77  , m(m)
78  , r(r)
79  , fname(fname)
80  , dirichlet_function_left(-1.0)
81  , dirichlet_function_right(1.0)
82  , dirichlet_function_in(1.0)
83  , circle_left(Point<2>(-x0, 0.0))
84  , circle_right(Point<2>(x0, 0.0))
85  , circle_infty(Point<2>(0.0, 0.0))
86  {
87  TimerOutput::OutputFrequency tf =
88  (print_time_tables) ? TimerOutput::summary : TimerOutput::never;
89 
90  TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
91 
92  {
93  TMR("Solver run");
95  }
96  }
97 
98  ~SolverABC() = default;
99 
100 private:
101  const unsigned int m;
102  const unsigned int r;
103  const std::string fname;
104 
105  const ExactSolutionABC_PHI<dim> exact_solution;
106  const Functions::ZeroFunction<dim> dirichlet_function_infty;
107  const Functions::ConstantFunction<dim> dirichlet_function_left;
108  const Functions::ConstantFunction<dim> dirichlet_function_right;
109  const Functions::ConstantFunction<dim> dirichlet_function_in;
110 
111  const SphericalManifold<2> circle_left;
112  const SphericalManifold<2> circle_right;
113  const SphericalManifold<2> circle_infty;
114 
115  SphericalManifold<3> sphere;
116 
117  // Sets the ID of the outer boundary to 0, 2, 5 depending on the value of
118  // BC_TYPE__ set CMakeLists.txt. The ID of the outer boundary in the geo
119  // file equals 2.
120  void renumber_boundaries();
121  virtual void make_mesh() override final;
122  virtual void fill_dirichlet_stack() override final;
123  virtual void solve() override final;
124 };
125 
126 template<int dim>
127 void
128 SolverABC<dim>::renumber_boundaries()
129 {
130  for (auto cell : Solver<dim>::triangulation.active_cell_iterators()) {
131  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f) {
132  if (cell->face(f)->at_boundary()) {
133  if (cell->face(f)->boundary_id() == 2)
134  cell->face(f)->set_all_boundary_ids(bid_infty);
135  }
136  }
137  }
138 }
139 
140 template<int dim>
141 void
143 {
144  SolverControl control(Solver<dim>::system_rhs.size(),
145  1e-8 * Solver<dim>::system_rhs.l2_norm(),
146  false,
147  false);
148 
149  if (log_cg_convergence)
150  control.enable_history_data();
151 
152  GrowingVectorMemory<Vector<double>> memory;
153  SolverCG<Vector<double>> cg(control, memory);
154 
155  PreconditionJacobi<SparseMatrix<double>> preconditioner;
156  preconditioner.initialize(Solver<dim>::system_matrix, 1.0);
157 
161  preconditioner);
162 
164 
165  if (log_cg_convergence) {
166  const std::vector<double> history_data = control.get_history_data();
167 
168  std::ofstream ofs(fname + "_cg_convergence.csv");
169 
170  unsigned int i = 1;
171  for (auto item : history_data) {
172  ofs << i << ", " << item << "\n";
173  i++;
174  }
175  ofs.close();
176  }
177 }
178 
179 #endif
Describes exact solution, , of the Asymptotic boundary condition (abc/) numerical experiment.
Global settings for the Asymptotic boundary conditions (abc/) numerical experiment.
Definition: settings.hpp:26
Implements the solver of the Asymptotic boundary condition (abc/) numerical experiment.
Definition: solver.hpp:41
SolverABC(unsigned int m, unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Definition: solver.hpp:62
Solves static scalar boundary value problem.
void run()
Runs the simulation.