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 SolverSSOLIIIAXI_H__
13 #define SolverSSOLIIIAXI_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 
38  : public SettingsSSOLIIIAXI
39  , public Solver<2>
40 {
41 public:
42  SolverSSOLIIIAXI() = delete;
43 
59  SolverSSOLIIIAXI(unsigned int p,
60  unsigned int mapping_degree,
61  unsigned int r,
62  std::string fname)
63  : Solver<2>(p,
64  mapping_degree,
65  1,
66  fname,
67  nullptr,
68  true,
69  true,
70  SettingsSSOLIIIAXI::print_time_tables,
71  SettingsSSOLIIIAXI::project_exact_solution,
72  true)
73  , r(r)
74  , fname(fname)
75  , fe_slice(1)
76  {
77  TimerOutput::OutputFrequency tf = (SettingsSSOLIIIAXI::print_time_tables)
78  ? TimerOutput::summary
79  : TimerOutput::never;
80 
81  TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
82 
83  {
84  TMR("Solver run");
86  }
87  {
88  TMR("Data slice");
89  data_slice(fname);
90  }
91  }
92 
93  ~SolverSSOLIIIAXI() = default;
94 
95 private:
96  const unsigned int r;
97  const std::string fname;
98 
99  const dealii::Functions::ZeroFunction<2> 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. r 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. r coordinate.
107  Triangulation<1, 2> triangulation_slice;
108  FE_Q<1, 2> fe_slice;
109  DoFHandler<1, 2> dof_handler_slice;
110  Vector<double> solution_slice;
111 
112  SphericalManifold<2> 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. r coordinate.
121  void data_slice(std::string fname);
122 };
123 
124 #endif
Global settings for the Axisymmetric - thick spherical coil with magnetic core (ssol-iii-axi/) numeri...
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:141
Implements the Axisymmetric - thick spherical coil with magnetic core (ssol-iii-axi/) numerical exper...
Definition: solver.hpp:40
SolverSSOLIIIAXI(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.