BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
observers.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 
18 #include "observers.hpp"
19 #include "euclidean_action.hpp"
20 #include "field_profiles.hpp"
21 
22 #include <Eigen/Core>
23 
24 #include <iomanip>
25 #include <iostream>
26 
27 namespace BubbleProfiler {
28 
29 Profile_observer::Profile_observer(const std::vector<std::string>& fields_,
30  const std::string& output_path_,
31  const Potential& potential_)
32  : fields(fields_)
33  , output_path(output_path_)
34  , potential(potential_)
35 {
36  action_file.open(output_path + "/action.txt");
37  profiles_file.open(output_path + "/field_profiles.txt");
38  perturbations_file.open(output_path + "/perturbations.txt");
39 }
40 
41 
43  const Field_profiles& profiles)
44 {
45  action_file << "#" << ' '
46  << std::setw(12) << "perturbation" << ' '
47  << std::setw(16) << "action" << std::endl;
48 
49  write_action_to_file(profiles);
50 }
51 
53  const Field_profiles& profiles)
54 {
55  logger.log_message(logging::Log_level::Trace, "Writing ansatz to file");
56 
57  profiles_file << "#" << ' '
58  << std::setw(12) << "perturbation" << ' '
59  << std::setw(16) << "rho";
60  for (const std::string& name: fields) {
61  profiles_file << ' ' << std::setw(16) << name;
62  }
63  profiles_file << std::endl;
64 
66 }
67 
69 {
70  perturbations_file << "#" << ' '
71  << std::setw(12) << "perturbation" << ' '
72  << std::setw(16) << "rho";
73  for (const std::string& name: fields) {
74  perturbations_file << ' ' << std::setw(16) << "eps_" + name;
75  }
76  perturbations_file << std::endl;
77 }
78 
80  std::ofstream& stream, const Field_profiles& profiles)
81 {
82  const auto coord_values = profiles.get_spatial_grid();
83  const auto field_values = profiles.get_field_profiles();
84 
85  const int n_grid_points = coord_values.size();
86  const int n_fields = fields.size();
87 
88  for (int i = 0; i < n_grid_points; ++i) {
89  stream << " "
90  << std::setw(12) << iteration_count << ' '
91  << std::setw(16) << std::setprecision(8)
92  << std::scientific << coord_values(i);
93  for (int j = 0; j < n_fields; ++j) {
94  stream << ' '
95  << std::setw(16) << std::setprecision(8)
96  << std::scientific << field_values(i, j);
97  }
98  stream << std::endl;
99  }
100 }
101 
103 {
104  const double action = calculate_action(potential, profiles);
105 
107  + std::to_string(action));
108 
109  action_file << " "
110  << std::setw(12) << iteration_count << ' '
111  << std::setw(16) << std::setprecision(8)
112  << std::scientific << action << std::endl;
113 }
114 
116  const Field_profiles& profiles, const Field_profiles& perturbations)
117 {
118  if (iteration_count == 0) {
122  } else {
125  write_action_to_file(profiles);
126  }
127  iteration_count++;
128 }
129 
130 } // namespace BubbleProfiler
contains the definition of the Field_profiles clas
contains helper functions for calculating the Euclidean action
const Eigen::VectorXd & get_spatial_grid() const
Get a vector of the grid point coordinates.
void write_action_to_file(const Field_profiles &profiles)
Definition: observers.cpp:102
void write_initial_action_to_file(const Field_profiles &profiles)
Definition: observers.cpp:42
void write_initial_profiles_to_file(const Field_profiles &profiles)
Definition: observers.cpp:52
void write_field_profiles_to_file(std::ofstream &stream, const Field_profiles &profiles)
Definition: observers.cpp:79
void log_message(Log_level level, const std::string &msg) const
void operator()(const Field_profiles &profile, const Field_profiles &perturbation)
Definition: observers.cpp:115
Profile_observer(const std::vector< std::string > &fields_, const std::string &output_path_, const Potential &potential_)
Definition: observers.cpp:29
std::vector< std::string > fields
Definition: observers.hpp:62
logging::Basic_logger logger
Definition: observers.hpp:68
Abstract base class for a generic potential.
Definition: potential.hpp:36
double calculate_action(const Potential &potential, const Field_profiles &profiles, std::size_t max_intervals=1000, double rel_tol=1.e-4, double abs_tol=1.e-4, Integration_rule rule=Integration_rule::GK31, bool use_kinetic=true)
Calculates the action using either the kinetic or potential terms.
Discretized set of field profiles.
const Eigen::MatrixXd & get_field_profiles() const
Get the field profile data in matrix form.