BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
kink_profile_guesser.hpp
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 #ifndef BUBBLEPROFILER_KINK_PROFILE_GUESSER_HPP_INCLUDED
19 #define BUBBLEPROFILER_KINK_PROFILE_GUESSER_HPP_INCLUDED
20 
21 #include "basic_logger.hpp"
22 #include "profile_guesser.hpp"
23 
24 #include <Eigen/Core>
25 
26 #include <array>
27 
28 namespace BubbleProfiler {
29 
35 public:
36  virtual ~Kink_profile_guesser() = default;
37 
58  const Eigen::VectorXd&,
59  int, double, double,
60  double, double) override;
61 
73  void set_trial_distance(double dist) { trial_dist = dist; }
74 
86  void set_max_domain_start(double r);
87 
99  void set_min_domain_end(double r);
100 
101 private:
103  const static std::array<double,46> alpha_grid_3d;
105  const static std::array<double, 50> alpha_grid_thin_3d;
107  const static std::array<double,46> phi0_grid_3d;
109  const static std::array<double,46> delta_grid_3d;
111  const static std::array<double, 50> delta_grid_thin_3d;
113  const static std::array<double,46> lw_grid_3d;
114 
116  const static std::array<double,44> alpha_grid_4d;
118  const static std::array<double, 50> alpha_grid_thin_4d;
120  const static std::array<double,44> phi0_grid_4d;
122  const static std::array<double,44> delta_grid_4d;
124  const static std::array<double, 50> delta_grid_thin_4d;
126  const static std::array<double,44> lw_grid_4d;
127 
129  double trial_dist{0.5};
131  double max_domain_start{1.e-8};
133  double min_domain_end{4.};
135  double alpha{0.5};
137  double aE{1.};
139  double phi0{0.};
141  double delta{0.};
143  double lw{0.};
145  double dist_true_vacuum{0.};
147  int num_fields{0};
148 
153  Eigen::MatrixXd cob_matrix{};
155 
156  void compute_vacuum_distance(const Potential&, const Eigen::VectorXd&);
158  const Eigen::VectorXd&);
159  void fit_ansatz_parameters(const Potential&,
160  const Eigen::VectorXd&, int);
161  double evaluate_ansatz_at(double) const;
162  double evaluate_ansatz_deriv_at(double) const;
163  double guess_domain_start() const;
164  double guess_domain_end() const;
165 
167  double alpha_threshold = 0.514;
168 
169  Field_profiles calculate_field_profiles(int, double, double, double, double);
170 };
171 
172 } // namespace BubbleProfiler
173 
174 #endif
void set_min_domain_end(double r)
Sets the minimum required distance if guessing the domain end.
void set_max_domain_start(double r)
Sets the maximum allowed distance if guessing the domain start.
static const std::array< double, 44 > delta_grid_4d
static const std::array< double, 46 > delta_grid_3d
static const std::array< double, 46 > phi0_grid_3d
static const std::array< double, 44 > lw_grid_4d
static const std::array< double, 46 > alpha_grid_3d
void compute_vacuum_distance(const Potential &, const Eigen::VectorXd &)
void calculate_potential_parameters(const Potential &, const Eigen::VectorXd &)
double alpha_threshold
Minimum value of alpha before throwing a Thin_wall_error.
static const std::array< double, 50 > alpha_grid_thin_3d
void fit_ansatz_parameters(const Potential &, const Eigen::VectorXd &, int)
Abstract class to represent ansatz generators.
static const std::array< double, 46 > lw_grid_3d
void set_trial_distance(double dist)
Sets the trial distance used in estimating the potential parameters.
static const std::array< double, 50 > delta_grid_thin_4d
virtual Field_profiles get_profile_guess(const Potential &, const Eigen::VectorXd &, int, double, double, double, double) override
Calculate an initial guess for the bubble profile.
Field_profiles calculate_field_profiles(int, double, double, double, double)
static const std::array< double, 50 > delta_grid_thin_3d
static const std::array< double, 44 > alpha_grid_4d
Computes an initial guess for the bubble profile using a kink ansatz.
static const std::array< double, 50 > alpha_grid_thin_4d
static const std::array< double, 44 > phi0_grid_4d
Abstract base class for a generic potential.
Definition: potential.hpp:36
Discretized set of field profiles.