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_Tz_to_Jxy.hpp"
22 #include "solver.hpp"
23 
24 #include "settings.hpp"
25 
26 using namespace Misc;
27 using namespace StaticScalarSolver;
28 
34 class BatchCVPII : public SettingsCVPII
35 {
36 public:
37  void run()
38  {
39  if (nr_threads_max > 0)
40  MultithreadInfo::set_thread_limit(nr_threads_max);
41 
42  std::string fname;
43  std::string dir = "Data/circle-linear/";
44 
45  std::cout << "Program: cvp-ii\n"
46  << "Dimensions: 2\n"
47  << "Writing to: " << dir << "\n";
48 
49  MainOutputTable table_T(2);
50  MainOutputTable table_J(2);
51 
52  for (unsigned int p = 1; p < 4; p++) {
53  table_T.clear();
54  table_J.clear();
55 
56  for (unsigned int r = 12; r < 16; r++) {
57  std::cout << "Solving for T ...\n";
58 
59  fname =
60  dir + "solution_T_p" + std::to_string(p) + "_r" + std::to_string(r);
61 
62  table_T.add_value("r", r);
63  table_T.add_value("p", p);
64 
66  std::cout << "Time table\n";
67 
68  SolverCVPII problem(p, 2, r, fname);
69  problem.clear();
70 
71  table_T.add_value("ndofs", problem.get_n_dofs());
72  table_T.add_value("ncells", problem.get_n_cells());
73  table_T.add_value("L2", problem.get_L2_norm());
74  table_T.add_value("H1", problem.get_H1_norm());
75 
76  std::cout << "Converting T to J. \n";
77 
78  fname =
79  dir + "solution_J_p" + std::to_string(p) + "_r" + std::to_string(r);
80 
81  table_J.add_value("r", r);
82  table_J.add_value("p", p);
83 
84  ExactSolutionCVPII_Jf exact_solution;
85 
87  std::cout << "Time table\n";
88 
89  ProjectTzToJxy projector(p - 1,
90  2,
91  problem.get_tria(),
92  problem.get_dof_handler(),
93  problem.get_solution(),
94  fname,
95  &exact_solution,
99  true);
100 
101  table_J.add_value("ndofs", projector.get_n_dofs());
102  table_J.add_value("ncells", projector.get_n_cells());
103  table_J.add_value("L2", projector.get_L2_norm());
104  table_J.add_value("H1", 0.0);
105  }
106  // Saving convergence table
107  std::cout << "Table T\n";
108  table_T.save(dir + "table_T_p" + std::to_string(p));
109 
110  std::cout << "Table J\n";
111  table_J.save(dir + "table_J_p" + std::to_string(p));
112  }
113  }
114 };
115 
116 int
117 main()
118 {
119  try {
120  BatchCVPII batch;
121  batch.run();
122  } catch (std::exception& exc) {
123  std::cerr << std::endl
124  << std::endl
125  << "----------------------------------------------------"
126  << std::endl;
127  std::cerr << "Exception on processing: " << std::endl
128  << exc.what() << std::endl
129  << "Aborting!" << std::endl
130  << "----------------------------------------------------"
131  << std::endl;
132  return 1;
133  } catch (...) {
134  std::cerr << std::endl
135  << std::endl
136  << "----------------------------------------------------"
137  << std::endl;
138  std::cerr << "Unknown exception!" << std::endl
139  << "Aborting!" << std::endl
140  << "----------------------------------------------------"
141  << std::endl;
142  return 1;
143  }
144  return 0;
145 }
This is a wrap-around class. It contains the main loop of the program that implements the Current vec...
Definition: main.cpp:35
Describes the given volume free-current density, , in the Current vector potential (cvp-ii/) numerica...
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 Current vector potential (cvp-ii/) 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:88
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 Current vector potential (cvp-ii/) numerical experiment.
Definition: solver.hpp:45
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.
Computes the two-dimensional free-current density as a vector curl of the current vector potential,...
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.