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