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_T_to_J.hpp"
22 #include "solver.hpp"
23 
24 #include "settings.hpp"
25 
26 using namespace Misc;
27 using namespace StaticVectorSolver;
28 
34 class BatchCVPI : public SettingsCVPI
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/sphere-linear/";
44 
45  std::cout << "Program: cvp-i\n"
46  << "Dimensions: "
47  << "3\n"
48  << "Writing to: " << dir << "\n";
49 
50  MainOutputTable table_J(3);
51 
52  for (unsigned int p = 0; p < 3; p++) {
53  table_J.clear();
54 
55  for (unsigned int r = 12; r < 16; r++) {
56  std::cout << "Solving for T. \n";
57 
58  fname =
59  dir + "solution_T_p" + std::to_string(p) + "_r" + std::to_string(r);
60 
62  std::cout << "Time table\n";
63 
64  SolverCVPI problem(p, 2, r, fname);
65  problem.clear();
66 
67  std::cout << "Converting T to J. \n";
68 
69  fname =
70  dir + "solution_J_p" + std::to_string(p) + "_r" + std::to_string(r);
71 
72  table_J.add_value("r", r);
73  table_J.add_value("p", p);
74 
75  ExactSolutionCVPI_Jf exact_solution;
76 
78  std::cout << "Time table\n";
79 
80  ProjectTtoJ projector(p,
81  2,
82  problem.get_tria(),
83  problem.get_dof_handler(),
84  problem.get_solution(),
85  fname,
86  &exact_solution,
90  true);
91 
92  table_J.add_value("ndofs", projector.get_n_dofs());
93  table_J.add_value("ncells", projector.get_n_cells());
94  table_J.add_value("L2", projector.get_L2_norm());
95  table_J.add_value("H1", 0.0);
96  }
97  // Saving convergence table
98 
99  std::cout << "Table J\n";
100  table_J.save(dir + "table_J_p" + std::to_string(p));
101  }
102  }
103 };
104 
105 int
106 main()
107 {
108  try {
109  BatchCVPI batch;
110  batch.run();
111  } catch (std::exception& exc) {
112  std::cerr << std::endl
113  << std::endl
114  << "----------------------------------------------------"
115  << std::endl;
116  std::cerr << "Exception on processing: " << std::endl
117  << exc.what() << std::endl
118  << "Aborting!" << std::endl
119  << "----------------------------------------------------"
120  << std::endl;
121  return 1;
122  } catch (...) {
123  std::cerr << std::endl
124  << std::endl
125  << "----------------------------------------------------"
126  << std::endl;
127  std::cerr << "Unknown exception!" << std::endl
128  << "Aborting!" << std::endl
129  << "----------------------------------------------------"
130  << std::endl;
131  return 1;
132  }
133  return 0;
134 }
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-i/) 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 Current vector potential (cvp-i/) 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-i/) numerical experiment.
Definition: solver.hpp:45
unsigned int get_n_cells() const
Returns the number of active cells in the mesh.
unsigned int get_n_dofs() const
Returns the total amount of the degrees of freedom.
Computes the volume free-current density as a curl of the current vector potential,...
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.