BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
algebraic_potential.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_ALGEBRAIC_POTENTIAL_HPP_INCLUDED
19 #define BUBBLEPROFILER_ALGEBRAIC_POTENTIAL_HPP_INCLUDED
20 
21 #include "basic_logger.hpp"
22 #include "potential.hpp"
23 
24 #include <ginac/ginac.h>
25 
26 namespace BubbleProfiler {
27 
29 
37 public:
39 
46  Algebraic_potential(const std::vector<std::string>& fields,
47  const std::string& expr);
48 
49  virtual ~Algebraic_potential() = default;
50 
53  return new Algebraic_potential(*this);
54  }
55 
56  virtual double operator()(const Eigen::VectorXd& coords) const override;
57  virtual double partial(const Eigen::VectorXd& coords, int i) const override;
58  virtual double partial(const Eigen::VectorXd& coords, int i, int j) const override;
59  virtual std::size_t get_number_of_fields() const override {
60  return fields.size(); }
61 
62  virtual void translate_origin(const Eigen::VectorXd&) override;
63  virtual void apply_basis_change(const Eigen::MatrixXd&) override;
64  virtual void add_constant_term(double) override;
65 
66  // Printable
67  inline friend std::ostream& operator<<(
68  std::ostream& os, const Algebraic_potential& p) {
69  os << p.v;
70  return os;
71  }
72 
73 private:
74  std::vector<std::string> fields{};
75  GiNaC::ex v{};
76  std::vector<GiNaC::ex> first_partials{};
77  std::vector<std::vector<GiNaC::ex> > second_partials{};
78  std::vector<GiNaC::symbol> syms{};
80 
81  // Evaluate a GiNaCs expression
82  double eval(const GiNaC::ex&, const Eigen::VectorXd&) const;
83 
84  // Rebuild algebraic derivatives
85  void build_derivatives();
86 };
87 
88 } // namespace BubbleProfiler
89 
90 #endif
virtual void add_constant_term(double) override
Add a constant offset to the potential.
virtual void apply_basis_change(const Eigen::MatrixXd &) override
Apply a change of basis matrix.
virtual std::size_t get_number_of_fields() const override
virtual double partial(const Eigen::VectorXd &coords, int i) const override
Partial derivative WRT coordinate i at a point.
virtual ~Algebraic_potential()=default
Algebraic_potential * clone() const
Get a fresh copy of this potential.
GiNaC::ex v
GiNaC expression to hold potential.
std::vector< GiNaC::symbol > syms
GiNaC symbols matched to fields vec.
Algebraic_potential(const std::vector< std::string > &fields, const std::string &expr)
Construct from list of fields and GiNaC potential string.
std::vector< std::vector< GiNaC::ex > > second_partials
Second partials.
virtual double operator()(const Eigen::VectorXd &coords) const override
Evaluate potential at point.
double eval(const GiNaC::ex &, const Eigen::VectorXd &) const
virtual void translate_origin(const Eigen::VectorXd &) override
Shift the location of the origin by a specified vector.
std::vector< GiNaC::ex > first_partials
First partials.
friend std::ostream & operator<<(std::ostream &os, const Algebraic_potential &p)
Implements a multi dimensional scalar field potential, which may be expressed as almost any algerbrai...
Abstract base class for a generic potential.
Definition: potential.hpp:36