Green’s Functions

We will here describe the inheritance hierarchy for generating Green’s functions, in order to use and extend it properly. The runtime creation of Green’s functions objects relies on the Factory Method pattern [GHJV94][Ale01], implemented through the generic Factory class.

../_images/green.png

IGreensFunction

class IGreensFunction

Interface for Green’s function classes.

The Non-Virtual Interface (NVI) idiom is used. Notice also that some of the return types are in some cases “auto”, meaning that we will let the compiler deduce them.

Author
Luca Frediani and Roberto Di Remigio
Date
2012-2015

Subclassed by GreensFunction< DerivativeTraits, IntegratorPolicy, ProfilePolicy, Derived >, GreensFunction< DerivativeTraits, IntegratorPolicy, Anisotropic, AnisotropicLiquid< DerivativeTraits, IntegratorPolicy > >, GreensFunction< DerivativeTraits, IntegratorPolicy, Sharp, SphericalSharp< DerivativeTraits, IntegratorPolicy > >, GreensFunction< DerivativeTraits, IntegratorPolicy, Uniform, UniformDielectric< DerivativeTraits, IntegratorPolicy > >, GreensFunction< DerivativeTraits, IntegratorPolicy, Uniform, Vacuum< DerivativeTraits, IntegratorPolicy > >, GreensFunction< DerivativeTraits, IntegratorPolicy, Yukawa, IonicLiquid< DerivativeTraits, IntegratorPolicy > >, GreensFunction< Numerical, IntegratorPolicy, ProfilePolicy, Derived >, GreensFunction< Numerical, IntegratorPolicy, ProfilePolicy, SphericalDiffuse< IntegratorPolicy, ProfilePolicy > >

Public Functions

double kernelS(const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the kernel of the \(\mathcal{S}\) integral operator, i.e. the value of the Greens’s function for the pair of points p1, p2: \( G(\mathbf{p}_1, \mathbf{p}_2)\)

Parameters
  • p1 -

    first point

  • p2 -

    second point

double kernelD(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the kernel of the \(\mathcal{D}\) integral operator for the pair of points p1, p2: \( [\boldsymbol{\varepsilon}\nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)]\cdot \mathbf{n}_{\mathbf{p}_2}\) To obtain the kernel of the \(\mathcal{D}^\dagger\) operator call this methods with \(\mathbf{p}_1\) and \(\mathbf{p}_2\) exchanged and with \(\mathbf{n}_{\mathbf{p}_2} = \mathbf{n}_{\mathbf{p}_1}\)

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

virtual bool uniform() const = 0

Whether the Green’s function describes a uniform environment

virtual Permittivity permittivity() const = 0

Returns a dielectric permittivity profile

virtual Eigen::MatrixXd singleLayer(const std::vector<Element> &e) const = 0

Calculates the matrix representation of the S operator

Parameters
  • e -

    list of finite elements

virtual Eigen::MatrixXd doubleLayer(const std::vector<Element> &e) const = 0

Calculates the matrix representation of the D operator

Parameters
  • e -

    list of finite elements

Protected Functions

virtual double kernelS_impl(const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const = 0

Returns value of the kernel of the \(\mathcal{S}\) integral operator, i.e. the value of the Greens’s function for the pair of points p1, p2: \( G(\mathbf{p}_1, \mathbf{p}_2)\)

Parameters
  • p1 -

    first point

  • p2 -

    second point

virtual double kernelD_impl(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const = 0

Returns value of the kernel of the \(\mathcal{D}\) integral operator for the pair of points p1, p2: \( [\boldsymbol{\varepsilon}\nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)]\cdot \mathbf{n}_{\mathbf{p}_2}\) To obtain the kernel of the \(\mathcal{D}^\dagger\) operator call this methods with \(\mathbf{p}_1\) and \(\mathbf{p}_2\) exchanged and with \(\mathbf{n}_{\mathbf{p}_2} = \mathbf{n}_{\mathbf{p}_1}\)

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

GreensFunction

template <typename DerivativeTraits, typename IntegratorPolicy, typename ProfilePolicy, typename Derived>
class GreensFunction

Templated interface for Green’s functions.

Author
Luca Frediani and Roberto Di Remigio
Date
2012-2014
Template Parameters
  • DerivativeTraits -

    evaluation strategy for the function and its derivatives

  • IntegratorPolicy -

    policy for the calculation of the matrix represenation of S and D

  • ProfilePolicy -

    dielectric profile type

Inherits from IGreensFunction

Public Functions

virtual double derivativeSource(const Eigen::Vector3d &normal_p1, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the directional derivative of the Greens’s function for the pair of points p1, p2: \( \nabla_{\mathbf{p_1}}G(\mathbf{p}_1, \mathbf{p}_2)\cdot \mathbf{n}_{\mathbf{p}_1}\) Notice that this method returns the directional derivative with respect to the source point.

Parameters
  • normal_p1 -

    the normal vector to p1

  • p1 -

    first point

  • p2 -

    second point

virtual double derivativeProbe(const Eigen::Vector3d &normal_p2, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the directional derivative of the Greens’s function for the pair of points p1, p2: \( \nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)\cdot \mathbf{n}_{\mathbf{p}_2}\) Notice that this method returns the directional derivative with respect to the probe point.

Parameters
  • normal_p2 -

    the normal vector to p2

  • p1 -

    first point

  • p2 -

    second point

Eigen::Vector3d gradientSource(const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns full gradient of Greens’s function for the pair of points p1, p2: \( \nabla_{\mathbf{p_1}}G(\mathbf{p}_1, \mathbf{p}_2)\) Notice that this method returns the gradient with respect to the source point.

Parameters
  • p1 -

    first point

  • p2 -

    second point

Eigen::Vector3d gradientProbe(const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns full gradient of Greens’s function for the pair of points p1, p2: \( \nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)\) Notice that this method returns the gradient with respect to the probe point.

Parameters
  • p1 -

    first point

  • p2 -

    second point

virtual bool uniform() const

Whether the Green’s function describes a uniform environment

virtual Permittivity permittivity() const

Returns a dielectric permittivity profile

Protected Functions

virtual DerivativeTraits operator()(DerivativeTraits *source, DerivativeTraits *probe) const = 0

Evaluates the Green’s function given a pair of points

Parameters
  • source -

    the source point

  • probe -

    the probe point

virtual double kernelS_impl(const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the kernel of the \(\mathcal{S}\) integral operator, i.e. the value of the Greens’s function for the pair of points p1, p2: \( G(\mathbf{p}_1, \mathbf{p}_2)\)

Note
Relies on the implementation of operator() in the subclasses and that is all subclasses need to implement. Thus this method is marked __final.
Parameters
  • p1 -

    first point

  • p2 -

    second point

Vacuum

template <typename DerivativeTraits = AD_directional, typename IntegratorPolicy = CollocationIntegrator>
class Vacuum

Green’s function for vacuum.

Author
Luca Frediani and Roberto Di Remigio
Date
2012-2014
Template Parameters
  • DerivativeTraits -

    evaluation strategy for the function and its derivatives

  • IntegratorPolicy -

    policy for the calculation of the matrix represenation of S and D

Inherits from GreensFunction< DerivativeTraits, IntegratorPolicy, Uniform, Vacuum< DerivativeTraits, IntegratorPolicy > >

Public Functions

virtual Eigen::MatrixXd singleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the S operator

Parameters
  • e -

    list of finite elements

virtual Eigen::MatrixXd doubleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the D operator

Parameters
  • e -

    list of finite elements

Private Functions

virtual DerivativeTraits operator()(DerivativeTraits *sp, DerivativeTraits *pp) const

Evaluates the Green’s function given a pair of points

Parameters
  • sp -

    the source point

  • pp -

    the probe point

virtual double kernelD_impl(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the kernel of the \(\mathcal{D}\) integral operator for the pair of points p1, p2: \( [\boldsymbol{\varepsilon}\nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)]\cdot \mathbf{n}_{\mathbf{p}_2}\) To obtain the kernel of the \(\mathcal{D}^\dagger\) operator call this methods with \(\mathbf{p}_1\) and \(\mathbf{p}_2\) exchanged and with \(\mathbf{n}_{\mathbf{p}_2} = \mathbf{n}_{\mathbf{p}_1}\)

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

UniformDielectric

template <typename DerivativeTraits = AD_directional, typename IntegratorPolicy = CollocationIntegrator>
class UniformDielectric

Green’s function for uniform dielectric.

Author
Luca Frediani and Roberto Di Remigio
Date
2012-2015
Template Parameters
  • DerivativeTraits -

    evaluation strategy for the function and its derivatives

  • IntegratorPolicy -

    policy for the calculation of the matrix represenation of S and D

Inherits from GreensFunction< DerivativeTraits, IntegratorPolicy, Uniform, UniformDielectric< DerivativeTraits, IntegratorPolicy > >

Public Functions

virtual Eigen::MatrixXd singleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the S operator

Parameters
  • e -

    list of finite elements

virtual Eigen::MatrixXd doubleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the D operator

Parameters
  • e -

    list of finite elements

Private Functions

virtual DerivativeTraits operator()(DerivativeTraits *sp, DerivativeTraits *pp) const

Evaluates the Green’s function given a pair of points

Parameters
  • sp -

    the source point

  • pp -

    the probe point

virtual double kernelD_impl(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the kernel of the \(\mathcal{D}\) integral operator for the pair of points p1, p2: \( [\boldsymbol{\varepsilon}\nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)]\cdot \mathbf{n}_{\mathbf{p}_2}\) To obtain the kernel of the \(\mathcal{D}^\dagger\) operator call this methods with \(\mathbf{p}_1\) and \(\mathbf{p}_2\) exchanged and with \(\mathbf{n}_{\mathbf{p}_2} = \mathbf{n}_{\mathbf{p}_1}\)

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

IonicLiquid

template <typename DerivativeTraits = AD_directional, typename IntegratorPolicy = CollocationIntegrator>
class IonicLiquid

Green’s functions for ionic liquid, described by the linearized Poisson-Boltzmann equation.

Author
Luca Frediani, Roberto Di Remigio
Date
2013-2015
Template Parameters
  • DerivativeTraits -

    evaluation strategy for the function and its derivatives

  • IntegratorPolicy -

    policy for the calculation of the matrix represenation of S and D

Inherits from GreensFunction< DerivativeTraits, IntegratorPolicy, Yukawa, IonicLiquid< DerivativeTraits, IntegratorPolicy > >

Public Functions

virtual Eigen::MatrixXd singleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the S operator

Parameters
  • e -

    list of finite elements

virtual Eigen::MatrixXd doubleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the D operator

Parameters
  • e -

    list of finite elements

Private Functions

virtual DerivativeTraits operator()(DerivativeTraits *sp, DerivativeTraits *pp) const

Evaluates the Green’s function given a pair of points

Parameters
  • sp -

    the source point

  • pp -

    the probe point

virtual double kernelD_impl(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the directional derivative of the Greens’s function for the pair of points p1, p2: \( \nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)\cdot \mathbf{n}_{\mathbf{p}_2}\) Notice that this method returns the directional derivative with respect to the probe point, thus assuming that the direction is relative to that point.

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

AnisotropicLiquid

template <typename DerivativeTraits = AD_directional, typename IntegratorPolicy = CollocationIntegrator>
class AnisotropicLiquid

Green’s functions for anisotropic liquid, described by a tensorial permittivity.

Author
Roberto Di Remigio
Date
2015
Template Parameters
  • DerivativeTraits -

    evaluation strategy for the function and its derivatives

  • IntegratorPolicy -

    policy for the calculation of the matrix represenation of S and D

Inherits from GreensFunction< DerivativeTraits, IntegratorPolicy, Anisotropic, AnisotropicLiquid< DerivativeTraits, IntegratorPolicy > >

Public Functions

AnisotropicLiquid(const Eigen::Vector3d &eigen_eps, const Eigen::Vector3d &euler_ang)

Parameters
  • eigen_eps -

    eigenvalues of the permittivity tensors

  • euler_ang -

    Euler angles in degrees

AnisotropicLiquid(const Eigen::Vector3d &eigen_eps, const Eigen::Vector3d &euler_ang, double f)

Parameters
  • eigen_eps -

    eigenvalues of the permittivity tensors

  • euler_ang -

    Euler angles in degrees

virtual Eigen::MatrixXd singleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the S operator

Parameters
  • e -

    list of finite elements

virtual Eigen::MatrixXd doubleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the D operator

Parameters
  • e -

    list of finite elements

Private Functions

virtual DerivativeTraits operator()(DerivativeTraits *source, DerivativeTraits *probe) const

Evaluates the Green’s function given a pair of points

Parameters
  • source -

    the source point

  • probe -

    the probe point

virtual double kernelD_impl(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the kernel for the calculation of the \(\mathcal{D}\) integral operator for the pair of points p1, p2: \( [\boldsymbol{\varepsilon}\nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)]\cdot \mathbf{n}_{\mathbf{p}_2}\) To obtain the kernel of the \(\mathcal{D}^\dagger\) operator call this methods with \(\mathbf{p}_1\) and \(\mathbf{p}_2\) exchanged and with \(\mathbf{n}_{\mathbf{p}_2} = \mathbf{n}_{\mathbf{p}_1}\)

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

SphericalDiffuse

template <typename IntegratorPolicy = CollocationIntegrator, typename ProfilePolicy = OneLayerTanh>
class SphericalDiffuse

Green’s function for a diffuse interface with spherical symmetry.

This class is general, in the sense that no specific dielectric profile has been set in its definition. In principle any profile that can be described by:

  1. a left-side dielectric constant;
  2. a right-side dielectric constant;
  3. an interface layer width;
  4. an interface layer center can be used to define a new diffuse interface with spherical symmetry. The origin of the dielectric sphere can be changed by means of the constructor. The solution of the differential equation defining the Green’s function is always performed assuming that the dielectric sphere is centered in the origin of the coordinate system. Whenever the public methods are invoked to “sample” the Green’s function at a pair of points, a translation of the sampling points is performed first.
Author
Hui Cao, Ville Weijo, Luca Frediani and Roberto Di Remigio
Date
2010-2015
Template Parameters
  • IntegratorPolicy -

    policy for the calculation of the matrix represenation of S and D

  • ProfilePolicy -

    functional form of the diffuse layer

Inherits from GreensFunction< Numerical, IntegratorPolicy, ProfilePolicy, SphericalDiffuse< IntegratorPolicy, ProfilePolicy > >

Unnamed Group

int maxLGreen_

Parameters and functions for the calculation of the Green’s function, including Coulomb singularity

Maximum angular momentum in the __final summation over Legendre polynomials to obtain G

std::vector<RadialFunction<interfaces::StateType, interfaces::LnTransformedRadial, Zeta>> zeta_

First independent radial solution, used to build Green’s function.

Note
The vector has dimension maxLGreen_ and has r^l behavior

std::vector<RadialFunction<interfaces::StateType, interfaces::LnTransformedRadial, Omega>> omega_

Second independent radial solution, used to build Green’s function.

Note
The vector has dimension maxLGreen_ and has r^(-l-1) behavior

double imagePotentialComponent_impl(int L, const Eigen::Vector3d &sp, const Eigen::Vector3d &pp, double Cr12) const

Returns L-th component of the radial part of the Green’s function.

Note
This function shifts the given source and probe points by the location of the dielectric sphere.
Parameters
  • L -

    angular momentum

  • sp -

    source point

  • pp -

    probe point

  • Cr12 -

    Coulomb singularity separation coefficient

Unnamed Group

int maxLC_

Parameters and functions for the calculation of the Coulomb singularity separation coefficient

Maximum angular momentum to obtain C(r, r’), needed to separate the Coulomb singularity

RadialFunction<interfaces::StateType, interfaces::LnTransformedRadial, Zeta> zetaC_

First independent radial solution, used to build coefficient.

Note
This is needed to separate the Coulomb singularity and has r^l behavior

RadialFunction<interfaces::StateType, interfaces::LnTransformedRadial, Omega> omegaC_

Second independent radial solution, used to build coefficient.

Note
This is needed to separate the Coulomb singularity and has r^(-l-1) behavior

double coefficient_impl(const Eigen::Vector3d &sp, const Eigen::Vector3d &pp) const

Returns coefficient for the separation of the Coulomb singularity.

Note
This function shifts the given source and probe points by the location of the dielectric sphere.
Parameters
  • sp -

    first point

  • pp -

    second point

Public Functions

SphericalDiffuse(double e1, double e2, double w, double c, const Eigen::Vector3d &o, int l)

Constructor for a one-layer interface

Parameters
  • e1 -

    left-side dielectric constant

  • e2 -

    right-side dielectric constant

  • w -

    width of the interface layer

  • c -

    center of the diffuse layer

  • o -

    center of the sphere

SphericalDiffuse(double e1, double e2, double w, double c, const Eigen::Vector3d &o, int l, double f)

Constructor for a one-layer interface

Parameters
  • e1 -

    left-side dielectric constant

  • e2 -

    right-side dielectric constant

  • w -

    width of the interface layer

  • c -

    center of the diffuse layer

  • o -

    center of the sphere

virtual Eigen::MatrixXd singleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the S operator

Parameters
  • e -

    list of finite elements

virtual Eigen::MatrixXd doubleLayer(const std::vector<Element> &e) const

Calculates the matrix representation of the D operator

Parameters
  • e -

    list of finite elements

double coefficientCoulomb(const Eigen::Vector3d &source, const Eigen::Vector3d &probe) const

Returns Coulomb singularity separation coefficient.

Parameters
  • source -

    location of the source charge

  • probe -

    location of the probe charge

double Coulomb(const Eigen::Vector3d &source, const Eigen::Vector3d &probe) const

Returns singular part of the Green’s function.

Parameters
  • source -

    location of the source charge

  • probe -

    location of the probe charge

double imagePotential(const Eigen::Vector3d &source, const Eigen::Vector3d &probe) const

Returns non-singular part of the Green’s function (image potential)

Parameters
  • source -

    location of the source charge

  • probe -

    location of the probe charge

double coefficientCoulombDerivative(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the directional derivative of the Coulomb singularity separation coefficient for the pair of points p1, p2: \( \nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)\cdot \mathbf{n}_{\mathbf{p}_2}\) Notice that this method returns the directional derivative with respect to the probe point, thus assuming that the direction is relative to that point.

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

double CoulombDerivative(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the directional derivative of the singular part of the Greens’s function for the pair of points p1, p2: \( \nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)\cdot \mathbf{n}_{\mathbf{p}_2}\) Notice that this method returns the directional derivative with respect to the probe point, thus assuming that the direction is relative to that point.

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

double imagePotentialDerivative(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the directional derivative of the non-singular part (image potential) of the Greens’s function for the pair of points p1, p2: \( \nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)\cdot \mathbf{n}_{\mathbf{p}_2}\) Notice that this method returns the directional derivative with respect to the probe point, thus assuming that the direction is relative to that point.

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

pcm::tuple<double, double> epsilon(const Eigen::Vector3d &point) const

Handle to the dielectric profile evaluation

Private Functions

virtual Numerical operator()(Numerical *sp, Numerical *pp) const

Evaluates the Green’s function given a pair of points

Note
This takes care of the origin shift
Parameters
  • sp -

    the source point

  • pp -

    the probe point

virtual double kernelD_impl(const Eigen::Vector3d &direction, const Eigen::Vector3d &p1, const Eigen::Vector3d &p2) const

Returns value of the kernel of the \(\mathcal{D}\) integral operator for the pair of points p1, p2: \( [\boldsymbol{\varepsilon}\nabla_{\mathbf{p_2}}G(\mathbf{p}_1, \mathbf{p}_2)]\cdot \mathbf{n}_{\mathbf{p}_2}\) To obtain the kernel of the \(\mathcal{D}^\dagger\) operator call this methods with \(\mathbf{p}_1\) and \(\mathbf{p}_2\) exchanged and with \(\mathbf{n}_{\mathbf{p}_2} = \mathbf{n}_{\mathbf{p}_1}\)

Parameters
  • direction -

    the direction

  • p1 -

    first point

  • p2 -

    second point

void initProfilePolicy(double e1, double e2, double w, double c)

Initializes a one-layer profile

Parameters
  • e1 -

    left-side dielectric constant

  • e2 -

    right-side dielectric constant

  • w -

    width of the interface layer

  • c -

    center of the diffuse layer

void initSphericalDiffuse()

This calculates all the components needed to evaluate the Green’s function

Private Members

Eigen::Vector3d origin_

Center of the dielectric sphere