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;
60 unsigned int mapping_degree,
71 project_exact_solution,
77 TimerOutput::OutputFrequency tf =
78 (print_time_tables) ? TimerOutput::summary : TimerOutput::never;
80 TimerOutput timer(std::cout, tf, TimerOutput::cpu_and_wall_times_grouped);
96 const std::string fname;
99 const dealii::Functions::ZeroFunction<dim> dirichlet_function;
103 const unsigned int nr_slice_global_refs = 10;
107 Triangulation<1, dim> triangulation_slice;
108 FE_Q<1, dim> fe_slice;
109 DoFHandler<1, dim> dof_handler_slice;
110 Vector<double> solution_slice;
112 SphericalManifold<dim> sphere;
114 virtual void make_mesh() override final;
115 virtual
void fill_dirichlet_stack() override final;
116 virtual
void solve() override final;
118 void mark_materials();
121 void data_slice(std::
string fname);
138 if (cell->center().norm() < a) {
139 cell->set_material_id(mid_2);
141 cell->set_material_id(mid_1);
144 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; f++) {
145 double dif_norm = 0.0;
146 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; v++)
147 dif_norm += std::abs(cell->face(f)->vertex(0).norm() -
148 cell->face(f)->vertex(v).norm());
150 if ((dif_norm < eps) && (cell->center().norm() > rd1))
151 cell->face(f)->set_all_manifold_ids(1);
162 GridGenerator::hyper_cube(triangulation_slice, 0.0 + eps, b - eps);
163 triangulation_slice.refine_global(nr_slice_global_refs);
165 dof_handler_slice.reinit(triangulation_slice);
166 dof_handler_slice.distribute_dofs(fe_slice);
167 solution_slice.reinit(dof_handler_slice.n_dofs());
172 VectorTools::interpolate(dof_handler_slice, potential, solution_slice);
175 DataOut<1, dim> data_out;
177 data_out.attach_dof_handler(dof_handler_slice);
178 data_out.add_data_vector(solution_slice,
"solution_slice");
179 data_out.build_patches();
181 std::ofstream out(fname +
"_slice" +
".gpi");
183 data_out.write_gnuplot(out);
196 if (log_cg_convergence)
197 control.enable_history_data();
199 GrowingVectorMemory<Vector<double>> memory;
200 SolverCG<Vector<double>> cg(control, memory);
202 PreconditionJacobi<SparseMatrix<double>> preconditioner;
212 if (log_cg_convergence) {
213 const std::vector<double> history_data = control.get_history_data();
215 std::ofstream ofs(fname +
"_cg_convergence.csv");
218 for (
auto item : history_data) {
219 ofs << i <<
", " << item <<
"\n";
Describes exact solution, , of the Volume charge (rho/) numerical experiment in two and three dimensi...
Global settings for the Volume charge (rho/) numerical experiment.
Implements the Volume charge (rho/) numerical experiment.
SolverRHO(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Solves static scalar boundary value problem.
void run()
Runs the simulation.