Logbook  (07-04-2025)
Static problems
main.cpp
1 /******************************************************************************
2  * Copyright (C) Siarhei Uzunbajakau, 2023.
3  *
4  * This program is free software. You can use, modify, and redistribute it under
5  * the terms of the GNU Lesser General Public License as published by the Free
6  * Software Foundation, either version 3 or (at your option) any later version.
7  * This program is distributed without any warranty.
8  *
9  * Refer to COPYING.LESSER for more details.
10  ******************************************************************************/
11 
12 #define BOOST_ALLOW_DEPRECATED_HEADERS
13 
14 #include "deal.II/base/multithread_info.h"
15 
16 #include <iostream>
17 #include <string>
18 
19 #include "misc.hpp"
20 #include "project_Az_to_Bxy.hpp"
21 #include "project_Az_to_Hxy.hpp"
22 #include "solver.hpp"
23 
24 using namespace Misc;
25 
31 class BatchMWR : public SettingsMWR
32 {
33 public:
34  void run()
35  {
36  if (nr_threads_max > 0)
37  MultithreadInfo::set_thread_limit(nr_threads_max);
38 
39  std::string dir = "Data/circle/";
40  std::string fname;
41 
42  std::cout << "Program: mwr\n"
43  << "Dimensions: " << 2 << "\n"
44  << "Writing to: " << dir << "\n";
45 
46  MainOutputTable table_A(2);
47  MainOutputTable table_H(2);
48  MainOutputTable table_B(2);
49 
50  for (unsigned int p = 1; p < 4; p++) {
51  table_A.clear();
52  table_H.clear();
53  table_B.clear();
54 
55  for (unsigned int r = 5; r < 9;
56  r++) { // Calculate magnetic vector potential.
57  fname =
58  dir + "solution_A_p" + std::to_string(p) + "_r" + std::to_string(r);
59 
60  table_A.add_value("r", r);
61  table_A.add_value("p", p);
62 
64  std::cout << "Time table A \n";
65 
66  SolverMWR problem(p, 2, r, fname);
67  table_A.add_value("ndofs", problem.get_n_dofs());
68  table_A.add_value("ncells", problem.get_n_cells());
69  table_A.add_value("L2", problem.get_L2_norm() / mu_0);
70  table_A.add_value("H1", problem.get_H1_norm() / mu_0);
71 
72  problem.clear();
73  { // Calculate auxiliary H-field.
74  fname =
75  dir + "solution_H_p" + std::to_string(p) + "_r" + std::to_string(r);
76 
77  table_H.add_value("r", r);
78  table_H.add_value("p", p);
79 
81  std::cout << "Time table H \n";
82 
83  ExactSolutionMWR_H exact_solution;
84 
85  ProjectAzToHxy projector(p - 1,
86  problem.get_mapping_degree(),
87  problem.get_tria(),
88  problem.get_dof_handler(),
89  problem.get_solution(),
90  fname,
91  &exact_solution,
95  true);
96 
97  table_H.add_value("ndofs", projector.get_n_dofs());
98  table_H.add_value("ncells", projector.get_n_cells());
99  table_H.add_value("L2", projector.get_L2_norm());
100  table_H.add_value("H1", 0.0);
101  }
102  { // Calculate the magnetic field.
103  fname =
104  dir + "solution_B_p" + std::to_string(p) + "_r" + std::to_string(r);
105 
106  table_B.add_value("r", r);
107  table_B.add_value("p", p);
108 
110  std::cout << "Time table D \n";
111 
112  ExactSolutionMWR_B exact_solution;
113 
114  ProjectAzToBxy projector(p - 1,
115  problem.get_mapping_degree(),
116  problem.get_tria(),
117  problem.get_dof_handler(),
118  problem.get_solution(),
119  fname,
120  &exact_solution,
124  true);
125 
126  table_B.add_value("ndofs", projector.get_n_dofs());
127  table_B.add_value("ncells", projector.get_n_cells());
128  table_B.add_value("L2", projector.get_L2_norm() / mu_0);
129  table_B.add_value("H1", 0.0);
130  }
131  }
132  // Saving convergence tables
133  std::cout << "Table A\n";
134  table_A.save(dir + "table_A_p" + std::to_string(p));
135 
136  std::cout << "Table H\n";
137  table_H.save(dir + "table_H_p" + std::to_string(p));
138 
139  std::cout << "Table B\n";
140  table_B.save(dir + "table_B_p" + std::to_string(p));
141  }
142  }
143 };
144 
145 int
146 main()
147 {
148  try {
149  BatchMWR batch;
150  batch.run();
151  } catch (std::exception& exc) {
152  std::cerr << std::endl
153  << std::endl
154  << "----------------------------------------------------"
155  << std::endl;
156  std::cerr << "Exception on processing: " << std::endl
157  << exc.what() << std::endl
158  << "Aborting!" << std::endl
159  << "----------------------------------------------------"
160  << std::endl;
161  return 1;
162  } catch (...) {
163  std::cerr << std::endl
164  << std::endl
165  << "----------------------------------------------------"
166  << std::endl;
167  std::cerr << "Unknown exception!" << std::endl
168  << "Aborting!" << std::endl
169  << "----------------------------------------------------"
170  << std::endl;
171  return 1;
172  }
173 
174  return 0;
175 }
This is a wrap-around class. It contains the main loop of the program that implements the Magnetic wi...
Definition: main.cpp:32
Describes exact solution, , of the Magnetic wire (mwr/) numerical experiment.
Describes exact solution, , of the Magnetic wire (mwr/) numerical experiment.
The convergence table used in multiple numerical experiments.
Definition: misc.hpp:25
void save(std::string fname)
Saves the data in text and tex formats, and prints the data on screen.
Definition: misc.cpp:47
const bool project_exact_solution
If set to true, the program will project the exact solution.
Definition: settings.hpp:80
const bool log_cg_convergence
If set to true, saves the residual at each iteration of the CG solver. The names of the files fit the...
Definition: settings.hpp:89
const bool print_time_tables
If set to true, the program will print time tables on the screen.
Definition: settings.hpp:70
Global settings for the Magnetic wire (mwr/) numerical experiment.
Definition: settings.hpp:26
const bool print_time_tables
If set to true, the program will print time tables on the screen.
Definition: settings.hpp:109
Implements the Magnetic wire (mwr/) numerical experiment.
Definition: solver.hpp:39
const Triangulation< dim > & get_tria() const
Returns a reference to triangulation.
const Vector< double > & get_solution() const
Returns a reference to the solution.
double get_L2_norm() const
Returns error norm.
void clear()
Releases computer memory associated with the system matrix and right-hand side.
const DoFHandler< dim > & get_dof_handler() const
Returns a reference to dof handler.
double get_H1_norm() const
Returns error norm.
unsigned int get_n_dofs() const
Returns the total amount of the degrees of freedom.
unsigned int get_n_cells() const
Returns the number of active cells in the mesh.
unsigned int get_mapping_degree() const
Returns degree of the interpolating Lagrange polynomials used for mapping from the reference cell to ...