AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
alma::analytic1D::SPR_calculator_FourierLaplace Class Reference

Class for computing the exact analytical RTA single pulse energy density response of the infinite bulk medium in Fourier-Laplace domain. More...

#include <analytic1d.hpp>

Public Member Functions

 SPR_calculator_FourierLaplace (const alma::Crystal_structure *poscar, const alma::Gamma_grid *grid, const Eigen::ArrayXXd *w, double T)
 Constructor: initialise internal variables.
 
void setDirection (const Eigen::Vector3d unitvector)
 Function setting the thermal transport axis.
 
void setLinSpatialGrid (double ximin, double ximax, int Nxi)
 Construct a linear grid of spatial frequencies from ximin to ximax with Nxi elements.
 
void setLogSpatialGrid (double ximin, double ximax, int Nxi)
 Construct a logarithmic grid of spatial frequencies from ximin to ximax with Nxi elements.
 
void setLinTemporalGrid (double fmin, double fmax, int Nf)
 Construct a linear grid of temporal frequencies from fmin to fmax with Nf elements. More...
 
void setLogTemporalGrid (double fmin, double fmax, int Nf)
 Construct a logarithmic grid of temporal frequencies from fmin to fmax with Nf elements. More...
 
Eigen::VectorXd getSpatialFrequencies ()
 Retrieve the spatial frequency grid.
 
Eigen::VectorXd getTemporalFrequencies ()
 Retrieve the temporal frequency grid.
 
Eigen::MatrixXcd getSPR ()
 Compute SPR and retrieve it.
 

Detailed Description

Class for computing the exact analytical RTA single pulse energy density response of the infinite bulk medium in Fourier-Laplace domain.

[see PRB 91 085202 (2015)]

Member Function Documentation

◆ setLinTemporalGrid()

void alma::analytic1D::SPR_calculator_FourierLaplace::setLinTemporalGrid ( double  fmin,
double  fmax,
int  Nf 
)

Construct a linear grid of temporal frequencies from fmin to fmax with Nf elements.

From this a Laplace grid s = 2*pi*1i*f is created.

◆ setLogTemporalGrid()

void alma::analytic1D::SPR_calculator_FourierLaplace::setLogTemporalGrid ( double  fmin,
double  fmax,
int  Nf 
)

Construct a logarithmic grid of temporal frequencies from fmin to fmax with Nf elements.

From this a Laplace grid s = 2*pi*1i*f is created.


The documentation for this class was generated from the following files: