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 #if DOMAIN__ == 0
44  const unsigned int mapping_degree = 1;
45  const std::string dir = "Data/cube/";
46 #endif
47 
48 #if DOMAIN__ == 1
49  const unsigned int mapping_degree = 2;
50  const std::string dir = "Data/sphere/";
51 #endif
52 
53  std::cout << "Program: mss-vt-i\n"
54  << "Dimensions: "
55  << "3\n"
56  << "Domain type: " << DOMAIN__ << "\n"
57  << "Writing to: " << dir << "\n";
58 
59  MainOutputTable table_J(3);
60  MainOutputTable table_B(3);
61 
62  for (unsigned int p = 0; p < 3; p++) {
63  table_J.clear();
64  table_B.clear();
65 
66  for (unsigned int r = 5; r < 9; r++) {
67  table_J.add_value("r", r);
68  table_J.add_value("p", p);
69 
70  table_B.add_value("r", r);
71  table_B.add_value("p", p);
72 
73  // Stage 0 -------------------------------------------------------------
74  std::cout << "Stage 0: solving for T ...\n";
75 
76  SolverMMSVTI_T stage0(p,
77  mapping_degree,
78  r,
79  dir + "solution_T_p" + std::to_string(p) + "_r" +
80  std::to_string(r));
81 
82  stage0.clear();
83 
84  // Stage 1 -------------------------------------------------------------
85  std::cout << "Stage 1: projecting T in H(curl) to Jf in H(div) ...\n";
86 
87  ExactSolutionMMSVTI_Jf exact_solution_Jf;
88 
89  ProjectTtoJ stage1(p,
90  mapping_degree,
91  stage0.get_tria(),
92  stage0.get_dof_handler(),
93  stage0.get_solution(),
94  dir + "solution_J_p" + std::to_string(p) + "_r" +
95  std::to_string(r),
96  &exact_solution_Jf,
100  true);
101 
102  table_J.add_value("ndofs", stage0.get_n_dofs());
103  table_J.add_value("ncells", stage0.get_n_cells());
104  table_J.add_value("L2", stage1.get_L2_norm());
105  table_J.add_value("H1", 0.0);
106 
107  stage1.clear();
108 
109  // Stage 2 -------------------------------------------------------------
110  std::cout << "Stage 2: solving for A ...\n";
111 
112  SolverMMSVTI_A stage2(p,
113  mapping_degree,
114  stage0.get_tria(),
115  stage0.get_dof_handler(),
116  stage0.get_solution(),
117  r,
118  dir + "solution_A_p" + std::to_string(p) + "_r" +
119  std::to_string(r));
120 
121  stage2.clear();
122 
123  // Stage 3 -------------------------------------------------------------
124  std::cout << "Stage 3: projecting A in H(curl) to B in H(div) ...\n";
125 
126  ExactSolutionMMSVTI_B exact_solution_B;
127 
128  ProjectAtoB stage3(p,
129  mapping_degree,
130  stage0.get_tria(),
131  stage2.get_dof_handler(),
132  stage2.get_solution(),
133  dir + "solution_B_p" + std::to_string(p) + "_r" +
134  std::to_string(r),
135  &exact_solution_B,
139  true);
140 
141  table_B.add_value("ndofs", stage0.get_n_dofs());
142  table_B.add_value("ncells", stage0.get_n_cells());
143  table_B.add_value("L2", stage3.get_L2_norm());
144  table_B.add_value("H1", 0.0);
145  }
146 
147  std::cout << "Table J \n";
148  table_J.save(dir + "main_table_J_p" + std::to_string(p));
149 
150  std::cout << "Table B \n";
151  table_B.save(dir + "main_table_B_p" + std::to_string(p));
152  }
153  }
154 };
155 
156 int
157 main()
158 {
159  try {
160  BatchMMSVTI batch;
161  batch.run();
162  } catch (std::exception& exc) {
163  std::cerr << std::endl
164  << std::endl
165  << "----------------------------------------------------"
166  << std::endl;
167  std::cerr << "Exception on processing: " << std::endl
168  << exc.what() << std::endl
169  << "Aborting!" << std::endl
170  << "----------------------------------------------------"
171  << std::endl;
172  return 1;
173  } catch (...) {
174  std::cerr << std::endl
175  << std::endl
176  << "----------------------------------------------------"
177  << std::endl;
178  std::cerr << "Unknown exception!" << std::endl
179  << "Aborting!" << std::endl
180  << "----------------------------------------------------"
181  << std::endl;
182  return 1;
183  }
184  return 0;
185 }
This is a wrap-around class. It contains the main loop of the program that implements the Method of m...
Definition: main.cpp:36
Describes exact solution, , of the Method of manufactured solutions, vector potential (mms-vt-i/) num...
Describes exact solution, , of the Method of manufactured solutions, vector potential (mms-vt-i/) num...
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 Method of manufactured solutions, vector potential (mms-vt-i/) numerical expe...
Definition: settings.hpp:26
Implements the solver for magnetic vector potential, , in the Method of manufactured solutions,...
Definition: solver.hpp:104
Implements the solver for current vector potential, , in the Method of manufactured solutions,...
Definition: solver.hpp:46
Computes the magnetic field as a curl of the magnetic vector potential, , i.e., .
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,...
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.
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.