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-i-axi\n"
44  << "Dimensions: " << 2 << "\n"
45  << "Writing to: " << dir << "\n";
46 
47  MainOutputTable table_A(2);
48  MainOutputTable table_B(2);
49 
50  for (unsigned int p = 1; p < 4; p++) {
51  table_A.clear();
52  table_B.clear();
53 
54  for (unsigned int r = 15; r < 19;
55  r++) { // Calculate magnetic vector potential.
56  fname =
57  dir + "solution_A_p" + std::to_string(p) + "_r" + std::to_string(r);
58 
59  table_A.add_value("r", r);
60  table_A.add_value("p", p);
61 
63  std::cout << "Time table A \n";
64 
65  SolverSSOLIAXI problem(p, 2, r, fname);
66  table_A.add_value("ndofs", problem.get_n_dofs());
67  table_A.add_value("ncells", problem.get_n_cells());
68  table_A.add_value("L2", problem.get_L2_norm() / mu_0);
69  table_A.add_value("H1", problem.get_H1_norm() / mu_0);
70 
71  problem.clear();
72  { // Calculate the magnetic field.
73  fname =
74  dir + "solution_B_p" + std::to_string(p) + "_r" + std::to_string(r);
75 
76  table_B.add_value("r", r);
77  table_B.add_value("p", p);
78 
80  std::cout << "Time table B \n";
81 
82  ExactSolutionSSOLIAXI_B exact_solution;
83 
84  ProjectAphiToBrz projector(p - 1,
85  problem.get_mapping_degree(),
86  problem.get_tria(),
87  problem.get_dof_handler(),
88  problem.get_solution(),
89  fname,
90  &exact_solution,
94  true);
95 
96  table_B.add_value("ndofs", projector.get_n_dofs());
97  table_B.add_value("ncells", projector.get_n_cells());
98  table_B.add_value("L2", projector.get_L2_norm() / mu_0);
99  table_B.add_value("H1", 0.0);
100  }
101  }
102  // Saving convergence tables
103  std::cout << "Table A\n";
104  table_A.save(dir + "table_A_p" + std::to_string(p));
105 
106  std::cout << "Table B\n";
107  table_B.save(dir + "table_B_p" + std::to_string(p));
108  }
109  }
110 };
111 
112 int
113 main()
114 {
115  try {
116  BatchSSOLIAXI batch;
117  batch.run();
118  } catch (std::exception& exc) {
119  std::cerr << std::endl
120  << std::endl
121  << "----------------------------------------------------"
122  << std::endl;
123  std::cerr << "Exception on processing: " << std::endl
124  << exc.what() << std::endl
125  << "Aborting!" << std::endl
126  << "----------------------------------------------------"
127  << std::endl;
128  return 1;
129  } catch (...) {
130  std::cerr << std::endl
131  << std::endl
132  << "----------------------------------------------------"
133  << std::endl;
134  std::cerr << "Unknown exception!" << std::endl
135  << "Aborting!" << std::endl
136  << "----------------------------------------------------"
137  << std::endl;
138  return 1;
139  }
140 
141  return 0;
142 }
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 - thin spherical coil (ssol-i-axi) numerical experime...
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 - thin spherical coil (ssol-i-axi) numerical experiment.
Definition: settings.hpp:32
const bool print_time_tables
If set to true, the program will print the time tables on the screen.
Definition: settings.hpp:103
Implements the Axisymmetric - thin spherical coil (ssol-i-axi) numerical experiment.
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.
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 ...