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  std::string dir =
40  (CYLINDER__ == 1) ? "Data/cylinder-axi/" : "Data/sphere-axi/";
41  std::string fname;
42 
43  unsigned int mapping_degree = (CYLINDER__ == 1) ? 1 : 2;
44 
45  std::cout << "Program: int-axi\n"
46  << "Dimensions: " << 2 << "\n"
47  << "Condition: Cylinder = " << CYLINDER__ << "\n"
48  << "Mapping degree = " << mapping_degree << "\n"
49  << "Writing to: " << dir << "\n";
50 
51  MainOutputTable table_PHI(2);
52  MainOutputTable table_E(2);
53  MainOutputTable table_D(2);
54 
55  for (unsigned int p = 1; p < 4; p++) {
56  table_PHI.clear();
57  table_E.clear();
58  table_D.clear();
59 
60 #if CYLINDER__ == 1
61  for (unsigned int r = 9; r < 13; r++)
62 #endif
63 #if CYLINDER__ == 0
64  for (unsigned int r = 5; r < 9; r++)
65 #endif
66  {
67  fname = dir + "solution_PHI_p" + std::to_string(p) + "_r" +
68  std::to_string(r);
69 
70  table_PHI.add_value("r", r);
71  table_PHI.add_value("p", p);
72 
74  std::cout << "Time table PHI \n";
75 
76  SolverINTAXI<static_cast<bool>(CYLINDER__)> problem(
77  p, mapping_degree, r, fname);
78  table_PHI.add_value("ndofs", problem.get_n_dofs());
79  table_PHI.add_value("ncells", problem.get_n_cells());
80  table_PHI.add_value("L2", problem.get_L2_norm());
81  table_PHI.add_value("H1", problem.get_H1_norm());
82 
83  problem.clear();
84 
85  { // Calculating electrostatic field
86  fname = dir + "solution_E" + "_p" + std::to_string(p) + "_r" +
87  std::to_string(r);
88 
89  table_E.add_value("r", r);
90  table_E.add_value("p", p);
91 
93  std::cout << "Time table E \n";
94 
95  ExactSolutionINTAXI_E<CYLINDER__> exact_solution;
96 
97  ProjectPHItoE projector(p - 1,
98  problem.get_mapping_degree(),
99  problem.get_tria(),
100  problem.get_dof_handler(),
101  problem.get_solution(),
102  fname,
103  &exact_solution,
104  true,
108  true);
109 
110  table_E.add_value("ndofs", projector.get_n_dofs());
111  table_E.add_value("ncells", projector.get_n_cells());
112  table_E.add_value("L2", projector.get_L2_norm());
113  table_E.add_value("H1", 0.0);
114  }
115  { // Calculating displacement
116  fname = dir + "solution_D" + "_p" + std::to_string(p) + "_r" +
117  std::to_string(r);
118 
119  table_D.add_value("r", r);
120  table_D.add_value("p", p);
121 
123  std::cout << "Time table D \n";
124 
125  ExactSolutionINTAXI_D<CYLINDER__> exact_solution;
126 
127  ProjectPHItoD projector(p - 1,
128  problem.get_mapping_degree(),
129  problem.get_tria(),
130  problem.get_dof_handler(),
131  problem.get_solution(),
132  fname,
133  &exact_solution,
134  true,
138  true);
139 
140  table_D.add_value("ndofs", projector.get_n_dofs());
141  table_D.add_value("ncells", projector.get_n_cells());
142  table_D.add_value("L2", projector.get_L2_norm() / ep_0);
143  table_D.add_value("H1", 0.0);
144  }
145  }
146  // Saving convergence tables
147  std::cout << "Table PHI\n";
148  table_PHI.save(dir + "table_PHI_p" + std::to_string(p));
149 
150  std::cout << "Table E\n";
151  table_E.save(dir + "table_E_p" + std::to_string(p));
152 
153  std::cout << "Table D\n";
154  table_D.save(dir + "table_D_p" + std::to_string(p));
155  }
156  }
157 };
158 
159 int
160 main()
161 {
162  try {
163  BatchINTAXI batch;
164  batch.run();
165  } catch (std::exception& exc) {
166  std::cerr << std::endl
167  << std::endl
168  << "----------------------------------------------------"
169  << std::endl;
170  std::cerr << "Exception on processing: " << std::endl
171  << exc.what() << std::endl
172  << "Aborting!" << std::endl
173  << "----------------------------------------------------"
174  << std::endl;
175  return 1;
176  } catch (...) {
177  std::cerr << std::endl
178  << std::endl
179  << "----------------------------------------------------"
180  << std::endl;
181  std::cerr << "Unknown exception!" << std::endl
182  << "Aborting!" << std::endl
183  << "----------------------------------------------------"
184  << std::endl;
185  return 1;
186  }
187 
188  return 0;
189 }
This is a wrap-around class. It contains the main loop of the program that implements the Axisymmetri...
Definition: main.cpp:32
Describes exact solution, , of the Axisymmetric - interface between dielectrics (int-axi/) numerical ...
Describes exact solution, , of the Axisymmetric - interface between dielectrics (int-axi/) numerical ...
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 - interface between dielectrics (int-axi/) 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:96
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 - interface between dielectrics (int-axi/) numerical experiment.
Definition: solver.hpp:42