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_Aphi_to_Brz.hpp"
21 #include "project_Aphi_to_Hrz.hpp"
22 #include "solver.hpp"
23 
24 using namespace Misc;
25 
33 {
34 public:
35  void run()
36  {
37  if (nr_threads_max > 0)
38  MultithreadInfo::set_thread_limit(nr_threads_max);
39 
40  std::string dir = "Data/circle/";
41  std::string fname;
42 
43  std::cout << "Program: ssol-iii-axi\n"
44  << "Dimensions: " << 2 << "\n"
45  << "Writing to: " << dir << "\n";
46 
47  MainOutputTable table_H(2);
48  MainOutputTable table_B(2);
49 
50  for (unsigned int p = 1; p < 4; p++) {
51  table_H.clear();
52  table_B.clear();
53 
54  for (unsigned int r = 15; r < 19; r++) { // Calculate A'.
55  fname =
56  dir + "solution_A_p" + std::to_string(p) + "_r" + std::to_string(r);
57 
59  std::cout << "Time table A \n";
60 
61  SolverSSOLIIIAXI problem(p, 2, r, fname);
62 
63  problem.clear();
64  { // Calculate H'.
65  fname =
66  dir + "solution_H_p" + std::to_string(p) + "_r" + std::to_string(r);
67 
68  table_H.add_value("r", r);
69  table_H.add_value("p", p);
70 
72  std::cout << "Time table H \n";
73 
74  ExactSolutionSSOLIIIAXI_H exact_solution;
75 
76  ProjectAphiToHrz projector(p - 1,
77  problem.get_mapping_degree(),
78  problem.get_tria(),
79  problem.get_dof_handler(),
80  problem.get_solution(),
81  fname,
82  &exact_solution,
86  true);
87 
88  table_H.add_value("ndofs", projector.get_n_dofs());
89  table_H.add_value("ncells", projector.get_n_cells());
90  table_H.add_value("L2", projector.get_L2_norm());
91  table_H.add_value("H1", 0.0);
92  }
93  { // Calculate B'.
94  fname =
95  dir + "solution_B_p" + std::to_string(p) + "_r" + std::to_string(r);
96 
97  table_B.add_value("r", r);
98  table_B.add_value("p", p);
99 
101  std::cout << "Time table D \n";
102 
103  ExactSolutionSSOLIIIAXI_B exact_solution;
104 
105  ProjectAphiToBrz projector(p - 1,
106  problem.get_mapping_degree(),
107  problem.get_tria(),
108  problem.get_dof_handler(),
109  problem.get_solution(),
110  fname,
111  &exact_solution,
115  true);
116 
117  table_B.add_value("ndofs", projector.get_n_dofs());
118  table_B.add_value("ncells", projector.get_n_cells());
119  table_B.add_value("L2", projector.get_L2_norm() / mu_0);
120  table_B.add_value("H1", 0.0);
121  }
122  }
123  // Save convergence tables
124  std::cout << "Table H\n";
125  table_H.save(dir + "table_H_p" + std::to_string(p));
126 
127  std::cout << "Table B\n";
128  table_B.save(dir + "table_B_p" + std::to_string(p));
129  }
130  }
131 };
132 
133 int
134 main()
135 {
136  try {
137  BatchSSOLIIIAXI batch;
138  batch.run();
139  } catch (std::exception& exc) {
140  std::cerr << std::endl
141  << std::endl
142  << "----------------------------------------------------"
143  << std::endl;
144  std::cerr << "Exception on processing: " << std::endl
145  << exc.what() << std::endl
146  << "Aborting!" << std::endl
147  << "----------------------------------------------------"
148  << std::endl;
149  return 1;
150  } catch (...) {
151  std::cerr << std::endl
152  << std::endl
153  << "----------------------------------------------------"
154  << std::endl;
155  std::cerr << "Unknown exception!" << std::endl
156  << "Aborting!" << std::endl
157  << "----------------------------------------------------"
158  << std::endl;
159  return 1;
160  }
161 
162  return 0;
163 }
This is a wrap-around class. It contains the main loop of the program that implements the Axisymmetri...
Definition: main.cpp:33
Describes exact solution, , of the Axisymmetric - thick spherical coil with magnetic core (ssol-iii-a...
Describes exact solution, , of the Axisymmetric - thick spherical coil with magnetic core (ssol-iii-a...
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 Axisymmetric - thick spherical coil with magnetic core (ssol-iii-axi/) numeri...
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:141
Implements the Axisymmetric - thick spherical coil with magnetic core (ssol-iii-axi/) numerical exper...
Definition: solver.hpp:40
const Triangulation< dim > & get_tria() const
Returns a reference to triangulation.
const Vector< double > & get_solution() const
Returns a reference to the solution.
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.
unsigned int get_mapping_degree() const
Returns degree of the interpolating Lagrange polynomials used for mapping from the reference cell to ...