AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
steady_montecarlo1d.cpp File Reference

VRMC solver based on the Peraud-Hadjiconstantinou method [APL 101, 153114 (2012)]. More...

#include <string>
#include <vector>
#include <map>
#include <iostream>
#include <fstream>
#include <algorithm>
#include <Eigen/Dense>
#include <boost/format.hpp>
#include <boost/filesystem.hpp>
#include <boost/serialization/map.hpp>
#include <boost/mpi.hpp>
#include <boost/property_tree/ptree.hpp>
#include <boost/property_tree/xml_parser.hpp>
#include <boost/math/distributions/gamma.hpp>
#include <randutils.hpp>
#include <cmakevars.hpp>
#include <constants.hpp>
#include <utilities.hpp>
#include <io_utils.hpp>
#include <structures.hpp>
#include <bulk_hdf5.hpp>
#include <sampling.hpp>
#include <bulk_properties.hpp>
#include <deviational_particle.hpp>
#include <isotopic_scattering.hpp>
#include <analytic1d.hpp>
Include dependency graph for steady_montecarlo1d.cpp:

Classes

class  Steady_1d_simulator
 Class that helps run simulations of general 1D structures sandwiched between two thermal reservoirs, in the steady-state regime. More...
 

Functions

double gamma_pdf (double k, double theta, double x)
 Probability density function for a Gamma distribution. More...
 
int get_bin_index (const Eigen::Ref< const Eigen::VectorXd > &xgrid, double xtarget)
 Helper function that determines in which position bin of the a particle is situated. More...
 
double calc_dotEeff_wall (const alma::Crystal_structure &poscar, const alma::Gamma_grid &grid, double Twall, double Tref, const Eigen::Ref< const Eigen::Vector3d > &normal)
 Compute the deviational power per unit area at an isothermal surface. More...
 
int main (int argc, char **argv)
 

Variables

boost::filesystem::path launch_path
 Working directory of the user.
 
std::string target_directory = "AUTO"
 Target directory for writing output.
 

Detailed Description

VRMC solver based on the Peraud-Hadjiconstantinou method [APL 101, 153114 (2012)].

Internally the solver runs multiple partial simulations. This allows us to provide the stochastic uncertainty on the reported output variables.

Function Documentation

◆ calc_dotEeff_wall()

double calc_dotEeff_wall ( const alma::Crystal_structure poscar,
const alma::Gamma_grid grid,
double  Twall,
double  Tref,
const Eigen::Ref< const Eigen::Vector3d > &  normal 
)
inline

Compute the deviational power per unit area at an isothermal surface.

Parameters
[in]poscar- a description of the unit cell
[in]grid- phonon spectrum on a regular grid
[in]Twall- the temperature of the wall in K
[in]Tref- the simulation temperature in K
[in]normal- a normal vector pointing out of the wall
Returns
the deviational power per unit area in units of kB / ps / nm ** 2

◆ gamma_pdf()

double gamma_pdf ( double  k,
double  theta,
double  x 
)

Probability density function for a Gamma distribution.

Parameters
[in]k- shape parameter
[in]theta- scale parameter
[in]x- point at which to evaluate the function
Returns
the value of the pdf

◆ get_bin_index()

int get_bin_index ( const Eigen::Ref< const Eigen::VectorXd > &  xgrid,
double  xtarget 
)
inline

Helper function that determines in which position bin of the a particle is situated.

Parameters
[in]xgrid- sorted grid of x values defining the bins [in] xtarget - position of the particle
Returns
the index of the bin where the particle is