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;
62 unsigned int mapping_degree,
73 project_exact_solution,
77 , dirichlet_function_in(1.0)
80 if (DIMENSION__ == 2) {
81 fname_mesh =
"../../gmsh/data/ring_r" + std::to_string(r) +
".msh";
83 fname_mesh =
"../../gmsh/data/shell_r" + std::to_string(r) +
".msh";
86 TimerOutput::OutputFrequency tf =
87 (print_time_tables) ? TimerOutput::summary : TimerOutput::never;
89 TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
92 TMR(
"Solver<dim>::run");
104 const unsigned int r;
105 const std::string fname;
106 std::string fname_mesh;
109 const dealii::Functions::ZeroFunction<dim> dirichlet_function_out;
110 const dealii::Functions::ConstantFunction<dim> dirichlet_function_in;
114 const unsigned int nr_slice_global_refs = 10;
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;
123 SphericalManifold<dim> sphere;
125 virtual void make_mesh() override final;
126 virtual
void fill_dirichlet_stack() override final;
127 virtual
void solve() override final;
129 void mark_materials();
132 void data_slice(std::
string fname);
140 { bid_out, &dirichlet_function_out } };
150 std::ifstream ifs(fname_mesh);
151 gridin.read_msh(ifs);
154 if (cell->center().norm() < d)
155 cell->set_material_id(mid_1);
157 if (cell->center().norm() > d)
158 cell->set_material_id(mid_2);
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);
165 cell->face(f)->set_all_boundary_ids(bid_in);
179 GridGenerator::hyper_cube(triangulation_slice, a + eps, b - eps);
180 triangulation_slice.refine_global(nr_slice_global_refs);
182 dof_handler_slice.reinit(triangulation_slice);
183 dof_handler_slice.distribute_dofs(fe_slice);
184 solution_slice.reinit(dof_handler_slice.n_dofs());
189 VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
192 DataOut<1, dim> data_out;
194 data_out.attach_dof_handler(dof_handler_slice);
195 data_out.add_data_vector(solution_slice,
"solution_slice");
196 data_out.build_patches();
198 std::ofstream out(fname +
"_slice" +
".gpi");
200 data_out.write_gnuplot(out);
213 if (log_cg_convergence)
214 control.enable_history_data();
216 GrowingVectorMemory<Vector<double>> memory;
217 SolverCG<Vector<double>> cg(control, memory);
219 PreconditionJacobi<SparseMatrix<double>> preconditioner;
229 if (log_cg_convergence) {
230 const std::vector<double> history_data = control.get_history_data();
232 std::ofstream ofs(fname +
"_cg_convergence.csv");
235 for (
auto item : history_data) {
236 ofs << i <<
", " << item <<
"\n";
Describes the exact solution, , of the Interface between dielectrics (int/) numerical experiment.
Global settings for the Interface between dielectrics (int/) numerical experiment.
Implements the solver of the Interface between dielectrics (int/) numerical experiment.
SolverINT(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Solves static scalar boundary value problem.
void run()
Runs the simulation.