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_PHI_to_D.hpp"
21 #include "project_PHI_to_E.hpp"
22 #include "solver.hpp"
23 
24 using namespace Misc;
25 
32 {
33 public:
34  void run()
35  {
36  if (nr_threads_max > 0)
37  MultithreadInfo::set_thread_limit(nr_threads_max);
38 
39  Assert(DIMENSION__ > 1, ExcInternalError());
40  Assert(DIMENSION__ < 4, ExcInternalError());
41 
42  std::string dir =
43  (DIMENSION__ == 2) ? "Data/cylinder2d/" : "Data/cylinder3d/";
44  std::string fname;
45 
46  std::cout << "Program: mms-axi\n"
47  << "Dimensions: " << DIMENSION__ << "\n"
48  << "Writing to: " << dir << "\n";
49 
50  MainOutputTable table_PHI(DIMENSION__);
51  MainOutputTable table_E(DIMENSION__);
52  MainOutputTable table_D(DIMENSION__);
53 
54  bool axisymmetric;
55  for (unsigned int p = 1; p < 4; p++) {
56  table_PHI.clear();
57  table_E.clear();
58  table_D.clear();
59 
60 #if DIMENSION__ == 2
61  axisymmetric = true;
62  for (unsigned int r = 27; r < 31; r++)
63 #endif
64 #if DIMENSION__ == 3
65  axisymmetric = false;
66  for (unsigned int r = 5; r < 9; r++)
67 #endif
68  {
69  fname =
70  dir + "solution_PHI_p" + std::to_string(p) + "_r" + std::to_string(r);
71 
72  table_PHI.add_value("r", r);
73  table_PHI.add_value("p", p);
74 
75  SolverMMSAXI<DIMENSION__> problem(p, r, fname);
76 
77  table_PHI.add_value("ndofs", problem.get_n_dofs());
78  table_PHI.add_value("ncells", problem.get_n_cells());
79  table_PHI.add_value("L2", problem.get_L2_norm());
80  table_PHI.add_value("H1", problem.get_H1_norm());
81 
82  problem.clear();
83 
84  {
85  fname = dir + "solution_E" + "_p" + std::to_string(p) + "_r" +
86  std::to_string(r);
87 
88  table_E.add_value("r", r);
89  table_E.add_value("p", p);
90 
92 
93  ProjectPHItoE<DIMENSION__> projector(p - 1,
94  1,
95  problem.get_tria(),
96  problem.get_dof_handler(),
97  problem.get_solution(),
98  fname,
99  &exact_solution,
100  axisymmetric,
104 
105  table_E.add_value("ndofs", problem.get_n_dofs());
106  table_E.add_value("ncells", problem.get_n_cells());
107  table_E.add_value("L2", projector.get_L2_norm());
108  table_E.add_value("H1", 0.0);
109  }
110  {
111  fname = dir + "solution_D" + "_p" + std::to_string(p) + "_r" +
112  std::to_string(r);
113 
114  table_D.add_value("r", r);
115  table_D.add_value("p", p);
116 
117  ExactSolutionMMSAXI_D<DIMENSION__> exact_solution;
118 
119  ProjectPHItoD<DIMENSION__> projector(p - 1,
120  1,
121  problem.get_tria(),
122  problem.get_dof_handler(),
123  problem.get_solution(),
124  fname,
125  &exact_solution,
126  axisymmetric,
130  table_D.add_value("ndofs", problem.get_n_dofs());
131  table_D.add_value("ncells", problem.get_n_cells());
132  table_D.add_value("L2", projector.get_L2_norm() / ep_0);
133  table_D.add_value("H1", 0.0);
134  }
135  }
136 
137  std::cout << "Table PHI\n";
138  table_PHI.save(dir + "table_PHI_p" + std::to_string(p));
139 
140  std::cout << "Table E\n";
141  table_E.save(dir + "table_E_p" + std::to_string(p));
142 
143  std::cout << "Table D\n";
144  table_D.save(dir + "table_D_p" + std::to_string(p));
145  }
146  }
147 };
148 
149 int
150 main()
151 {
152  try {
153  BatchMMSAXI batch;
154  batch.run();
155  } catch (std::exception& exc) {
156  std::cerr << std::endl
157  << std::endl
158  << "----------------------------------------------------"
159  << std::endl;
160  std::cerr << "Exception on processing: " << std::endl
161  << exc.what() << std::endl
162  << "Aborting!" << std::endl
163  << "----------------------------------------------------"
164  << std::endl;
165  return 1;
166  } catch (...) {
167  std::cerr << std::endl
168  << std::endl
169  << "----------------------------------------------------"
170  << std::endl;
171  std::cerr << "Unknown exception!" << std::endl
172  << "Aborting!" << std::endl
173  << "----------------------------------------------------"
174  << std::endl;
175  return 1;
176  }
177  return 0;
178 }
This is a wrap-around class. It contains the main loop of the program that implements the Axisymmetri...
Definition: main.cpp:32
Describes the exact solution, , of the Axisymmetric - method of manufactured solutions (mms-axi/) num...
Describes the exact solution, , of the Axisymmetric - method of manufactured solutions (mms-axi/) 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
Global settings for the Axisymmetric - method of manufactured solutions (mms-axi/) numerical experime...
Definition: settings.hpp:26
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
Implements the Axisymmetric - method of manufactured solutions (mms-axi/) numerical experiment.
Definition: solver.hpp:41
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.