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

Builds HDF5 files describing binary superlattices, based on XML input. More...

#include <map>
#include <set>
#include <string>
#include <vector>
#include <cstdint>
#include <iostream>
#include <fstream>
#include <algorithm>
#include <boost/format.hpp>
#include <boost/filesystem.hpp>
#include <boost/log/core.hpp>
#include <boost/log/trivial.hpp>
#include <boost/log/expressions.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/iostreams/stream_buffer.hpp>
#include <boost/iostreams/device/null.hpp>
#include <boost/uuid/sha1.hpp>
#include <boost/endian/conversion.hpp>
#include <Eigen/Dense>
#include <basen.hpp>
#include <randutils.hpp>
#include <cmakevars.hpp>
#include <constants.hpp>
#include <utilities.hpp>
#include <io_utils.hpp>
#include <structures.hpp>
#include <vc.hpp>
#include <vasp_io.hpp>
#include <bulk_hdf5.hpp>
#include <sampling.hpp>
#include <bulk_properties.hpp>
#include <deviational_particle.hpp>
#include <isotopic_scattering.hpp>
#include <superlattices.hpp>
Include dependency graph for superlattice_builder.cpp:

Functions

std::string get_strprofile (const std::vector< double > &input)
 Return a single string describing a composition profile. More...
 
std::string get_uid (const std::vector< double > &input)
 Return a reproducible short identifier obtained from a set of doubles. More...
 
int main (int argc, char **argv)
 

Variables

boost::iostreams::stream_buffer< boost::iostreams::null_sink > null_buf
 
std::vector< std::string > materialBase
 
int gridDensityA
 
int gridDensityB
 
int gridDensityC
 
std::string materials_repository = "."
 
std::string target_directory = "AUTO"
 
std::string target_filename = "AUTO"
 
int materialType = -1
 
std::vector< double > mixfractions
 
Eigen::Vector3i normal
 
std::size_t nqline
 
bool do3ph = true
 
bool overwrite = false
 

Detailed Description

Builds HDF5 files describing binary superlattices, based on XML input.

Function Documentation

◆ get_strprofile()

std::string get_strprofile ( const std::vector< double > &  input)

Return a single string describing a composition profile.

Parameters
[in]input- vector of doubles
Returns
the string built by representing each concentration as a string with %.6f and separating them by commas.

◆ get_uid()

std::string get_uid ( const std::vector< double > &  input)

Return a reproducible short identifier obtained from a set of doubles.

The identifier is defined as the first 8 characters of the base32-encoded SHA1 hash of the string built by representing each concentration as a string with %.6f and separating them by commas.

Parameters
[in]input- vector of doubles
Returns
the identifier as a string

Variable Documentation

◆ null_buf

boost::iostreams::stream_buffer<boost::iostreams::null_sink> null_buf
Initial value:
{
boost::iostreams::null_sink()}