12 #ifndef SolverSCHAXI_H__
13 #define SolverSCHAXI_H__
15 #define BOOST_ALLOW_DEPRECATED_HEADERS
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>
22 #include <deal.II/numerics/fe_field_function.h>
24 #include "exact_solution.hpp"
25 #include "settings.hpp"
26 #include "static_scalar_solver.hpp"
28 #define TMR(__name) TimerOutput::Scope timer_section(timer, __name)
30 using namespace StaticScalarSolver;
37 template<
bool is_cylinder>
61 unsigned int mapping_degree,
79 ? TimerOutput::summary
82 TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
98 const std::string fname;
101 const dealii::Functions::ZeroFunction<2> dirichlet_function;
105 const unsigned int nr_slice_global_refs = 3;
109 Triangulation<1, 2> triangulation_slice;
111 DoFHandler<1, 2> dof_handler_slice;
112 Vector<double> solution_slice;
114 SphericalManifold<2> sphere;
116 virtual void make_mesh() override final;
117 virtual
void fill_dirichlet_stack() override final;
118 virtual
void solve() override final;
121 void data_slice(std::
string fname);
124 template<
bool is_cylinder>
126 SolverSCHAXI<is_cylinder>::data_slice(std::
string fname)
128 triangulation_slice.refine_global(nr_slice_global_refs);
130 dof_handler_slice.reinit(triangulation_slice);
131 dof_handler_slice.distribute_dofs(fe_slice);
132 solution_slice.reinit(dof_handler_slice.n_dofs());
137 VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
140 DataOut<1, 2> data_out;
142 data_out.attach_dof_handler(dof_handler_slice);
143 data_out.add_data_vector(solution_slice,
"solution_slice");
144 data_out.build_patches();
146 std::ofstream out(fname +
"_slice" +
".gpi");
148 data_out.write_gnuplot(out);
152 template<
bool is_cylinder>
159 template<
bool is_cylinder>
168 if (log_cg_convergence)
169 control.enable_history_data();
171 GrowingVectorMemory<Vector<double>> memory;
172 SolverCG<Vector<double>> cg(control, memory);
174 PreconditionJacobi<SparseMatrix<double>> preconditioner;
184 if (log_cg_convergence) {
185 const std::vector<double> history_data = control.get_history_data();
187 std::ofstream ofs(fname +
"_cg_convergence.csv");
190 for (
auto item : history_data) {
191 ofs << i <<
", " << item <<
"\n";
Describes exact solution, , of the Axisymmetric - surface charge (sch-axi/) numerical experiment.
Global settings for the Axisymmetric - surface charge (sch-axi/) numerical experiment.
const bool print_time_tables
If set to true, the program will print time tables on the screen.
Implements the Axisymmetric - surface charge (sch-axi/) numerical experiment.
SolverSCHAXI(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Solves static scalar boundary value problem.
void run()
Runs the simulation.