|
| 1 | +/* |
| 2 | + * PCMSolver, an API for the Polarizable Continuum Model |
| 3 | + * Copyright (C) 2019 Roberto Di Remigio, Luca Frediani and contributors. |
| 4 | + * |
| 5 | + * This file is part of PCMSolver. |
| 6 | + * |
| 7 | + * PCMSolver is free software: you can redistribute it and/or modify |
| 8 | + * it under the terms of the GNU Lesser General Public License as published by |
| 9 | + * the Free Software Foundation, either version 3 of the License, or |
| 10 | + * (at your option) any later version. |
| 11 | + * |
| 12 | + * PCMSolver is distributed in the hope that it will be useful, |
| 13 | + * but WITHOUT ANY WARRANTY; without even the implied warranty of |
| 14 | + * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
| 15 | + * GNU Lesser General Public License for more details. |
| 16 | + * |
| 17 | + * You should have received a copy of the GNU Lesser General Public License |
| 18 | + * along with PCMSolver. If not, see <http://www.gnu.org/licenses/>. |
| 19 | + * |
| 20 | + * For information on the complete list of contributors to the |
| 21 | + * PCMSolver API, see: <http://pcmsolver.readthedocs.io/> |
| 22 | + */ |
| 23 | + |
| 24 | +#pragma once |
| 25 | + |
| 26 | +#include <array> |
| 27 | +#include <cmath> |
| 28 | +#include <fstream> |
| 29 | +#include <functional> |
| 30 | +#include <iomanip> |
| 31 | +#include <tuple> |
| 32 | +#include <vector> |
| 33 | + |
| 34 | +#include <Eigen/Core> |
| 35 | + |
| 36 | +#include "utils/MathUtils.hpp" |
| 37 | +#include "utils/RungeKutta.hpp" |
| 38 | +#include "utils/SplineFunction.hpp" |
| 39 | + |
| 40 | +/*! \file InterfacesImpl.hpp */ |
| 41 | + |
| 42 | +namespace pcm { |
| 43 | +namespace green { |
| 44 | +namespace detail { |
| 45 | +/*! \brief Abstract class for a system of ordinary differential equations |
| 46 | + * \tparam Order The order of the ordinary differential equation |
| 47 | + * \author Roberto Di Remigio |
| 48 | + * \date 2018 |
| 49 | + */ |
| 50 | +template <size_t Order = 1> class ODESystem { |
| 51 | +public: |
| 52 | + typedef std::array<double, Order> StateType; |
| 53 | + size_t ODEorder() const { return Order; } |
| 54 | + void operator()(const StateType & f, StateType & dfdx, const double t) const { |
| 55 | + RHS(f, dfdx, t); |
| 56 | + } |
| 57 | + virtual ~ODESystem() {} |
| 58 | + |
| 59 | +private: |
| 60 | + virtual void RHS(const StateType & f, StateType & dfdx, const double t) const = 0; |
| 61 | +}; |
| 62 | + |
| 63 | +/*! \typedef ProfileEvaluator |
| 64 | + * \brief sort of a function pointer to the dielectric profile evaluation function |
| 65 | + */ |
| 66 | +typedef std::function<std::tuple<double, double>(const double)> ProfileEvaluator; |
| 67 | + |
| 68 | +/*! \class LnTransformedRadial |
| 69 | + * \brief system of ln-transformed first-order radial differential equations |
| 70 | + * \author Roberto Di Remigio |
| 71 | + * \date 2015 |
| 72 | + * |
| 73 | + * Provides a handle to the system of differential equations for the integrator. |
| 74 | + * The dielectric profile comes in as a boost::function object. |
| 75 | + */ |
| 76 | +class LnTransformedRadial final : public pcm::green::detail::ODESystem<2> { |
| 77 | +public: |
| 78 | + /*! Type of the state vector of the ODE */ |
| 79 | + typedef pcm::green::detail::ODESystem<2>::StateType StateType; |
| 80 | + /*! Constructor from profile evaluator and angular momentum */ |
| 81 | + LnTransformedRadial(const ProfileEvaluator & e, int lval) : eval_(e), l_(lval) {} |
| 82 | + |
| 83 | +private: |
| 84 | + /*! Dielectric profile function and derivative evaluation */ |
| 85 | + const ProfileEvaluator eval_; |
| 86 | + /*! Angular momentum */ |
| 87 | + const int l_; |
| 88 | + /*! \brief Provides a functor for the evaluation of the system of first-order ODEs. |
| 89 | + * \param[in] rho state vector holding the function and its first derivative |
| 90 | + * \param[out] drhodr state vector holding the first and second derivative |
| 91 | + * \param[in] y logarithmic position on the integration grid |
| 92 | + * |
| 93 | + * The second-order ODE and the system of first-order ODEs |
| 94 | + * are reported in the manuscript. |
| 95 | + */ |
| 96 | + virtual void RHS(const StateType & rho, StateType & drhodr, const double y) const { |
| 97 | + // Evaluate the dielectric profile |
| 98 | + double eps = 0.0, epsPrime = 0.0; |
| 99 | + std::tie(eps, epsPrime) = eval_(std::exp(y)); |
| 100 | + if (utils::numericalZero(eps)) |
| 101 | + PCMSOLVER_ERROR("Division by zero!"); |
| 102 | + double gamma_epsilon = std::exp(y) * epsPrime / eps; |
| 103 | + // System of equations is defined here |
| 104 | + drhodr[0] = rho[1]; |
| 105 | + drhodr[1] = -rho[1] * (rho[1] + 1.0 + gamma_epsilon) + l_ * (l_ + 1); |
| 106 | + } |
| 107 | +}; |
| 108 | +} // namespace detail |
| 109 | + |
| 110 | +using detail::ProfileEvaluator; |
| 111 | + |
| 112 | +/*! \class RadialFunction |
| 113 | + * \brief represents solutions to the radial 2nd order ODE |
| 114 | + * \author Roberto Di Remigio |
| 115 | + * \date 2018 |
| 116 | + * \tparam ODE system of 1st order ODEs replacing the 2nd order ODE |
| 117 | + * |
| 118 | + * The sign of the integration interval, i.e. the sign of the difference |
| 119 | + * between the minimum and the maximum of the interval, is used to discriminate |
| 120 | + * between the two independent solutions (zeta-type and omega-type) of the |
| 121 | + * differential equation: |
| 122 | + * - When the sign is positive (y_min_ < y_max_), we are integrating **forwards** |
| 123 | + * for the zeta-type solution. |
| 124 | + * - When the sign is negative (y_max_ < y_min_), we are integrating **backwards** |
| 125 | + * for the omega-type solution. |
| 126 | + * |
| 127 | + * This is based on an idea Luca had in Wuhan/Chongqing for the planar |
| 128 | + * diffuse interface. |
| 129 | + * With respect to the previous, much more heavily template-d implementation, |
| 130 | + * it introduces a bit more if-else if-else branching in the function_impl and |
| 131 | + * derivative_impl functions. |
| 132 | + * This is largely offset by the better readability and reduced code duplication. |
| 133 | + */ |
| 134 | +template <typename ODE> class RadialFunction final { |
| 135 | +public: |
| 136 | + RadialFunction() : L_(0), y_min_(0.0), y_max_(0.0), y_sign_(1) {} |
| 137 | + RadialFunction(int l, |
| 138 | + double ymin, |
| 139 | + double ymax, |
| 140 | + double ystep, |
| 141 | + const ProfileEvaluator & eval) |
| 142 | + : L_(l), |
| 143 | + y_min_(ymin), |
| 144 | + y_max_(ymax), |
| 145 | + y_sign_(pcm::utils::sign(y_max_ - y_min_)) { |
| 146 | + compute(ystep, eval); |
| 147 | + } |
| 148 | + std::tuple<double, double> operator()(double point) const { |
| 149 | + return std::make_tuple(function_impl(point), derivative_impl(point)); |
| 150 | + } |
| 151 | + friend std::ostream & operator<<(std::ostream & os, RadialFunction & obj) { |
| 152 | + for (size_t i = 0; i < obj.function_[0].size(); ++i) { |
| 153 | + // clang-format off |
| 154 | + os << std::fixed << std::left << std::setprecision(14) |
| 155 | + << obj.function_[0][i] << " " |
| 156 | + << obj.function_[1][i] << " " |
| 157 | + << obj.function_[2][i] << std::endl; |
| 158 | + // clang-format on |
| 159 | + } |
| 160 | + return os; |
| 161 | + } |
| 162 | + |
| 163 | +private: |
| 164 | + typedef typename ODE::StateType StateType; |
| 165 | + /*! Angular momentum of the solution */ |
| 166 | + int L_; |
| 167 | + /*! Lower bound of the integration interval */ |
| 168 | + double y_min_; |
| 169 | + /*! Upper bound of the integration interval */ |
| 170 | + double y_max_; |
| 171 | + /*! The sign of the integration interval */ |
| 172 | + double y_sign_; |
| 173 | + /*! The actual data: grid, function value and first derivative values */ |
| 174 | + std::array<std::vector<double>, 3> function_; |
| 175 | + /*! Reports progress of differential equation integrator */ |
| 176 | + void push_back(const StateType & x, double y) { |
| 177 | + function_[0].push_back(y); |
| 178 | + function_[1].push_back(x[0]); |
| 179 | + function_[2].push_back(x[1]); |
| 180 | + } |
| 181 | + /*! \brief Calculates radial solution |
| 182 | + * \param[in] step ODE integrator step |
| 183 | + * \param[in] eval dielectric profile evaluator function object |
| 184 | + * \return the number of integration steps |
| 185 | + * |
| 186 | + * This function discriminates between the first (zeta-type), i.e. the one |
| 187 | + * with r^l behavior, and the second (omega-type) radial solution, i.e. the |
| 188 | + * one with r^(-l-1) behavior, based on the sign of the integration interval |
| 189 | + * y_sign_. |
| 190 | + */ |
| 191 | + size_t compute(const double step, const ProfileEvaluator & eval) { |
| 192 | + ODE system(eval, L_); |
| 193 | + // Set initial conditions |
| 194 | + StateType init; |
| 195 | + if (y_sign_ > 0.0) { // zeta-type solution |
| 196 | + init[0] = y_sign_ * L_ * y_min_; |
| 197 | + init[1] = y_sign_ * L_; |
| 198 | + } else { // omega-type solution |
| 199 | + init[0] = y_sign_ * (L_ + 1) * y_min_; |
| 200 | + init[1] = y_sign_ * (L_ + 1); |
| 201 | + } |
| 202 | + pcm::utils::RungeKutta4<StateType> stepper; |
| 203 | + size_t nSteps = |
| 204 | + pcm::utils::integrate_const(stepper, |
| 205 | + system, |
| 206 | + init, |
| 207 | + y_min_, |
| 208 | + y_max_, |
| 209 | + y_sign_ * step, |
| 210 | + std::bind(&RadialFunction<ODE>::push_back, |
| 211 | + this, |
| 212 | + std::placeholders::_1, |
| 213 | + std::placeholders::_2)); |
| 214 | + |
| 215 | + // clang-format off |
| 216 | + // Reverse order of function_ if omega-type solution was computed |
| 217 | + // this ensures that they are in ascending order, as later expected by |
| 218 | + // function_impl and derivative_impl |
| 219 | + if (y_sign_ < 0.0) { |
| 220 | + for (auto & comp : function_) { |
| 221 | + std::reverse(comp.begin(), comp.end()); |
| 222 | + } |
| 223 | + // clang-format on |
| 224 | + } |
| 225 | + return nSteps; |
| 226 | + } |
| 227 | + /*! \brief Returns value of function at given point |
| 228 | + * \param[in] point evaluation point |
| 229 | + * |
| 230 | + * We first check if point is below y_min_, if yes we use |
| 231 | + * the asymptotic form. |
| 232 | + */ |
| 233 | + double function_impl(double point) const { |
| 234 | + double val = 0.0; |
| 235 | + if (point <= y_min_ && y_sign_ > 0.0) { // Asymptotic zeta-type solution |
| 236 | + val = y_sign_ * L_ * point; |
| 237 | + } else if (point >= y_min_ && y_sign_ < 0.0) { // Asymptotic omega-type solution |
| 238 | + val = y_sign_ * (L_ + 1) * point; |
| 239 | + } else { |
| 240 | + val = utils::splineInterpolation(point, function_[0], function_[1]); |
| 241 | + } |
| 242 | + return val; |
| 243 | + } |
| 244 | + /*! \brief Returns value of 1st derivative of function at given point |
| 245 | + * \param[in] point evaluation point |
| 246 | + * |
| 247 | + * Below y_min_, the asymptotic form is used. Otherwise we interpolate. |
| 248 | + */ |
| 249 | + double derivative_impl(double point) const { |
| 250 | + double val = 0.0; |
| 251 | + if (point <= y_min_ && y_sign_ > 0.0) { // Asymptotic zeta-type solution |
| 252 | + val = y_sign_ * L_; |
| 253 | + } else if (point >= y_min_ && y_sign_ < 0.0) { // Asymptotic omega-type solution |
| 254 | + val = y_sign_ * (L_ + 1); |
| 255 | + } else { |
| 256 | + val = utils::splineInterpolation(point, function_[0], function_[2]); |
| 257 | + } |
| 258 | + return val; |
| 259 | + } |
| 260 | +}; |
| 261 | + |
| 262 | +/*! \brief Write contents of a RadialFunction to file |
| 263 | + * \param[in] f solution to the radial 2nd order ODE to print |
| 264 | + * \param[in] fname name of the file |
| 265 | + * \author Roberto Di Remigio |
| 266 | + * \date 2018 |
| 267 | + * \tparam ODE system of 1st order ODEs replacing the 2nd order ODE |
| 268 | + */ |
| 269 | +template <typename ODE> |
| 270 | +void writeToFile(RadialFunction<ODE> & f, const std::string & fname) { |
| 271 | + std::ofstream fout; |
| 272 | + fout.open(fname.c_str()); |
| 273 | + fout << f << std::endl; |
| 274 | + fout.close(); |
| 275 | +} |
| 276 | +} // namespace green |
| 277 | +} // namespace pcm |
0 commit comments