BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
shooting_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_SHOOTING_PROFILE_GUESSER_HPP_INCLUDED
19 #define BUBBLEPROFILER_SHOOTING_PROFILE_GUESSER_HPP_INCLUDED
20 
21 #include "basic_logger.hpp"
22 #include "profile_guesser.hpp"
23 #include "shooting.hpp"
24 
25 #include <Eigen/Core>
26 
27 namespace BubbleProfiler {
28 
40 public:
41  virtual ~Shooting_profile_guesser() = default;
42 
44  const Eigen::VectorXd&,
45  int, double, double,
46  double, double) override;
47 
48  void set_trial_distance(double dist) { trial_dist = dist; }
49  void set_solver(const Shooting& s) { solver = s; }
50 
51  const Shooting& get_solver() const { return solver; }
52  Shooting& get_solver() { return solver; }
53 
54 private:
55  double trial_dist{0.5};
56  double alpha{0.5};
57  double E{1.};
58  double dist_true_vacuum{0.};
59  int num_fields{0};
61  // Change of (field space) basis matrix, takes vectors
62  // in the original basis to those in the ansatz basis.
63  Eigen::MatrixXd cob_matrix{};
65 
66  void compute_vacuum_distance(const Potential&, const Eigen::VectorXd&);
68  const Eigen::VectorXd&);
69  double calculate_false_min_location() const;
70  double calculate_barrier_location() const;
71  double calculate_true_min_location() const;
72  double get_large_distance_solution(double, const Field_profiles&, int) const;
73 
74  void solve_one_dimensional(int);
75  Field_profiles calculate_field_profiles(int, double, double, double, double);
76  Eigen::VectorXd update_spatial_grid(const Eigen::VectorXd&, double, double, int) const;
77 };
78 
79 } // namespace BubbleProfiler
80 
81 #endif
Constructs a solution ansatz using the overshoot/undershoot method.
Field_profiles calculate_field_profiles(int, double, double, double, double)
Abstract class to represent ansatz generators.
Eigen::VectorXd update_spatial_grid(const Eigen::VectorXd &, double, double, int) const
double get_large_distance_solution(double, const Field_profiles &, int) const
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.
Solve the one-dimensional problem using the shooting method.
Definition: shooting.hpp:44
One-dimensional shooting method.
void compute_vacuum_distance(const Potential &, const Eigen::VectorXd &)
Abstract base class for a generic potential.
Definition: potential.hpp:36
void calculate_potential_parameters(const Potential &, const Eigen::VectorXd &)
Discretized set of field profiles.