AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
All Classes Files Functions Variables Typedefs Enumerations Friends Pages
steady_montecarlo1d_powersource.cpp File Reference

1D Monte Carlo solver for multilayer structures that uses a heat source with prescribed power density at the top and heat sink at the bottom instead of isothermal reservoirs. 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_powersource.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...
 

Typedefs

typedef std::pair< std::size_t, std::size_t > idx_pair
 
typedef std::array< double, 3 > emissionID
 

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 > &grid, double target)
 Helper function that determines in which bin a particle is situated. 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

1D Monte Carlo solver for multilayer structures that uses a heat source with prescribed power density at the top and heat sink at the bottom instead of isothermal reservoirs.

Function Documentation

◆ 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 > &  grid,
double  target 
)
inline

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

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