12 #ifndef SolverSLDI_H__
13 #define SolverSLDI_H__
15 #define BOOST_ALLOW_DEPRECATED_HEADERS
17 #include <deal.II/base/function.h>
18 #include <deal.II/grid/grid_in.h>
19 #include <deal.II/grid/grid_out.h>
20 #include <deal.II/grid/grid_tools.h>
21 #include <deal.II/grid/manifold_lib.h>
23 #include <deal.II/numerics/fe_field_function.h>
25 #include "exact_solution.hpp"
26 #include "settings.hpp"
27 #include "static_scalar_solver.hpp"
29 #define TMR(__name) TimerOutput::Scope timer_section(timer, __name)
31 using namespace StaticScalarSolver;
33 #pragma GCC diagnostic push
34 #pragma GCC diagnostic ignored "-Wunused-parameter"
44 ,
public Function<dim>
49 virtual double value(
const Point<dim>& r,
50 const unsigned int component = 0)
const override final;
53 #pragma GCC diagnostic pop
83 unsigned int mapping_degree,
94 project_exact_solution,
105 const unsigned int r;
106 const std::string fname;
111 SphericalManifold<dim> sphere;
113 virtual void make_mesh() override final;
114 virtual
void fill_dirichlet_stack() override final;
115 virtual
void solve() override final;
122 switch (type_of_bc) {
126 case DirichletNeumann:
142 std::string fname_mesh =
143 (dim == 2) ?
"../../gmsh/data/square_r" + std::to_string(r) +
".msh"
144 :
"../../gmsh/data/cube_r" + std::to_string(r) +
".msh";
146 std::ifstream ifs(fname_mesh);
147 gridin.read_msh(ifs);
153 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; f++) {
154 if (((type_of_bc == DirichletOnly) || (type_of_bc == Exact)) &&
155 cell->face(f)->at_boundary())
156 cell->face(f)->set_boundary_id(bid);
159 if (cell->center().norm() < a) {
161 cell->set_material_id(mid_1);
162 }
else if (cell->center().norm() < b) {
164 cell->set_material_id(mid_2);
167 cell->set_material_id(mid_3);
170 if ((cell->center().norm() > rd1) && (cell->center().norm() < b))
171 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; f++) {
172 double dif_norm = 0.0;
173 for (
unsigned int v = 1; v < GeometryInfo<dim>::vertices_per_face; v++)
174 dif_norm += std::abs(cell->face(f)->vertex(0).norm() -
175 cell->face(f)->vertex(v).norm());
178 cell->face(f)->set_all_manifold_ids(1);
194 if (log_cg_convergence)
195 control.enable_history_data();
197 GrowingVectorMemory<Vector<double>> memory;
198 SolverCG<Vector<double>> cg(control, memory);
200 PreconditionJacobi<SparseMatrix<double>> preconditioner;
210 if (log_cg_convergence) {
211 const std::vector<double> history_data = control.get_history_data();
213 std::ofstream ofs(fname +
"_cg_convergence.csv");
216 for (
auto item : history_data) {
217 ofs << i <<
", " << item <<
"\n";
The Dirichlet boundary condition of the Magnetostatic shield - 1 (sld-i/) numerical experiment....
Describes exact solution, , of the Magnetostatic shield - 1 (sld-i/) numerical experiment in two and ...
Global settings for the Magnetostatic shield - 1 (sld-i/) numerical experiment.
Implements the solver that solves for in the Magnetostatic shield - 1 (sld-i/) numerical experiment.
SolverSLDI(unsigned int p, unsigned int mapping_degree, unsigned int r, std::string fname)
Solves static scalar boundary value problem.
void run()
Runs the simulation.