24 #include <Eigen/Dense> 28 namespace analytic1D {
41 const Eigen::ArrayXXd* w,
63 void setLogRTbins(
double taumin,
double taumax,
int Nbins);
98 const Eigen::Vector3d normal,
99 double specularity = 0.0);
100 void setCrossPlaneFilm(
double filmthickness);
158 const Eigen::ArrayXXd* w;
163 Eigen::Vector3d unitvector;
165 Eigen::Vector3d normalvector;
167 Eigen::VectorXd MFPbins;
169 Eigen::VectorXd taubins;
171 Eigen::VectorXd omegabins;
173 double internal_MFPmax;
175 double internal_MFPprojmax;
177 double internal_taumax;
179 double internal_omegamax;
183 double filmthickness;
190 double dominantProjMFP;
197 Eigen::Matrix3d FuchsRotation;
199 void initFuchsRotation();
203 Eigen::VectorXd kappacumul;
205 Eigen::VectorXd Cvcumul;
214 int kappacumulIdentifier;
231 const Eigen::ArrayXXd* w,
239 void setLinGrid(
double ximin,
double ximax,
int Nxi);
243 void setLogGrid(
double ximin,
double ximax,
int Nxi);
246 void setXiGrid(
const Eigen::Ref<const Eigen::VectorXd> xigrid);
249 Eigen::VectorXd getSpatialFrequencies();
253 void normaliseOutput(
bool norm);
259 Eigen::VectorXd getPsi();
267 const Eigen::ArrayXXd* w;
272 Eigen::Vector3d unitvector;
279 void updateDiffusivity();
295 const Eigen::ArrayXXd* w,
303 void setLinSpatialGrid(
double ximin,
double ximax,
int Nxi);
307 void setLogSpatialGrid(
double ximin,
double ximax,
int Nxi);
312 void setLinTemporalGrid(
double fmin,
double fmax,
int Nf);
317 void setLogTemporalGrid(
double fmin,
double fmax,
int Nf);
320 Eigen::VectorXd getSpatialFrequencies();
323 Eigen::VectorXd getTemporalFrequencies();
326 Eigen::MatrixXcd getSPR();
334 const Eigen::ArrayXXd* w;
339 Eigen::Vector3d unitvector;
346 Eigen::MatrixXcd Pxis;
367 const Eigen::ArrayXXd* w,
375 void setLinGrid(
double xmin,
double xmax,
int Nx);
379 void setLogGrid(
double xmin,
double xmax,
int Nx);
385 void declareGridNormalised(
bool norm);
389 void normaliseOutput(
bool norm);
392 void setTime(
double t);
398 Eigen::VectorXd getGrid();
399 Eigen::VectorXd getNormalisedGrid();
411 Eigen::VectorXd getSourceTransient(
412 const Eigen::Ref<Eigen::VectorXd> timegrid);
415 Eigen::VectorXd getSPR();
418 Eigen::MatrixXd resolveSPRbyMFP();
426 const Eigen::ArrayXXd* w;
431 Eigen::Vector3d unitvector;
433 Eigen::Vector3d filmnormal;
436 bool gridIsNormalised;
438 bool normaliseoutput;
440 Eigen::VectorXd MFPbins;
445 void updateDiffusivity();
449 Eigen::MatrixXd Pmodes;
470 const Eigen::ArrayXXd* w,
478 void setLinGrid(
double fmin,
double fmax,
int Nf);
482 void setLogGrid(
double fmin,
double fmax,
int Nf);
485 void setLaplaceGrid(
const Eigen::Ref<const Eigen::VectorXcd> sgrid);
488 Eigen::VectorXd getGrid();
489 Eigen::VectorXcd getLaplaceGrid();
492 Eigen::VectorXcd getMSD();
500 const Eigen::ArrayXXd* w;
505 Eigen::Vector3d unitvector;
511 Eigen::VectorXcd MSD;
528 const Eigen::ArrayXXd* w,
536 void setLinGrid(
double tmin,
double tmax,
int Nt);
540 void setLogGrid(
double tmin,
double tmax,
int Nt);
543 void setTimeGrid(
const Eigen::Ref<const Eigen::VectorXd> tgrid);
546 Eigen::VectorXd getGrid();
550 void normaliseOutput(
bool norm);
556 Eigen::VectorXd getMSD();
564 const Eigen::ArrayXXd* w;
569 Eigen::Vector3d unitvector;
571 Eigen::Vector3d filmnormal;
576 void updateDiffusivity();
578 bool normaliseoutput;
584 Eigen::VectorXd GS_coeffs;
586 double factorial(
int n);
Eigen::VectorXd getDOSgrid()
Obtain the DOS energy grid (in meV)
Definition: analytic1d.cpp:434
void setDOSgridsize(int nbins)
Set the number of bins to be used in DOS evaluation.
Definition: analytic1d.cpp:430
Class for computing the RTA propagator function psi(xi) of a medium.
Definition: analytic1d.hpp:226
Definition: analytic1d.hpp:26
void setAutoProjMFPbins(int Nbins)
Construct a logarithmic grid of MFP bins witn Nbins elements.
Definition: analytic1d.cpp:72
Class for computing basic thermal properties (kappa, Cv, diffusivity) and cumulative functions (resol...
Definition: analytic1d.hpp:36
Class for computing the exact analytical RTA single pulse energy density response of the infinite bul...
Definition: analytic1d.hpp:290
Class for computing the exact analytical RTA solution for mean square thermal energy displacement...
Definition: analytic1d.hpp:465
double getConductivity()
Retrieve thermal conductivity.
Definition: analytic1d.cpp:492
void resolveByRT()
Choose resolving cumulative curves by relaxation time.
Definition: analytic1d.cpp:564
void setLogRTbins(double taumin, double taumax, int Nbins)
Construct a logarithmic grid of tau bins from taumin to taumax with Nbins elements.
Definition: analytic1d.cpp:88
double getDominantRT()
Retrieve the dominant phonon relaxation time.
Definition: analytic1d.cpp:541
void setAutoMFPbins(int Nbins)
Construct a logarithmic grid of MFP bins witn Nbins elements.
Definition: analytic1d.cpp:78
double getDominantProjMFP()
Retrieve the dominant projected MFP as defined below.
Definition: analytic1d.cpp:526
Eigen::VectorXd getCumulativeCapacity()
Compute the cumulative heat capacity curve with respect to the MFP bins and retrieve it...
Definition: analytic1d.cpp:744
Eigen::VectorXd getDOS()
Calculate and obtain the DOS.
Definition: analytic1d.cpp:442
Classes and functions used to manipulate grids in reciprocal space.
void setAutoOmegabins(int Nbins)
Construct a linear grid of omega bins witn Nbins elements.
Definition: analytic1d.cpp:100
Eigen::VectorXd getRTbins()
Retrieve the relaxation time bins.
Definition: analytic1d.cpp:109
Eigen::VectorXd getMFPbins()
Retrieve the MFP bins.
Definition: analytic1d.cpp:105
void setDirection(const Eigen::Vector3d unitvector)
Function setting the thermal transport axis.
Definition: analytic1d.cpp:60
void setLogMFPbins(double MFPmin, double MFPmax, int Nbins)
Construct a logarithmic grid of MFP bins from MFPmin to MFPmax with Nbins elements.
Definition: analytic1d.cpp:66
double getSpectralConductivity()
Obtain spectrally computed thermal conductivity.
Definition: analytic1d.cpp:312
Eigen::VectorXd getOmegabins()
Retrieve the relaxation time bins.
Definition: analytic1d.cpp:113
void resolveByProjMFP()
Choose resolving cumulative curves by projected MFP.
Definition: analytic1d.cpp:560
Class for computing the exact analytical RTA solution for mean square thermal energy displacement in ...
Definition: analytic1d.hpp:523
void setAutoRTbins(int Nbins)
Construct a logarithmic grid of tau bins witn Nbins elements.
Definition: analytic1d.cpp:83
Eigen::VectorXd getCumulativeConductivity()
Compute the cumulative thermal conductivity curve with respect to the MFP bins and retrieve it...
Definition: analytic1d.cpp:572
BasicProperties_calculator(const alma::Crystal_structure *poscar, const alma::Gamma_grid *grid, const Eigen::ArrayXXd *w, double T)
Constructor: initialise internal variables.
Definition: analytic1d.cpp:36
void resolveByOmega()
Choose resolving cumulative curves by angular frequency.
Definition: analytic1d.cpp:568
double getDiffusivity()
Compute the Fourier diffusivity and retrieve it.
Definition: analytic1d.cpp:511
Definitions of the basic data-handling classes in ALMA.
Objects of this class represent a regular grid with the Gamma point in one corner.
Definition: qpoint_grid.hpp:32
Hold information about a crystal structure.
Definition: structures.hpp:51
void setBulk()
Specify that the medium should be treated as infinite bulk (default option)
Definition: analytic1d.cpp:118
Physical, mathematical and miscellaneous constants used in alma.
Class for computing the approximate analytical RTA single pulse energy density response in real space...
Definition: analytic1d.hpp:362
void resolveByMFP()
Choose resolving cumulative curves by MFP.
Definition: analytic1d.cpp:556
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. ...
Definition: analytic1d.cpp:130
double getAnisotropyIndex()
Compute the anisotropy index kappa_max/kappa_min and retrieve it.
Definition: analytic1d.cpp:483
double getCapacity()
Compute the heat capacity and retrieve it.
Definition: analytic1d.cpp:507
void setLinOmegabins(double omegamin, double omegamax, int Nbins)
Construct a linear grid of omega bins from omegamin to omegamax with Nbins elements.
Definition: analytic1d.cpp:94