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 
37 {
38 public:
39  void run()
40  {
41  if (nr_threads_max > 0)
42  MultithreadInfo::set_thread_limit(nr_threads_max);
43 
44  std::string dir = "Data/sphere/";
45  std::string fname;
46 
47  std::cout << "Program: ssol-iii\n"
48  << "Dimensions: "
49  << "3\n"
50  << "Writing to: " << dir << "\n";
51 
52  MainOutputTable table_J(3);
53  MainOutputTable table_B(3);
54 
55  for (unsigned int p = 0; p < 3; p++) {
56  table_J.clear();
57  table_B.clear();
58 
59  for (unsigned int r = 6; r < 10; r++) {
60  table_J.add_value("r", r);
61  table_J.add_value("p", p);
62 
63  table_B.add_value("r", r);
64  table_B.add_value("p", p);
65 
66  // Stage 0 -------------------------------------------------------------
67  std::cout << "Stage 0: solving for T ...\n";
68 
69  fname =
70  dir + "solution_T_p" + std::to_string(p) + "_r" + std::to_string(r);
71 
73  std::cout << "Time table T \n";
74 
75  SolverSSOLIII_T stage0(p, 2, r, fname);
76 
77  stage0.clear();
78 
79  // Stage 1 -------------------------------------------------------------
80  std::cout << "Stage 1: projecting T in H(curl) to Jf in H(div) ...\n";
81 
82  fname =
83  dir + "solution_J_p" + std::to_string(p) + "_r" + std::to_string(r);
84 
86  std::cout << "Time table J \n";
87 
88  ExactSolutionSSOLIII_Jf exact_solution_Jf;
89 
90  ProjectTtoJ<1> stage1(p,
91  2,
92  stage0.get_tria(),
93  stage0.get_dof_handler(),
94  stage0.get_solution(),
95  fname,
96  &exact_solution_Jf,
100  true);
101 
102  stage1.clear();
103 
104  table_J.add_value("ndofs", stage1.get_n_dofs());
105  table_J.add_value("ncells", stage1.get_n_cells());
106  table_J.add_value("L2", stage1.get_L2_norm());
107  table_J.add_value("H1", 0.0);
108 
109  // Stage 2 -------------------------------------------------------------
110  std::cout << "Stage 2: solving for A ...\n";
111 
112  fname =
113  dir + "solution_A_p" + std::to_string(p) + "_r" + std::to_string(r);
114 
116  std::cout << "Time table A \n";
117 
118  SolverSSOLIII_A stage2(p,
119  2,
120  r,
121  stage0.get_tria(),
122  stage0.get_dof_handler(),
123  stage0.get_solution(),
124  fname);
125 
126  stage2.clear();
127 
128  // Stage 3 -------------------------------------------------------------
129  std::cout << "Stage 3: projecting A in H(curl) to B in H(div) ...\n";
130 
131  fname =
132  dir + "solution_B_p" + std::to_string(p) + "_r" + std::to_string(r);
133 
135  std::cout << "Time table B \n";
136 
137  ExactSolutionSSOLIII_B exact_solution_B;
138 
139  ProjectAtoB<3> stage3(p,
140  2,
141  stage0.get_tria(),
142  stage2.get_dof_handler(),
143  stage2.get_solution(),
144  fname,
145  &exact_solution_B,
149  true);
150 
151  table_B.add_value("ndofs", stage3.get_n_dofs());
152  table_B.add_value("ncells", stage3.get_n_cells());
153  table_B.add_value("L2", stage3.get_L2_norm());
154  table_B.add_value("H1", 0.0);
155  }
156  // Saving convergence tables
157  std::cout << "Table J\n";
158  table_J.save(dir + "table_J_p" + std::to_string(p));
159 
160  std::cout << "Table B\n";
161  table_B.save(dir + "table_B_p" + std::to_string(p));
162  }
163  }
164 };
165 
166 int
167 main()
168 {
169  try {
170  BatchSSOLIII batch;
171  batch.run();
172  } catch (std::exception& exc) {
173  std::cerr << std::endl
174  << std::endl
175  << "----------------------------------------------------"
176  << std::endl;
177  std::cerr << "Exception on processing: " << std::endl
178  << exc.what() << std::endl
179  << "Aborting!" << std::endl
180  << "----------------------------------------------------"
181  << std::endl;
182  return 1;
183  } catch (...) {
184  std::cerr << std::endl
185  << std::endl
186  << "----------------------------------------------------"
187  << std::endl;
188  std::cerr << "Unknown exception!" << std::endl
189  << "Aborting!" << std::endl
190  << "----------------------------------------------------"
191  << std::endl;
192  return 1;
193  }
194  return 0;
195 }
This is a wrap-around class. It contains the main loop of the program that implements the Spherical c...
Definition: main.cpp:37
Describes exact solution, , of the Spherical coil with magnetic core (ssol-iii/) numerical experiment...
Describes exact solution, , of the Spherical coil with magnetic core (ssol-iii/) 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 Spherical coil with magnetic core (ssol-iii/) 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:133
Solves for the magnetic vector potential, , in the Spherical coil with magnetic core (ssol-iii/) nume...
Definition: solver.hpp:104
Solves for the current vector potential, , in the Spherical coil with magnetic core (ssol-iii/) numer...
Definition: solver.hpp:46
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.