BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
generalized_fubini_observer.cpp
Go to the documentation of this file.
1 /*
2  * This file is part of BubbleProfiler.
3  *
4  * BubbleProfiler is free software: you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation, either version 3 of the License, or
7  * (at your option) any later version.
8  *
9  * BubbleProfiler is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with BubbleProfiler. If not, see <http://www.gnu.org/licenses/>.
16  */
17 
19 #include "error.hpp"
20 #include "field_profiles.hpp"
21 
22 #include <iomanip>
23 
24 namespace BubbleProfiler {
25 
27  const Generalized_fubini_potential& potential_,
28  const std::string& output_file_)
29  : potential(potential_)
30  , output_file(output_file_)
31 {
33 }
34 
36 {
37  if (!output_stream.is_open() || !output_stream.good()) {
38  throw IO_error("unable to write to file '" + output_file + "'");
39  }
40  output_stream << "# "
41  << std::setw(12) << "iteration" << ' '
42  << std::setw(16) << "rho" << ' '
43  << std::setw(16) << "perturbation" << ' '
44  << std::setw(16) << "profile" << ' '
45  << std::setw(16) << "exact"
46  << std::endl;
47 }
48 
50  const Field_profiles& perturbation,
51  const Field_profiles& exact)
52 {
53  if (!output_stream.is_open() || !output_stream.good()) {
54  throw IO_error("unable to write to file '" + output_file + "'");
55  }
56 
57  const auto coord_values = profile.get_spatial_grid();
58  const auto numerical_values = profile.get_field_profiles();
59  const auto perturbation_values = perturbation.get_field_profiles();
60  const auto exact_values = exact.get_field_profiles();
61 
62  const int n_grid_points = coord_values.size();
63  for (int i = 0; i < n_grid_points; ++i) {
64  output_stream << " "
65  << std::setw(12) << iteration_count << ' '
66  << std::setw(16) << std::setprecision(8)
67  << std::scientific << coord_values(i) << ' '
68  << std::setw(16) << std::setprecision(8)
69  << std::scientific << perturbation_values(i, 0) << ' '
70  << std::setw(16) << std::setprecision(8)
71  << std::scientific << numerical_values(i, 0) << ' '
72  << std::setw(16) << std::setprecision(8)
73  << std::scientific << exact_values(i, 0) << std::endl;
74  }
75 }
76 
78  const Field_profiles& perturbation)
79 {
80  if (iteration_count == 0) {
82  }
83 
84  const auto exact = potential.get_profile(profile.get_spatial_grid());
85  write_data(profile, perturbation, exact);
87 }
88 
89 } // namespace BubbleProfiler
contains the definition of the Field_profiles clas
const Eigen::VectorXd & get_spatial_grid() const
Get a vector of the grid point coordinates.
Generalized_fubini_observer(const Generalized_fubini_potential &potential_, const std::string &output_file_)
void write_data(const Field_profiles &, const Field_profiles &, const Field_profiles &)
void operator()(const Field_profiles &profile, const Field_profiles &perturbation)
Exception indicating generic input/output error.
Definition: error.hpp:122
Discretized set of field profiles.
const Eigen::MatrixXd & get_field_profiles() const
Get the field profile data in matrix form.
Field_profiles get_profile(const Eigen::VectorXd &) const