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 #include <deal.II/base/timer.h>
16 
17 #include <iostream>
18 #include <string>
19 
20 #include "misc.hpp"
21 #include "project_A_to_B.hpp"
22 #include "project_T_to_J.hpp"
23 #include "solver.hpp"
24 
25 #include "settings.hpp"
26 
27 using namespace Misc;
28 using namespace StaticVectorSolver;
29 
36 {
37 public:
38  void run()
39  {
40  if (nr_threads_max > 0)
41  MultithreadInfo::set_thread_limit(nr_threads_max);
42 
43  std::string dir = "Data/sphere/";
44  std::string fname;
45 
46  std::cout << "Program: ssol-ii\n"
47  << "Dimensions: "
48  << "3\n"
49  << "Writing to: " << dir << "\n";
50 
51  MainOutputTable table_J(3);
52  MainOutputTable table_B(3);
53 
54  for (unsigned int p = 0; p < 3; p++) {
55  table_J.clear();
56  table_B.clear();
57 
58  for (unsigned int r = 6; r < 10; r++) {
59  table_J.add_value("r", r);
60  table_J.add_value("p", p);
61 
62  table_B.add_value("r", r);
63  table_B.add_value("p", p);
64 
65  // Stage 0 -------------------------------------------------------------
66  std::cout << "Stage 0: solving for T ...\n";
67 
68  fname =
69  dir + "solution_T_p" + std::to_string(p) + "_r" + std::to_string(r);
70 
72  std::cout << "Time table T \n";
73 
74  SolverSSOLII_T stage0(p, 2, r, fname);
75 
76  stage0.clear();
77 
78  // Stage 1 -------------------------------------------------------------
79  std::cout << "Stage 1: projecting T in H(curl) to Jf in H(div) ...\n";
80 
81  fname =
82  dir + "solution_J_p" + std::to_string(p) + "_r" + std::to_string(r);
83 
85  std::cout << "Time table J \n";
86 
87  ExactSolutionSSOLII_Jf exact_solution_Jf;
88 
89  ProjectTtoJ<1> stage1(p,
90  2,
91  stage0.get_tria(),
92  stage0.get_dof_handler(),
93  stage0.get_solution(),
94  fname,
95  &exact_solution_Jf,
99  true);
100 
101  stage1.clear();
102 
103  table_J.add_value("ndofs", stage1.get_n_dofs());
104  table_J.add_value("ncells", stage1.get_n_cells());
105  table_J.add_value("L2", stage1.get_L2_norm());
106  table_J.add_value("H1", 0.0);
107 
108  // Stage 2 -------------------------------------------------------------
109  std::cout << "Stage 2: solving for A ...\n";
110 
111  fname =
112  dir + "solution_A_p" + std::to_string(p) + "_r" + std::to_string(r);
113 
115  std::cout << "Time table A \n";
116 
117  SolverSSOLII_A stage2(p,
118  2,
119  r,
120  stage0.get_tria(),
121  stage0.get_dof_handler(),
122  stage0.get_solution(),
123  fname);
124 
125  stage2.clear();
126 
127  // Stage 3 -------------------------------------------------------------
128  std::cout << "Stage 3: projecting A in H(curl) to B in H(div) ...\n";
129 
130  fname =
131  dir + "solution_B_p" + std::to_string(p) + "_r" + std::to_string(r);
132 
134  std::cout << "Time table B \n";
135 
136  ExactSolutionSSOLII_B exact_solution_B;
137 
138  ProjectAtoB<3> stage3(p,
139  2,
140  stage0.get_tria(),
141  stage2.get_dof_handler(),
142  stage2.get_solution(),
143  fname,
144  &exact_solution_B,
148  true);
149 
150  table_B.add_value("ndofs", stage3.get_n_dofs());
151  table_B.add_value("ncells", stage3.get_n_cells());
152  table_B.add_value("L2", stage3.get_L2_norm());
153  table_B.add_value("H1", 0.0);
154  }
155  // Saving convergence tables
156  std::cout << "Table J\n";
157  table_J.save(dir + "table_J_p" + std::to_string(p));
158 
159  std::cout << "Table B\n";
160  table_B.save(dir + "table_B_p" + std::to_string(p));
161  }
162  }
163 };
164 
165 int
166 main()
167 {
168  try {
169  BatchSSOLII batch;
170  batch.run();
171  } catch (std::exception& exc) {
172  std::cerr << std::endl
173  << std::endl
174  << "----------------------------------------------------"
175  << std::endl;
176  std::cerr << "Exception on processing: " << std::endl
177  << exc.what() << std::endl
178  << "Aborting!" << std::endl
179  << "----------------------------------------------------"
180  << std::endl;
181  return 1;
182  } catch (...) {
183  std::cerr << std::endl
184  << std::endl
185  << "----------------------------------------------------"
186  << std::endl;
187  std::cerr << "Unknown exception!" << std::endl
188  << "Aborting!" << std::endl
189  << "----------------------------------------------------"
190  << std::endl;
191  return 1;
192  }
193  return 0;
194 }
This is a wrap-around class. It contains the main loop of the program that implements the Thick spher...
Definition: main.cpp:36
Describes the exact solution, , of the Thick spherical coil (ssol-ii/) numerical experiment.
Describes the exact solution, , of the Thick spherical coil (ssol-ii/) 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 Thick spherical coil (ssol-ii/) numerical experiment.
Definition: settings.hpp:25
const bool print_time_tables
If set to true, the program will print time tables on the screen.
Definition: settings.hpp:99
Solves for the magnetic vector potential, , in the Thick spherical coil (ssol-ii/) numerical experime...
Definition: solver.hpp:101
Solves for the current vector potential, , in the Thick spherical coil (ssol-ii/) numerical experimen...
Definition: solver.hpp:44
Computes the magnetic field as a curl of the magnetic vector potential, , i.e., .
unsigned int get_n_cells() const
Returns the number of active cells in the mesh.
unsigned int get_n_dofs() const
Returns the total amount of the degrees of freedom.
void clear()
Releases computer memory associated with system matrix and right-hand side.
Computes the volume free-current density as a curl of the current vector potential,...
const DoFHandler< dim > & get_dof_handler() const
Returns a reference to dof handler.
const Vector< double > & get_solution() const
Returns a reference to solution.
const Triangulation< dim > & get_tria() const
Returns a reference to triangulation.
void clear()
Releases computer memory associated with system matrix and right-hand side.
void clear()
Releases computer memory associated with system matrix and right-hand side.
const Vector< double > & get_solution() const
Returns a reference to solution.
const DoFHandler< dim > & get_dof_handler() const
Returns a reference to a dof handler.