12 #ifndef SolverINTAXI_H__
13 #define SolverINTAXI_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/grid_out.h>
21 #include <deal.II/grid/grid_tools.h>
22 #include <deal.II/grid/manifold_lib.h>
24 #include <deal.II/numerics/fe_field_function.h>
26 #include "exact_solution.hpp"
27 #include "settings.hpp"
28 #include "static_scalar_solver.hpp"
30 #define TMR(__name) TimerOutput::Scope timer_section(timer, __name)
32 using namespace StaticScalarSolver;
38 template<
bool is_cylinder>
62 unsigned int mapping_degree,
77 , dirichlet_function_in(1.0)
80 if (CYLINDER__ == 1) {
81 fname_mesh =
"../../gmsh/data/cylinder_r" + std::to_string(r) +
".msh";
83 fname_mesh =
"../../gmsh/data/sphere_r" + std::to_string(r) +
".msh";
87 ? TimerOutput::summary
90 TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
105 const unsigned int r;
106 const std::string fname;
107 std::string fname_mesh;
110 const dealii::Functions::ZeroFunction<2> dirichlet_function_out;
111 const dealii::Functions::ConstantFunction<2> dirichlet_function_in;
115 const unsigned int nr_slice_global_refs = 3;
119 Triangulation<1, 2> triangulation_slice;
121 DoFHandler<1, 2> dof_handler_slice;
122 Vector<double> solution_slice;
124 SphericalManifold<2> sphere;
126 virtual void make_mesh() override final;
127 virtual
void fill_dirichlet_stack() override final;
128 virtual
void solve() override final;
131 void data_slice(std::
string fname);
134 template<
bool is_cylinder>
139 { bid_out, &dirichlet_function_out } };
142 template<
bool is_cylinder>
146 GridGenerator::hyper_cube(triangulation_slice, a + eps, b - eps);
147 triangulation_slice.refine_global(nr_slice_global_refs);
149 dof_handler_slice.reinit(triangulation_slice);
150 dof_handler_slice.distribute_dofs(fe_slice);
151 solution_slice.reinit(dof_handler_slice.n_dofs());
156 VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
158 DataOut<1, 2> data_out;
160 data_out.attach_dof_handler(dof_handler_slice);
161 data_out.add_data_vector(solution_slice,
"solution_slice");
162 data_out.build_patches();
164 std::ofstream out(fname +
"_slice" +
".gpi");
166 data_out.write_gnuplot(out);
170 template<
bool is_cylinder>
179 if (log_cg_convergence)
180 control.enable_history_data();
182 GrowingVectorMemory<Vector<double>> memory;
183 SolverCG<Vector<double>> cg(control, memory);
185 PreconditionJacobi<SparseMatrix<double>> preconditioner;
195 if (log_cg_convergence) {
196 const std::vector<double> history_data = control.get_history_data();
198 std::ofstream ofs(fname +
"_cg_convergence.csv");
201 for (
auto item : history_data) {
202 ofs << i <<
", " << item <<
"\n";
209 template<
bool is_cylinder>
216 std::ifstream ifs(fname_mesh);
217 gridin.read_msh(ifs);
221 if (cell->center()[0] < d) {
222 cell->set_material_id(mid_1);
224 cell->set_material_id(mid_2);
229 if (cell->center().norm() < d) {
230 cell->set_material_id(mid_1);
232 cell->set_material_id(mid_2);
Describes exact solution, , of the Axisymmetric - interface between dielectrics (int-axi/) numerical ...
Global settings for the Axisymmetric - interface between dielectrics (int-axi/) numerical experiment.
const bool print_time_tables
If set to true, the program will print time tables on the screen.
Implements the Axisymmetric - interface between dielectrics (int-axi/) numerical experiment.
SolverINTAXI(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Solves static scalar boundary value problem.
void run()
Runs the simulation.