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

Class for computing basic thermal properties (kappa, Cv, diffusivity) and cumulative functions (resolved for MFP, energy, etc.) along a given transport direction. More...

#include <analytic1d.hpp>

Public Member Functions

 BasicProperties_calculator (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 setLogMFPbins (double MFPmin, double MFPmax, int Nbins)
 Construct a logarithmic grid of MFP bins from MFPmin to MFPmax with Nbins elements.
 
void setAutoMFPbins (int Nbins)
 Construct a logarithmic grid of MFP bins witn Nbins elements. More...
 
void setAutoProjMFPbins (int Nbins)
 Construct a logarithmic grid of MFP bins witn Nbins elements. More...
 
void setLogRTbins (double taumin, double taumax, int Nbins)
 Construct a logarithmic grid of tau bins from taumin to taumax with Nbins elements.
 
void setAutoRTbins (int Nbins)
 Construct a logarithmic grid of tau bins witn Nbins elements. More...
 
void setLinOmegabins (double omegamin, double omegamax, int Nbins)
 Construct a linear grid of omega bins from omegamin to omegamax with Nbins elements.
 
void setAutoOmegabins (int Nbins)
 Construct a linear grid of omega bins witn Nbins elements. More...
 
Eigen::VectorXd getMFPbins ()
 Retrieve the MFP bins.
 
Eigen::VectorXd getRTbins ()
 Retrieve the relaxation time bins.
 
Eigen::VectorXd getOmegabins ()
 Retrieve the relaxation time bins.
 
void setBulk ()
 Specify that the medium should be treated as infinite bulk (default option)
 
void setInPlaneFilm (double filmthickness, const Eigen::Vector3d normal, double specularity=0.0)
 Specify that the medium should be treated as a thin film with provided thickness. More...
 
void setCrossPlaneFilm (double filmthickness)
 
double getConductivity ()
 Retrieve thermal conductivity.
 
double getSpectralConductivity ()
 Obtain spectrally computed thermal conductivity.
 
double getCapacity ()
 Compute the heat capacity and retrieve it.
 
double getDiffusivity ()
 Compute the Fourier diffusivity and retrieve it.
 
double getDominantProjMFP ()
 Retrieve the dominant projected MFP as defined below.
 
double getDominantRT ()
 Retrieve the dominant phonon relaxation time.
 
double getAnisotropyIndex ()
 Compute the anisotropy index kappa_max/kappa_min and retrieve it.
 
void resolveByMFP ()
 Choose resolving cumulative curves by MFP.
 
void resolveByProjMFP ()
 Choose resolving cumulative curves by projected MFP.
 
void resolveByRT ()
 Choose resolving cumulative curves by relaxation time.
 
void resolveByOmega ()
 Choose resolving cumulative curves by angular frequency.
 
Eigen::VectorXd getCumulativeConductivity ()
 Compute the cumulative thermal conductivity curve with respect to the MFP bins and retrieve it.
 
Eigen::VectorXd getCumulativeCapacity ()
 Compute the cumulative heat capacity curve with respect to the MFP bins and retrieve it.
 
void setDOSgridsize (int nbins)
 Set the number of bins to be used in DOS evaluation.
 
Eigen::VectorXd getDOSgrid ()
 Obtain the DOS energy grid (in meV)
 
Eigen::VectorXd getDOS ()
 Calculate and obtain the DOS.
 

Detailed Description

Class for computing basic thermal properties (kappa, Cv, diffusivity) and cumulative functions (resolved for MFP, energy, etc.) along a given transport direction.

Calculations involve appropriate thin film corrections for both cross-plane and in-plane transport.

Member Function Documentation

◆ setAutoMFPbins()

void alma::analytic1D::BasicProperties_calculator::setAutoMFPbins ( int  Nbins)

Construct a logarithmic grid of MFP bins witn Nbins elements.

MFPmax is assigned automatically to twice the largest MFP. MFPmin is assigned to 1e-6*largest MFP.

◆ setAutoOmegabins()

void alma::analytic1D::BasicProperties_calculator::setAutoOmegabins ( int  Nbins)

Construct a linear grid of omega bins witn Nbins elements.

omegamax is assigned automatically to 1.1*largest angular frequency. omegamin is automatically set to 0.

◆ setAutoProjMFPbins()

void alma::analytic1D::BasicProperties_calculator::setAutoProjMFPbins ( int  Nbins)

Construct a logarithmic grid of MFP bins witn Nbins elements.

MFPmax is assigned automatically to the twice the largest projected MFP. MFPmin is assigned to 1e-6*largest projected MFP.

◆ setAutoRTbins()

void alma::analytic1D::BasicProperties_calculator::setAutoRTbins ( int  Nbins)

Construct a logarithmic grid of tau bins witn Nbins elements.

taumax is assigned automatically to the twice the largest relaxation time. taumin is assigned to 1e-6*largest relaxation time.

◆ setInPlaneFilm()

void alma::analytic1D::BasicProperties_calculator::setInPlaneFilm ( double  filmthickness,
const Eigen::Vector3d  normal,
double  specularity = 0.0 
)

Specify that the medium should be treated as a thin film with provided thickness.

This option corrects MFPs and relaxation times in conductivity calculations.


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