High Performance Plasticity  0.5.0
hpp Namespace Reference

Classes

class  Crystal
 
class  CrystalError
 
struct  CrystalInitialConditions
 
class  CrystalOutputConfig
 
class  CrystalProperties
 
class  CrystalPropertiesCUDA
 
struct  CrystalSolverConfig
 
class  EulerAngles
 
class  Experiment
 
struct  FFTWConfigRealND
 
class  FreeDelete
 
class  GridOrientationGenerator
 Generates orientations on a fixed grid over an azimuthal angle of \( [0, 2\pi) \) and a zenithal angle of \( [0, \pi/2) \). Does not cover the area of a sphere equally. More...
 
class  GSHCoeffs
 
class  GSHCoeffsCUDA
 
class  HDF5Handler
 
struct  hdf_complex_t
 A complex number type. More...
 
struct  HDFReadWriteParams
 
struct  HDFReadWriteParamsC
 
class  HDFUtilsError
 
struct  idx2d
 A 2D index. More...
 
struct  idx4d
 A 4D index. More...
 
class  OrientationGenerator
 An abstract class for reproducibly generating a sequence of orientations. More...
 
struct  PolarDecomposition
 
class  Polycrystal
 
class  PolycrystalOutputConfig
 
class  RandomOrientationGenerator
 A class for reproducibly generating a sequence of random orientations. More...
 
class  SO3Discrete
 
struct  SpectralCoeffCUDA
 
class  SpectralCoordCUDA
 
class  SpectralCrystal
 
struct  SpectralCrystalCUDA
 
class  SpectralCrystalListCUDA
 
struct  SpectralCrystalSolverConfig
 
class  SpectralDatabase
 
class  SpectralDatabaseCUDA
 
class  SpectralDatabaseUnified
 
class  SpectralDatabaseUnifiedCUDA
 
struct  SpectralDataCUDA
 
class  SpectralDataset
 
class  SpectralDatasetCUDA
 
class  SpectralDatasetID
 
class  SpectralPolycrystal
 
class  SpectralPolycrystalCUDA
 
class  SpectralPolycrystalGSHCUDA
 
struct  StretchingTensorDecomposition
 
class  Tensor2
 A class for second order tensors. More...
 
class  Tensor2AsymmCUDA
 
class  Tensor2CUDA
 
class  Tensor2SymmCUDA
 
class  Tensor4
 
class  Tensor4CUDA
 
class  TensorError
 
class  Timer
 
class  VecCUDA
 

Enumerations

enum  CrystalType { CRYSTAL_TYPE_NONE, CRYSTAL_TYPE_FCC }
 
enum  HardeningLaw { HARDENING_LAW_BROWN, HARDENING_LAW_VOCE }
 
enum  CrystalDatasetIdx {
  SIGMA00, SIGMA11, SIGMA12, SIGMA02,
  SIGMA01, WP01, WP12, WP02,
  GAMMA
}
 
enum  HDFReadWrite : char { HDFReadWrite::Write ='w', HDFReadWrite::Read ='r' }
 
enum  SymmetryType { SYMMETRY_TYPE_NONE, SYMMETRY_TYPE_C4 }
 

Functions

std::vector< SpectralDatasetIDdefaultCrystalSpectralDatasetIDs ()
 
template<typename U >
hpp::Tensor2< U > tensorA (hpp::Tensor2< U > F_p, hpp::Tensor2< U > F_next)
 See kalidindi1992.TensorA. More...
 
template<typename U >
std::vector< hpp::Tensor2< U > > tensorC_alphas (const hpp::Tensor4< U > &L, const hpp::Tensor2< U > &A, const std::vector< hpp::Tensor2< U >> &S_0, const unsigned int n_alpha)
 See kalidindi1992.TensorC_alphas. More...
 
template<typename U >
void tensorC_alphasInPlace (const hpp::Tensor4< U > &L, const hpp::Tensor2< U > &A, const std::vector< hpp::Tensor2< U >> &S_0, const unsigned int n_alpha, hpp::Tensor2< U > &dumB_alpha, std::vector< hpp::Tensor2< U >> &C_alphas)
 In-place version of hpp.tensorC_alphas. More...
 
template<typename U >
plasticShearingRate (U tau_alpha, U s_alpha, U gammadot_0, U m)
 See kalidindi1992.PlasticShearingRate. More...
 
template<typename U >
std::vector< U > shearStrainRates (const CrystalProperties< U > &props, const hpp::Tensor2< U > &T, const std::vector< U > &s_alphas)
 See kalidindi1992.ShearStrainRates. More...
 
template<typename U >
void shearStrainRatesInPlace (const CrystalProperties< U > &props, const hpp::Tensor2< U > &T, const std::vector< U > &s_alphas, std::vector< U > &gammadot_alphas)
 In-place version of hpp.shearStrainRates. More...
 
template<typename U >
std::vector< U > shearStrainIncrements (const CrystalProperties< U > &props, const hpp::Tensor2< U > &T, const std::vector< U > &s_alphas, const U dt)
 See kalidindi1992.ShearStrainIncrements. More...
 
template<typename U >
void shearStrainIncrementsInPlace (const CrystalProperties< U > &props, const hpp::Tensor2< U > &T, const std::vector< U > &s_alphas, const U dt, std::vector< U > &Dgamma_alphas)
 In-place version of hpp.shearStrainIncrements. More...
 
template<typename U >
hpp::Tensor2< U > partialDGammaPartialT (U m, U gammadot_0, const hpp::Tensor2< U > &T, const hpp::Tensor2< U > &S_0_alpha, U s_alpha, U dt)
 See kalidindi1992.PartialDGammaPartialT. More...
 
template<typename U >
void partialDGammaPartialTInPlace (U m, U gammadot_0, const hpp::Tensor2< U > &T, const hpp::Tensor2< U > &S_0_alpha, U s_alpha, U dt, hpp::Tensor2< U > &pdgpt)
 In-place version of hpp.partialDGammaPartialT. More...
 
template<typename U >
hpp::Tensor4< U > tensorJ (const CrystalProperties< U > &props, const std::vector< hpp::Tensor2< U >> &C_alphas, const std::vector< U > &s_alphas, const hpp::Tensor2< U > &T, const U dt)
 See kalidindi1992.TensorJ. More...
 
template<typename U >
void tensorJInPlace (const CrystalProperties< U > &props, const std::vector< hpp::Tensor2< U >> &C_alphas, const std::vector< U > &s_alphas, const hpp::Tensor2< U > &T, const U dt, hpp::Tensor2< U > &dum2ndOrder, hpp::Tensor4< U > &dum4thOrder, hpp::Tensor4< U > &J)
 In-place version of hpp.tensorJ. More...
 
template<typename U >
hpp::Tensor2< U > tensorG (const hpp::Tensor4< U > &L, const hpp::Tensor2< U > &A, const hpp::Tensor2< U > &T_prev_iter, const std::vector< U > &Dgamma_alphas, const std::vector< hpp::Tensor2< U >> &C_alphas)
 See kalidindi1992.TensorG. More...
 
template<typename U >
void tensorGInPlace (const hpp::Tensor4< U > &L, const hpp::Tensor2< U > &A, const hpp::Tensor2< U > &T_prev_iter, const std::vector< U > &Dgamma_alphas, const std::vector< hpp::Tensor2< U >> &C_alphas, Tensor2< U > &dum2ndOrder, hpp::Tensor2< U > &dumT_tr, hpp::Tensor2< U > &G)
 In-place version of hpp.TensorG. More...
 
template<typename U >
hpp::Tensor2< U > newtonStressCorrection (const CrystalProperties< U > &props, const hpp::Tensor2< U > &A, const std::vector< hpp::Tensor2< U >> &C_alphas, const std::vector< U > &s_alphas, const hpp::Tensor2< U > &T_old, const U DT_max, const U dt)
 See kalidindi1992.NewtonStressCorrection. More...
 
template<typename U >
void newtonStressCorrectionInPlace (const CrystalProperties< U > &props, const hpp::Tensor2< U > &A, const std::vector< hpp::Tensor2< U >> &C_alphas, const std::vector< U > &s_alphas, const hpp::Tensor2< U > &T_old, const U DT_max, const U dt, std::vector< U > &dumDgamma_alphas, hpp::Tensor4< U > &dumJ, hpp::Tensor2< U > &dumT_tr, hpp::Tensor2< U > &dumG, hpp::Tensor2< U > &dum2ndOrder, hpp::Tensor4< U > &dum4thOrder, hpp::Tensor2< U > &DT)
 In-place version of hpp.newtonStressCorrection. More...
 
template<typename U >
slipHardeningRate (U h_0, U s_s, U a, U s_beta)
 
template<typename U >
hpp::Tensor2< U > strainHardeningRates (const CrystalProperties< U > &props, const std::vector< U > &s_alphas)
 
template<typename U >
strainHardeningRateVoce (const CrystalProperties< U > &props, const std::vector< U > &s_alphas)
 Voce hardening law. More...
 
template<typename U >
std::vector< U > slipDeformationResistanceUpdate (const CrystalProperties< U > &props, const hpp::Tensor2< U > &T_next, const std::vector< U > &s_alphas_current_time, const std::vector< U > &s_alphas_prev_iter, const U dt)
 
template<typename U >
slipDeformationResistanceStepSpectralSolver (const CrystalProperties< U > &props, const U s_alpha, const U gammaSum, const U dt)
 From kalidindi2006 Equation 5. More...
 
template<typename U >
hpp::Tensor2< U > plasticDeformationGradientUpdate (const CrystalProperties< U > &props, const hpp::Tensor2< U > &F_p_prev_time, const std::vector< U > &Dgamma_alphas)
 
template<typename U >
Tensor2< U > plasticSpinTensor (CrystalProperties< U > mprops, std::vector< U > gammadot_alphas, std::vector< std::vector< U >> m_alphas, std::vector< std::vector< U >> n_alphas)
 Mihaila 2014 equation 14. More...
 
template<typename U >
bool strainRateLowEnough (const std::vector< U > &Dgamma_alphas, U Dgamma_goal)
 
template<typename U >
setTimestepByShearRate (U dt_old, const std::vector< U > &Dgamma_alphas, U Dgamma_goal)
 The exact method in kalidindi1992. More...
 
template<typename T >
std::vector< T > cartesianToSpherical (const std::vector< T > &cartVec)
 
template<typename T >
Tensor2< T > histogramPoleEqualArea (const std::vector< EulerAngles< T >> &angles, const std::vector< T > &planeNormal)
 
template<typename T >
void writePoleHistogramHistoryHDF5 (H5::H5File &outfile, const std::string &dsetBaseName, std::vector< Tensor2< T >> &history, const std::vector< T > &pole)
 
void hdfSuppressErrors ()
 
void hdfRestoreErrors ()
 
hid_t getComplexHDFType ()
 
std::vector< hsize_t > getDatasetDims (hid_t dset_id)
 
bool MPIAllTrue (bool condition, MPI_Comm comm)
 Determine if condition is true for all processes. More...
 
FFTWConfigRealND prepareFFTWConfigRealND (const std::vector< ptrdiff_t > &realDims, MPI_Comm comm)
 
void destroyConfigRealND (FFTWConfigRealND &cfg)
 
unsigned int roundUpToMultiple (unsigned int val, unsigned int multiple)
 
bool operator< (const SpectralDatasetID &l, const SpectralDatasetID &r)
 
bool operator== (const SpectralDatasetID &l, const SpectralDatasetID &r)
 
void evaluateSpectralCompressionErrorFull (std::string rawDBName, std::string spectralDBName, std::string errorDBName, unsigned int nTermsMax, std::string outFilename, MPI_Comm comm)
 
void evaluateSpectralCompressionErrorFullUnified (std::string rawDBName, std::string spectralDBName, std::string errorDBName, unsigned int nTermsMax, std::string outFilename, MPI_Comm comm)
 
void evaluateSpectralCompressionErrorAxisSlice (std::string rawDBName, std::string spectralDBName, std::vector< std::string > dsetBasenames, unsigned int nTermsMax, unsigned int refineMult, std::vector< int > axisSlice, std::string outFilename, MPI_Comm comm)
 Evaluate the error incurred in the spectral reresentation of a database. More...
 
void evaluateSpectralCompressionErrorAxisSliceUnified (std::string rawDBName, std::string spectralDBName, std::vector< std::string > dsetBasenames, unsigned int nTermsMax, unsigned int refineMult, std::vector< int > axisSlice, std::string outFilename, MPI_Comm comm)
 Evaluate the error incurred in the spectral reresentation of a database. More...
 
template<typename U >
Tensor2< U > simpleShearDeformationGradient (U t, U shear_rate)
 The deformation gradient from a simple shear. More...
 
template<typename U >
Tensor2< U > simpleShearVelocityGradient (U t, U shear_rate)
 The velocity gradient from a simple shear. More...
 
template<typename U >
Tensor2< U > simpleCompressionDeformationGradient (U t, U comp_rate)
 The deformation gradient from a simple compression. More...
 
template<typename U >
Tensor2< U > simpleCompressionVelocityGradient (U t, U comp_rate)
 The velocity gradient from a simple compression. More...
 
template<typename U >
Tensor2< U > planeStrainCompressionDeformationGradient (U t, U comp_rate)
 The deformation gradient from a plane strain compression. More...
 
template<typename U >
Tensor2< U > planeStrainCompressionVelocityGradient (U t, U comp_rate)
 The velocity gradient from a plane strain compression. More...
 
template<typename U >
Tensor2< U > staticDeformationGradient (U t)
 The deformation gradient for no deformation. More...
 
template<typename U >
Tensor2< U > staticVelocityGradient (U t)
 The velocity gradient for no deformation. More...
 
template<typename U >
Tensor2< U > deformationGradFromVelGrad (U t_init, const Tensor2< U > &F_init, const Tensor2< U > &L, U t)
 
template<typename U >
Tensor2< U > stretchingVelocityGradient (U theta, U eDot)
 
template<typename T >
StretchingTensorDecomposition< T > getStretchingTensorDecomposition (const Tensor2< T > &L)
 
template<typename U >
Tensor4< U > isotropicElasticityTensor (U mu, U kappa)
 Construct a fourth order isotropic elasticity tensor. More...
 
template<typename U >
Tensor4< U > cubeSymmetricElasticityTensor (U c11, U c12, U c44)
 Construct a fourth order elasticity tensor accounting for cubic symmetry. More...
 
template<typename U >
Tensor4< U > volumeAverage (const std::vector< Tensor4< U >> &tVec, const std::vector< U > &vVec)
 Volume average fourth order tensors. More...
 
template<typename U >
Tensor4< U > getEshelbyTensorCubicMaterialSphericalInclusion (U c11, U c12, U c44, U I0, U I1, U I2)
 
template<typename U >
Tensor4< U > getHomogenizedStiffnessVolumeAverage (const std::vector< Tensor4< U >> &cVec, const std::vector< U > &vVec)
 Get the homogenized stiffness tensor using a volume average. More...
 
template<typename U >
Tensor4< U > getHomogenizedStiffnessELSC (const std::vector< Tensor4< U >> &cVec, const std::vector< U > &vVec, const Tensor4< U > &S)
 Get the homogenized stiffness tensor using the elastic self-consistent method. More...
 
constexpr int nSlipSystems (CrystalType crystalType)
 
template<typename U >
CrystalProperties< U > defaultCrystalProperties ()
 
template<typename U >
CrystalProperties< U > rotate (const CrystalProperties< U > &propsOld, hpp::Tensor2< U > rotTensor)
 
template<typename U >
CrystalSolverConfig< U > defaultCrystalSolverConfig ()
 
template<typename U >
CrystalSolverConfig< U > defaultConservativeCrystalSolverConfig ()
 
template<typename U >
CrystalInitialConditions< U > defaultCrystalInitialConditions ()
 
template<typename T >
EulerAngles< T > getEulerAnglesFromDeformationGradient (const hpp::Tensor2< T > &F_star)
 Get the corresponding Euler angles for the rotation induced by this deformation. More...
 
template<typename U >
bool operator== (const Crystal< U > &l, const Crystal< U > &r)
 
template<typename T >
bool operator== (const SpectralCrystalCUDA< T > &l, const SpectralCrystalCUDA< T > &r)
 
template<typename T , unsigned int N>
__device__ T slipDeformationResistanceStepSpectralSolver (const CrystalPropertiesCUDA< T, N > *props, const T s_alpha, const T gammaSum, const T dt)
 
template<typename T , unsigned int P>
__device__ void getSpectralCrystalDatabaseCoordinate (SpectralCrystalCUDA< T > &crystal, SpectralDatabaseUnifiedCUDA< T, 4, P > *db, Tensor2CUDA< T, 3, 3 > &RStretchingTensor, T theta, unsigned int *dbCoord)
 
template<typename T , unsigned int N>
__global__ void SPECTRAL_POLYCRYSTAL_STEP (unsigned int nCrystals, SpectralCrystalCUDA< T > *crystals, CrystalPropertiesCUDA< T, N > *props, Tensor2CUDA< T, 3, 3 > RStretchingTensor, Tensor2AsymmCUDA< T, 3 > WNext, T theta, T eDot, T dt, SpectralDatabaseCUDA< T, 4 > *db, Tensor2CUDA< T, 3, 3 > *TCauchyPerBlockSums)
 
template<typename T , unsigned int N, unsigned int P>
__global__ void SPECTRAL_POLYCRYSTAL_STEP_UNIFIED (unsigned int nCrystals, SpectralCrystalCUDA< T > *crystals, CrystalPropertiesCUDA< T, N > *props, Tensor2CUDA< T, 3, 3 > RStretchingTensor, Tensor2AsymmCUDA< T, 3 > WNext, T theta, T eDot, T dt, SpectralDatabaseUnifiedCUDA< T, 4, P > *db, Tensor2CUDA< T, 3, 3 > *TCauchyPerBlockSums)
 
template<typename T >
__global__ void GET_GSH_COEFFS (const SpectralCrystalCUDA< T > *crystals, unsigned int ncrystals, GSHCoeffsCUDA< T > *coeffsPerBlockSums)
 
template<typename T >
__global__ void GET_AVERAGE_TCAUCHY (unsigned int nCrystals, const SpectralCrystalCUDA< T > *crystals, Tensor2CUDA< T, 3, 3 > *TCauchyGlobal)
 
template<typename T >
__device__ Tensor2CUDA< T, 3, 3 > EulerZXZRotationMatrixCUDA (EulerAngles< T > angles)
 
template<typename T >
GSHCoeffs< T > operator+ (const GSHCoeffs< T > &coeffs1, const GSHCoeffs< T > &coeffs2)
 
template<typename T >
void operator+= (GSHCoeffs< T > &A, const GSHCoeffs< T > &B)
 
template<typename T >
GSHCoeffs< T > operator/ (const GSHCoeffs< T > &coeffs, const T val)
 
template<typename T >
void operator/= (GSHCoeffs< T > &A, const T B)
 
template<typename T >
__host__ __device__ GSHCoeffsCUDA< T > operator+ (const GSHCoeffsCUDA< T > &coeffs1, const GSHCoeffsCUDA< T > &coeffs2)
 
template<typename T >
__host__ __device__ void operator+= (GSHCoeffsCUDA< T > &A, const GSHCoeffsCUDA< T > &B)
 
template<typename T >
__host__ __device__ GSHCoeffsCUDA< T > operator/ (const GSHCoeffsCUDA< T > &coeffs, T val)
 
template<typename T >
__device__ GSHCoeffsCUDA< T > warpReduceSumGSHCoeffs (GSHCoeffsCUDA< T > coeffs)
 
template<typename T >
__device__ GSHCoeffsCUDA< T > blockReduceSumGSHCoeffs (GSHCoeffsCUDA< T > val)
 
template<typename T >
__global__ void BLOCK_REDUCE_KEPLER_GSH_COEFFS (GSHCoeffsCUDA< T > *in, GSHCoeffsCUDA< T > *out, int nTerms)
 
void hdfAssert (herr_t code, const char *file, int line)
 
template<typename T >
hid_t getHDF5TypeC ()
 Get HDF5 equivalent type. More...
 
template<typename T >
hid_t createHDF5Dataset (hid_t file_id, std::string datasetName, std::vector< hsize_t > dataDims)
 
template<typename T >
hid_t createHDF5GridOfArrays (hid_t file_id, std::string datasetName, std::vector< hsize_t > gridDims, std::vector< hsize_t > arrayDims)
 
template<typename T >
HDFReadWriteParamsC getReadWriteParametersForMultipleHDF5Arrays (hid_t dset_id, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims, std::vector< hsize_t > arrayCount)
 Gets the parameters and ensures that they're compatible with the actual dataset. More...
 
template<typename T >
HDFReadWriteParamsC getReadWriteParametersForSingleHDF5Array (hid_t dset_id, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims)
 Gets the parameters and ensures that they're compatible with the actual dataset. More...
 
template<typename T >
void readWriteHDF5SimpleArray (hid_t dset_id, hid_t plist_id, HDFReadWriteParamsC parms, T *output, HDFReadWrite mode)
 
template<typename T >
void readHDF5SimpleArray (hid_t dset_id, hid_t plist_id, HDFReadWriteParamsC parms, T *output)
 
template<typename T >
void writeHDF5SimpleArray (hid_t dset_id, hid_t plist_id, HDFReadWriteParamsC parms, T *output)
 
template<typename T >
void readSingleHDF5Array (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims, T *output)
 
template<typename T >
void readSingleHDF5Array (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > arrayDims, T *output)
 
template<typename T >
void writeSingleHDF5Array (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims, T *output)
 
template<typename T >
void writeSingleHDF5Array (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > arrayDims, T *output)
 
template<typename T >
void writeMultipleHDF5Arrays (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims, std::vector< hsize_t > arrayCount, T *output)
 
template<typename T >
void readMultipleHDF5Arrays (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims, std::vector< hsize_t > arrayCount, T *output)
 
template<typename T >
void readSingleHDF5Value (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > gridOffset, T *output)
 
template<typename T >
void writeSingleHDF5Value (hid_t dset_id, hid_t plist_id, std::vector< hsize_t > gridOffset, T *output)
 
template<typename T >
H5::DataType getHDF5Type ()
 Get HDF5 equivalent type. More...
 
template<typename T >
H5::DataSet createHDF5Dataset (H5::H5File &file, const H5std_string &datasetName, std::vector< hsize_t > dataDims)
 
template<typename T >
H5::DataSet createHDF5GridOfArrays (H5::H5File &file, const H5std_string &datasetName, std::vector< hsize_t > gridDims, std::vector< hsize_t > arrayDims)
 
template<typename T >
HDFReadWriteParams getReadWriteParametersForSingleHDF5Array (H5::DataSet &dataset, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims)
 
template<typename T >
void readWriteHDF5SimpleArray (H5::DataSet &dataset, HDFReadWriteParams parms, T *output, HDFReadWrite mode)
 
template<typename T >
void readHDF5SimpleArray (H5::DataSet &dataset, HDFReadWriteParams parms, T *output)
 
template<typename T >
void writeHDF5SimpleArray (H5::DataSet &dataset, HDFReadWriteParams parms, T *output)
 
template<typename T >
void readSingleHDF5Array (H5::DataSet &dataset, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims, T *output)
 
template<typename T >
void writeSingleHDF5Array (H5::DataSet &dataset, std::vector< hsize_t > gridOffset, std::vector< hsize_t > arrayDims, T *output)
 
template<typename T >
void writeVectorToHDF5Array (H5::H5File &file, const std::string &dsetName, std::vector< T > &vec)
 
template<typename T >
void addAttribute (H5::H5File &file, const std::string &attrName, T attrVal)
 
template<typename T >
MPI_Datatype MPIType ()
 Get the MPI datatype from the C++ type. More...
 
template<typename T >
MPIMin (const T &localVal, MPI_Comm comm)
 
template<typename T >
MPIMax (const T &localVal, MPI_Comm comm)
 
template<typename T >
MPISum (const T &localVal, MPI_Comm comm)
 
template<typename T >
std::vector< T > MPIConcatOnRoot (T localVal, MPI_Comm comm)
 
template<typename T >
std::vector< T > MPIConcatOnRoot (std::vector< T > localVec, MPI_Comm comm)
 
template<typename T >
MPIBroadcastFromRoot (T rootVal, MPI_Comm comm)
 
template<typename T >
std::vector< T > MPIBroadcastFromRoot (std::vector< T > rootVec, MPI_Comm comm)
 
template<typename T >
std::vector< T > MPISplitVectorEvenly (const std::vector< T > &rootVec, MPI_Comm comm)
 
template<typename T >
std::vector< T > MPISplitVectorEvenly (const std::vector< T > &rootVec, MPI_Comm comm, MPI_Datatype dtype)
 
template<class T >
boost::python::list toPythonList (const std::vector< T > &vec)
 
template<class T >
std::vector< T > toStdVector (const boost::python::list &list)
 
template<typename T >
MPI_Datatype getEulerAnglesTypeMPI ()
 
template<typename T >
EulerAngles< T > polarToEuler (T theta, T phi)
 Convert polar angles to Euler angles. More...
 
template<typename T >
std::ostream & operator<< (std::ostream &out, const EulerAngles< T > &angles)
 
template<typename T >
bool operator== (const EulerAngles< T > &l, const EulerAngles< T > &r)
 
template<typename T >
bool operator!= (const EulerAngles< T > &l, const EulerAngles< T > &r)
 
template<typename T >
Tensor2< T > EulerZXZRotationMatrix (T alpha, T beta, T gamma)
 
template<typename T >
Tensor2< T > EulerZXZRotationMatrix (EulerAngles< T > angle)
 
template<typename T >
EulerAngles< T > getEulerZXZAngles (Tensor2< T > R)
 Get Euler angles from rotation matrix. More...
 
std::mt19937 makeMt19937 (bool defaultSeed)
 
template<typename T >
void rotationTensorFrom3UniformRandomsArvo1992 (Tensor2< T > &A, T x1, T x2, T x3)
 
template<typename T >
void randomRotationTensorArvo1992 (std::mt19937 &gen, Tensor2< T > &A)
 
template<typename T >
void randomRotationTensorArvo1992 (Tensor2< T > &A, bool defaultSeed=false)
 
template<typename T >
void rotationTensorFrom3UniformRandoms (Tensor2< T > &A, T x1, T x2, T x3)
 
template<typename T >
Tensor2< T > randomRotationTensorStewart1980 (unsigned int dim, bool defaultSeed=false)
 Generate a random rotation tensor. More...
 
template<typename T >
void randomRotationTensorInPlace (unsigned int dim, Tensor2< T > &A, bool defaultSeed=false)
 Generate a random rotation tensor. More...
 
template<typename T >
Tensor2< T > randomRotationTensor (unsigned int dim, bool defaultSeed=false)
 
template<typename T >
EulerAngles< T > quaternionToEulerAngles (isoi::Quaternion &q)
 
int fftwFlat2 (int i0, int i1, int N0, int N1)
 
int fftwFlat4 (int i0, int i1, int i2, int i3, const std::vector< ptrdiff_t > &N)
 
double fftwComplexMag (fftw_complex comps)
 
std::string getCoordsName (std::string baseName)
 
std::string getCoeffsName (std::string baseName)
 
template<typename U >
std::string getComponentSuffix (const std::vector< U > &componentIdx)
 
std::string getDefaultCoeffsDatasetName (SpectralDatasetID dsetID)
 
template<typename U >
void loadSpectralDatabase (hpp::HDF5Handler &infile, std::string dsetNameCoords, std::string dsetNameCoeffs, std::vector< hsize_t > componentIdx, std::vector< std::vector< unsigned int >> &coordsList, std::vector< std::complex< U >> &coeffList, unsigned int nTerms, unsigned int refineMult=1)
 Load a spectral database. More...
 
template<typename U >
std::complex< U > TruncatedIFFTNDSingleSquare (const unsigned int spatialDim, const std::vector< std::complex< U >> &expTable, const std::vector< std::vector< unsigned int >> &spectralCoordsList, const std::vector< std::complex< U >> &coeffs, std::vector< unsigned int > spatialCoord, unsigned int nTerms)
 TruncatedIFFTND where all dimensions have the same number of points. More...
 
template<typename U >
TruncatedIFFTNDSingleSquareReal (const unsigned int spatialDim, const std::vector< std::complex< U >> &expTable, const std::vector< std::vector< unsigned int >> &spectralCoordsList, const std::vector< std::complex< U >> &coeffs, std::vector< unsigned int > spatialCoord, unsigned int nTerms)
 
template<typename U >
TruncatedIFFTNDSingleSquareRealScalar (const unsigned int nDims, const unsigned int spatialDim, const U expTable[], const int spectralCoordsList[], const U coeffs[], std::vector< unsigned int > spatialCoord, unsigned int nTerms)
 Scalar implementation of hpp.TruncatedIFFTNDSingleSquareReal. More...
 
template<typename U >
TruncatedIFFTNDSingleSquareRealScalar (const unsigned int nDims, const unsigned int spatialDim, const int spectralCoordsList[], const U coeffs[], std::vector< unsigned int > spatialCoord, unsigned int nTerms)
 
template<typename U >
TruncatedIFFT4DSingleSquareReal (const unsigned int spatialDim, const U expTable[], const int spectralCoordsList[], const U coeffs[], const std::vector< unsigned int > &spatialCoord, unsigned int nTerms)
 AVX implementation of hpp.TruncatedIFFT4DSingleSquareReal. More...
 
template<typename U >
TruncatedIFFTNDSingleSquareReal (const unsigned int nDims, const unsigned int spatialDim, const U expTable[], const int spectralCoordsList[], const U coeffs[], const std::vector< unsigned int > &spatialCoord, unsigned int nTerms)
 
template<typename U >
TruncatedIFFT4DSingleSquareReal (const unsigned int spatialDim, const int spectralCoordsList[], const U coeffs[], const std::vector< unsigned int > &spatialCoord, unsigned int nTerms)
 
template<typename U >
TruncatedIFFTNDSingleSquareReal (const unsigned int nDims, const unsigned int spatialDim, const int spectralCoordsList[], const U coeffs[], const std::vector< unsigned int > &spatialCoord, unsigned int nTerms)
 
template<typename T , unsigned int N>
__global__ void GET_IDFT_REAL (SpectralDatabaseCUDA< T, N > *db, unsigned int dsetIdx, unsigned int *spatialCoord, T *val)
 
template<typename T , unsigned int N, unsigned int P>
struct ALIGN (16) SpectralDataUnifiedCUDA
 
template<typename T , unsigned int N>
__device__ void getExpVal (unsigned int *spatialCoord, SpectralCoordCUDA< N > &coord, unsigned int gridDim, T expArgFactor, T *expValRe, T *expValIm)
 
template<typename T >
std::vector< T > operator/ (const std::vector< T > &vec, T scalar)
 Elementwise division. More...
 
template<typename T >
std::vector< std::vector< T > > operator/ (const std::vector< std::vector< T >> &veclist, T scalar)
 
template<typename T >
void operator*= (std::vector< T > &vec, const T scalar)
 
template<typename T >
void operator/= (std::vector< T > &vec, const T scalar)
 
template<typename T >
std::vector< T > operator* (const std::vector< T > &vec, const T scalar)
 
template<typename T >
std::vector< T > operator* (const T scalar, const std::vector< T > &vec)
 
template<typename T >
std::vector< std::vector< T > > operator* (T scalar, const std::vector< std::vector< T >> &veclist)
 
template<typename T >
std::vector< T > abs (const std::vector< T > &vec)
 
template<typename T >
min (const std::vector< T > &vec)
 
template<typename T >
max (const std::vector< T > &vec)
 
template<typename T >
std::vector< T > operator- (const std::vector< T > &vec1, const std::vector< T > &vec2)
 
template<typename T >
prod (const std::vector< T > &vec)
 Product of all the entries of a vector. More...
 
template<typename T >
std::ostream & operator<< (std::ostream &out, const std::vector< T > &vec)
 
template<typename T >
std::ostream & operator<< (std::ostream &out, const std::valarray< T > &vec)
 
idx2d unflat (unsigned int flatIdx, unsigned int n1, unsigned int n2)
 
idx4d unflat (unsigned int flatIdx, unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
 
template<typename T , typename U >
std::vector< T > unflatC (T flatIdx, std::vector< U > dims)
 
template<typename T , typename U >
flatC (std::vector< T > idx, std::vector< U > dims)
 
template<typename T >
std::vector< T > ones (unsigned int n)
 
template<typename T >
Tensor2< T > outer (const std::vector< T > &A, const std::vector< T > &B)
 
template<typename T >
void outerInPlace (const std::vector< T > &A, const std::vector< T > &B, Tensor2< T > &C)
 
template<typename T >
Tensor2< T > outer (const std::valarray< T > &A, const std::valarray< T > &B)
 
template<typename T >
Tensor2< T > identityTensor2 (unsigned int n)
 
template<typename T >
void identityTensor2InPlace (unsigned int n, Tensor2< T > &A)
 
template<typename T >
Tensor2< T > ones2 (unsigned int n)
 
template<typename T >
Tensor2< T > diag2 (const std::vector< T > &vals)
 
template<typename T >
Tensor2< T > diag2 (const std::valarray< T > &vals)
 
template<typename T >
bool areEqual (const Tensor2< T > &A, const Tensor2< T > &B, T tol)
 
template<typename T >
bool operator== (const Tensor2< T > &A, const Tensor2< T > &B)
 Test equality of two tensors. More...
 
template<typename T >
bool operator!= (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
bool areSameShape (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
void assertSameShape (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
Tensor2< T > operator+ (const Tensor2< T > &A, const Tensor2< T > &B)
 Addition. More...
 
template<typename T >
void operator+= (Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
Tensor2< T > operator+ (const Tensor2< T > &A, const T &B)
 
template<typename T >
Tensor2< T > operator+ (const T &A, const Tensor2< T > &B)
 
template<typename T >
void operator+= (Tensor2< T > &A, const T &B)
 
template<typename T >
Tensor2< T > MPISum (Tensor2< T > &local, MPI_Comm comm)
 
template<typename T >
Tensor2< T > operator- (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
void operator-= (Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
Tensor2< T > operator- (const Tensor2< T > &A, const T &B)
 
template<typename T >
Tensor2< T > operator- (const T &A, const Tensor2< T > &B)
 
template<typename T >
void operator-= (Tensor2< T > &A, const T &B)
 
template<typename T >
Tensor2< T > operator* (const Tensor2< T > &A, const T &B)
 
template<typename T >
Tensor2< T > operator* (const T &A, const Tensor2< T > &B)
 
template<typename T >
void operator*= (Tensor2< T > &A, const T &B)
 
template<typename T >
Tensor2< T > operator/ (const Tensor2< T > &A, const T &B)
 
template<typename T >
void operator/= (Tensor2< T > &A, const T &B)
 
template<typename T >
void assertCompatibleForMultiplication (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
Tensor2< T > operator* (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
void ABPlusBTransposeAInPlace (const hpp::Tensor2< T > &A, const hpp::Tensor2< T > &B, hpp::Tensor2< T > &C)
 
template<typename T >
contract (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
std::ostream & operator<< (std::ostream &out, const Tensor2< T > &A)
 
template<typename T >
std::ostream & operator<< (std::ostream &out, Tensor4< T > &A)
 
template<typename T >
Tensor4< T > identityTensor4 (unsigned int n)
 
template<typename T >
void identityTensor4InPlace (unsigned int n, Tensor4< T > &A)
 
template<typename T >
bool areSameShape (const Tensor4< T > &A, const Tensor4< T > &B)
 
template<typename T >
void assertSameShape (const Tensor4< T > &A, const Tensor4< T > &B)
 
template<typename T >
bool operator== (const Tensor4< T > &A, const Tensor4< T > &B)
 Test equality of two tensors. More...
 
template<typename T >
bool operator!= (const Tensor4< T > &A, const Tensor4< T > &B)
 
template<typename T >
Tensor4< T > operator+ (const Tensor4< T > &A, const Tensor4< T > &B)
 Addition. More...
 
template<typename T >
void operator+= (Tensor4< T > &A, const Tensor4< T > &B)
 
template<typename T >
Tensor4< T > operator+ (const Tensor4< T > &A, const T &B)
 
template<typename T >
Tensor4< T > operator+ (const T &A, const Tensor4< T > &B)
 
template<typename T >
void operator+= (Tensor4< T > &A, const T &B)
 
template<typename T >
Tensor4< T > operator- (const Tensor4< T > &A)
 
template<typename T >
Tensor4< T > operator- (const Tensor4< T > &A, const Tensor4< T > &B)
 
template<typename T >
void operator-= (Tensor4< T > &A, const Tensor4< T > &B)
 
template<typename T >
Tensor4< T > operator- (const Tensor4< T > &A, const T &B)
 
template<typename T >
Tensor4< T > operator- (const T &A, const Tensor4< T > &B)
 
template<typename T >
void operator-= (Tensor4< T > &A, const T &B)
 
template<typename T >
Tensor4< T > operator* (const Tensor4< T > &A, const T &B)
 
template<typename T >
Tensor4< T > operator* (const T &A, const Tensor4< T > &B)
 
template<typename T >
void operator*= (Tensor4< T > &A, const T &B)
 
template<typename T >
Tensor4< T > operator/ (const Tensor4< T > &A, const T &B)
 
template<typename T >
void operator/= (Tensor4< T > &A, const T &B)
 
template<typename T >
Tensor4< T > outer (const Tensor2< T > &A, const Tensor2< T > &B)
 
template<typename T >
void assertCompatibleForOuterProduct (const Tensor2< T > &A, const Tensor2< T > &B, const Tensor4< T > &C)
 
template<typename T >
void outerInPlace (const Tensor2< T > &A, const Tensor2< T > &B, Tensor4< T > &C)
 
template<typename T >
void assertCompatibleForContraction (const Tensor4< T > &A, const Tensor2< T > &B)
 
template<typename T >
void assertCompatibleForContraction (const Tensor4< T > &A, const Tensor2< T > &B, const Tensor2< T > &C)
 
template<typename T >
Tensor2< T > contract (const Tensor4< T > &A, const Tensor2< T > &B)
 
template<typename T >
void contractInPlace (const Tensor4< T > &A, const Tensor2< T > &B, Tensor2< T > &C)
 
template<typename U >
void assertCompatibleForContraction (const Tensor4< U > &A, const Tensor4< U > &B)
 
template<typename U >
Tensor4< U > contract (const Tensor4< U > &A, const Tensor4< U > &B)
 
template<typename T >
std::vector< T > operator* (const hpp::Tensor2< T > &A, const std::vector< T > &x)
 
template<typename T >
std::valarray< T > operator* (const hpp::Tensor2< T > &A, const std::valarray< T > &x)
 
template<typename T >
hpp::Tensor2< T > operator* (const std::vector< T > &x, const hpp::Tensor2< T > &A)
 Replicates the numpy behaviour of vec*mat. More...
 
template<typename T >
Tensor2< T > transformIntoFrame (const Tensor2< T > &A, const Tensor2< T > &Q)
 Transform tensor \( \mathbf{A} \) into the frame given by the columns of \( \mathbf{Q} \). More...
 
template<typename T >
Tensor2< T > transformOutOfFrame (const Tensor2< T > &A_star, const Tensor2< T > &Q)
 Transform tensor \( \mathbf{A}^* \) out of the frame given by the columns of \( \mathbf{Q} \). More...
 
template<typename T >
Tensor4< T > transformOutOfFrame (const Tensor4< T > &E_star, const Tensor2< T > &Q)
 Transform tensor \( \mathbf{E}^* \) out of the frame given by the columns of \( \mathbf{Q} \). More...
 
template<typename T , unsigned int N>
__device__ VecCUDA< T, N > operator/ (const VecCUDA< T, N > &inVec, T scalar)
 
template<typename T >
__device__ VecCUDA< T, 3 > cartesianToSpherical (const VecCUDA< T, 3 > &cartVec)
 Uses the "mathematics" convention. More...
 
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA< T, M, N > operator- (const Tensor2CUDA< T, M, N > &A, const Tensor2CUDA< T, M, N > &B)
 
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA< T, M, N > operator+ (const Tensor2CUDA< T, M, N > &A, const Tensor2CUDA< T, M, N > &B)
 
template<typename T , unsigned int M, unsigned N>
__host__ __device__ void operator+= (Tensor2CUDA< T, M, N > &A, const Tensor2CUDA< T, M, N > &B)
 
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA< T, M, N > operator* (const Tensor2CUDA< T, M, N > &A, T scalar)
 
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA< T, M, N > operator* (T scalar, const Tensor2CUDA< T, M, N > &A)
 
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA< T, M, N > operator/ (const Tensor2CUDA< T, M, N > &A, T scalar)
 
template<typename T , unsigned int M, unsigned N>
__host__ __device__ void operator/= (Tensor2CUDA< T, M, N > &A, T scalar)
 
template<typename T , unsigned int M, unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA< T, M, P > operator* (const Tensor2CUDA< T, M, N > &A, const Tensor2CUDA< T, N, P > &B)
 
template<typename T , unsigned int M, unsigned int N>
__host__ __device__ VecCUDA< T, M > operator* (const Tensor2CUDA< T, M, N > &A, const VecCUDA< T, N > &x)
 
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA< T, M, M > transformIntoFrame (const Tensor2CUDA< T, M, M > &A, const Tensor2CUDA< T, M, M > &Q)
 Transform tensor \( \mathbf{A} \) into the frame given by the columns of \( \mathbf{Q} \). More...
 
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA< T, M, M > transformOutOfFrame (const Tensor2CUDA< T, M, M > &A_star, const Tensor2CUDA< T, M, M > &Q)
 Transform tensor \( \mathbf{A}^* \) out of the frame given by the columns of \( \mathbf{Q} \). More...
 
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA< T, M, M > transformIntoFrame (const Tensor2AsymmCUDA< T, M > &A, const Tensor2CUDA< T, M, M > &Q)
 
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA< T, M, M > transformOutOfFrame (const Tensor2AsymmCUDA< T, M > &A_star, const Tensor2CUDA< T, M, M > &Q)
 
template<typename T >
__device__ Tensor2CUDA< T, 3, 3 > EulerZXZRotationMatrixCUDA (T alpha, T beta, T gamma)
 
template<typename T , unsigned int M, unsigned int N>
std::ostream & operator<< (std::ostream &out, const Tensor2CUDA< T, M, N > &A)
 
template<typename T , unsigned int M, unsigned N>
__device__ Tensor2CUDA< T, M, N > warpReduceSumTensor2 (Tensor2CUDA< T, M, N > A)
 
template<typename T , unsigned int M, unsigned N>
__device__ Tensor2CUDA< T, M, N > blockReduceSumTensor2 (Tensor2CUDA< T, M, N > val)
 
template<typename T , unsigned int M, unsigned N>
__global__ void BLOCK_REDUCE_KEPLER_TENSOR2 (Tensor2CUDA< T, M, N > *in, Tensor2CUDA< T, M, N > *out, int nTerms)
 
template<typename T , unsigned int M, unsigned int N>
__device__ EulerAngles< T > getEulerZXZAngles (const Tensor2CUDA< T, M, N > &R)
 Coresponds to the ZXZ Proper Euler Angles. More...
 
template<typename T , unsigned int N>
__host__ __device__ Tensor2SymmCUDA< T, N > operator- (const Tensor2SymmCUDA< T, N > &A, const Tensor2SymmCUDA< T, N > &B)
 
template<typename T , unsigned int N>
__host__ __device__ Tensor2SymmCUDA< T, N > operator+ (const Tensor2SymmCUDA< T, N > &A, const Tensor2SymmCUDA< T, N > &B)
 
template<typename T , unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA< T, N, P > operator* (const Tensor2SymmCUDA< T, N > &A, const Tensor2CUDA< T, N, P > &B)
 
template<typename T , unsigned int N>
__host__ __device__ Tensor2AsymmCUDA< T, N > operator- (const Tensor2AsymmCUDA< T, N > &A, const Tensor2AsymmCUDA< T, N > &B)
 
template<typename T , unsigned int N>
__host__ __device__ Tensor2CUDA< T, N, N > operator- (const Tensor2AsymmCUDA< T, N > &A, const Tensor2CUDA< T, N, N > &B)
 
template<typename T , unsigned int N>
__host__ __device__ Tensor2AsymmCUDA< T, N > operator+ (const Tensor2AsymmCUDA< T, N > &A, const Tensor2AsymmCUDA< T, N > &B)
 
template<typename T , unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA< T, N, P > operator* (const Tensor2AsymmCUDA< T, N > &A, const Tensor2CUDA< T, N, P > &B)
 
template<typename T , unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA< T, N, P > operator* (const Tensor2CUDA< T, N, P > &A, const Tensor2AsymmCUDA< T, P > &B)
 

Variables

H5E_auto2_t dummyHDFErrorHandler
 
void * dummyHDFClientData
 
hid_t hdf_complex_id = getComplexHDFType()
 

Enumeration Type Documentation

Enumerator
SIGMA00 
SIGMA11 
SIGMA12 
SIGMA02 
SIGMA01 
WP01 
WP12 
WP02 
GAMMA 
Enumerator
CRYSTAL_TYPE_NONE 
CRYSTAL_TYPE_FCC 
Enumerator
HARDENING_LAW_BROWN 
HARDENING_LAW_VOCE 
enum hpp::HDFReadWrite : char
strong
Enumerator
Write 
Read 
Enumerator
SYMMETRY_TYPE_NONE 
SYMMETRY_TYPE_C4 

Function Documentation

template<typename T >
void hpp::ABPlusBTransposeAInPlace ( const hpp::Tensor2< T > &  A,
const hpp::Tensor2< T > &  B,
hpp::Tensor2< T > &  C 
)
template<typename T >
std::vector<T> hpp::abs ( const std::vector< T > &  vec)
template<typename T >
void hpp::addAttribute ( H5::H5File &  file,
const std::string &  attrName,
attrVal 
)
template<typename T , unsigned int N, unsigned int P>
struct hpp::ALIGN ( 16  )
template<typename T >
bool hpp::areEqual ( const Tensor2< T > &  A,
const Tensor2< T > &  B,
tol 
)
template<typename T >
bool hpp::areSameShape ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
inline
template<typename T >
bool hpp::areSameShape ( const Tensor4< T > &  A,
const Tensor4< T > &  B 
)
inline
template<typename T >
void hpp::assertCompatibleForContraction ( const Tensor4< T > &  A,
const Tensor2< T > &  B 
)
inline
template<typename T >
void hpp::assertCompatibleForContraction ( const Tensor4< T > &  A,
const Tensor2< T > &  B,
const Tensor2< T > &  C 
)
inline
template<typename U >
void hpp::assertCompatibleForContraction ( const Tensor4< U > &  A,
const Tensor4< U > &  B 
)
inline
template<typename T >
void hpp::assertCompatibleForMultiplication ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T >
void hpp::assertCompatibleForOuterProduct ( const Tensor2< T > &  A,
const Tensor2< T > &  B,
const Tensor4< T > &  C 
)
template<typename T >
void hpp::assertSameShape ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T >
void hpp::assertSameShape ( const Tensor4< T > &  A,
const Tensor4< T > &  B 
)
template<typename T >
__global__ void hpp::BLOCK_REDUCE_KEPLER_GSH_COEFFS ( GSHCoeffsCUDA< T > *  in,
GSHCoeffsCUDA< T > *  out,
int  nTerms 
)
template<typename T , unsigned int M, unsigned N>
__global__ void hpp::BLOCK_REDUCE_KEPLER_TENSOR2 ( Tensor2CUDA< T, M, N > *  in,
Tensor2CUDA< T, M, N > *  out,
int  nTerms 
)
template<typename T >
__device__ GSHCoeffsCUDA<T> hpp::blockReduceSumGSHCoeffs ( GSHCoeffsCUDA< T >  val)
inline
template<typename T , unsigned int M, unsigned N>
__device__ Tensor2CUDA<T,M,N> hpp::blockReduceSumTensor2 ( Tensor2CUDA< T, M, N >  val)
inline
template<typename T >
__device__ VecCUDA<T,3> hpp::cartesianToSpherical ( const VecCUDA< T, 3 > &  cartVec)

Uses the "mathematics" convention.

Parameters
cartVecthe vector in cartesian coordinates
Returns
the vector in spherical coordinates
template<typename T >
std::vector<T> hpp::cartesianToSpherical ( const std::vector< T > &  cartVec)
template<typename T >
T hpp::contract ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T >
Tensor2<T> hpp::contract ( const Tensor4< T > &  A,
const Tensor2< T > &  B 
)
inline
template<typename U >
Tensor4<U> hpp::contract ( const Tensor4< U > &  A,
const Tensor4< U > &  B 
)
inline
template<typename T >
void hpp::contractInPlace ( const Tensor4< T > &  A,
const Tensor2< T > &  B,
Tensor2< T > &  C 
)
inline
template<typename T >
H5::DataSet hpp::createHDF5Dataset ( H5::H5File &  file,
const H5std_string &  datasetName,
std::vector< hsize_t >  dataDims 
)
template<typename T >
hid_t hpp::createHDF5Dataset ( hid_t  file_id,
std::string  datasetName,
std::vector< hsize_t >  dataDims 
)
template<typename T >
H5::DataSet hpp::createHDF5GridOfArrays ( H5::H5File &  file,
const H5std_string &  datasetName,
std::vector< hsize_t >  gridDims,
std::vector< hsize_t >  arrayDims 
)
template<typename T >
hid_t hpp::createHDF5GridOfArrays ( hid_t  file_id,
std::string  datasetName,
std::vector< hsize_t >  gridDims,
std::vector< hsize_t >  arrayDims 
)
template<typename U >
Tensor4<U> hpp::cubeSymmetricElasticityTensor ( c11,
c12,
c44 
)

Construct a fourth order elasticity tensor accounting for cubic symmetry.

Parameters
c11elastic constant \({C}_11\)
c12elastic constant \({C}_12\)
c44elastic constant \({C}_44\)
Returns
The elasticity tensor
template<typename U >
CrystalSolverConfig<U> hpp::defaultConservativeCrystalSolverConfig ( )
template<typename U >
CrystalInitialConditions<U> hpp::defaultCrystalInitialConditions ( )
template<typename U >
CrystalProperties<U> hpp::defaultCrystalProperties ( )
template<typename U >
CrystalSolverConfig<U> hpp::defaultCrystalSolverConfig ( )
std::vector< SpectralDatasetID > hpp::defaultCrystalSpectralDatasetIDs ( )
template<typename U >
Tensor2<U> hpp::deformationGradFromVelGrad ( t_init,
const Tensor2< U > &  F_init,
const Tensor2< U > &  L,
t 
)
void hpp::destroyConfigRealND ( FFTWConfigRealND cfg)
template<typename T >
Tensor2<T> hpp::diag2 ( const std::vector< T > &  vals)
template<typename T >
Tensor2<T> hpp::diag2 ( const std::valarray< T > &  vals)
template<typename T >
Tensor2<T> hpp::EulerZXZRotationMatrix ( alpha,
beta,
gamma 
)
template<typename T >
Tensor2<T> hpp::EulerZXZRotationMatrix ( EulerAngles< T >  angle)
template<typename T >
__device__ Tensor2CUDA<T,3,3> hpp::EulerZXZRotationMatrixCUDA ( EulerAngles< T >  angles)
template<typename T >
__device__ Tensor2CUDA<T,3,3> hpp::EulerZXZRotationMatrixCUDA ( alpha,
beta,
gamma 
)
void hpp::evaluateSpectralCompressionErrorAxisSlice ( std::string  rawDBName,
std::string  spectralDBName,
std::vector< std::string >  dsetBasenames,
unsigned int  nTermsMax,
unsigned int  refineMult,
std::vector< int >  axisSlice,
std::string  outFilename,
MPI_Comm  comm 
)

Evaluate the error incurred in the spectral reresentation of a database.

Parameters
rawDBName
spectralDBName
dsetBasenames
nTermsMax
refineMult
axisSlicea vector of length N, where N is the dimension of the dataset, containing [axis, otherCoord0, otherCoord1, ..., otherCoordN-1], where axis specifies the axis along which a slice is to be taken, and the values otherCoord0, ..., otherCoordN-1 specify in order what the fixed values of the coordinates should be in the other dimensions. For example, an axis slice of [2, 19, 3, 5] specifies that the slice is to be taken along axis/dimension "2", that is, the third axis. The following numbers then dictate that the fixed coordinates for the other axes should be:
  • 19 for axis 0
  • 3 for axis 1
  • 5 for axis 3
outFilename
comm
void hpp::evaluateSpectralCompressionErrorAxisSliceUnified ( std::string  rawDBName,
std::string  spectralDBName,
std::vector< std::string >  dsetBasenames,
unsigned int  nTermsMax,
unsigned int  refineMult,
std::vector< int >  axisSlice,
std::string  outFilename,
MPI_Comm  comm 
)

Evaluate the error incurred in the spectral reresentation of a database.

Parameters
rawDBName
spectralDBName
dsetBasenames
nTermsMax
refineMult
axisSlicea vector of length N, where N is the dimension of the dataset, containing [axis, otherCoord0, otherCoord1, ..., otherCoordN-1], where axis specifies the axis along which a slice is to be taken, and the values otherCoord0, ..., otherCoordN-1 specify in order what the fixed values of the coordinates should be in the other dimensions. For example, an axis slice of [2, 19, 3, 5] specifies that the slice is to be taken along axis/dimension "2", that is, the third axis. The following numbers then dictate that the fixed coordinates for the other axes should be:
  • 19 for axis 0
  • 3 for axis 1
  • 5 for axis 3
outFilename
comm
void hpp::evaluateSpectralCompressionErrorFull ( std::string  rawDBName,
std::string  spectralDBName,
std::string  errorDBName,
unsigned int  nTermsMax,
std::string  outFilename,
MPI_Comm  comm 
)
void hpp::evaluateSpectralCompressionErrorFullUnified ( std::string  rawDBName,
std::string  spectralDBName,
std::string  errorDBName,
unsigned int  nTermsMax,
std::string  outFilename,
MPI_Comm  comm 
)
double hpp::fftwComplexMag ( fftw_complex  comps)
inline
int hpp::fftwFlat2 ( int  i0,
int  i1,
int  N0,
int  N1 
)
inline

FFTW uses row-major flattening

int hpp::fftwFlat4 ( int  i0,
int  i1,
int  i2,
int  i3,
const std::vector< ptrdiff_t > &  N 
)
inline
template<typename T , typename U >
T hpp::flatC ( std::vector< T >  idx,
std::vector< U >  dims 
)
template<typename T >
__global__ void hpp::GET_AVERAGE_TCAUCHY ( unsigned int  nCrystals,
const SpectralCrystalCUDA< T > *  crystals,
Tensor2CUDA< T, 3, 3 > *  TCauchyGlobal 
)
template<typename T >
__global__ void hpp::GET_GSH_COEFFS ( const SpectralCrystalCUDA< T > *  crystals,
unsigned int  ncrystals,
GSHCoeffsCUDA< T > *  coeffsPerBlockSums 
)
template<typename T , unsigned int N>
__global__ void hpp::GET_IDFT_REAL ( SpectralDatabaseCUDA< T, N > *  db,
unsigned int  dsetIdx,
unsigned int *  spatialCoord,
T *  val 
)
std::string hpp::getCoeffsName ( std::string  baseName)
inline
hid_t hpp::getComplexHDFType ( )
template<typename U >
std::string hpp::getComponentSuffix ( const std::vector< U > &  componentIdx)
std::string hpp::getCoordsName ( std::string  baseName)
inline
std::vector< hsize_t > hpp::getDatasetDims ( hid_t  dset_id)
std::string hpp::getDefaultCoeffsDatasetName ( SpectralDatasetID  dsetID)
inline
template<typename U >
Tensor4<U> hpp::getEshelbyTensorCubicMaterialSphericalInclusion ( c11,
c12,
c44,
I0,
I1,
I2 
)
template<typename T >
EulerAngles<T> hpp::getEulerAnglesFromDeformationGradient ( const hpp::Tensor2< T > &  F_star)

Get the corresponding Euler angles for the rotation induced by this deformation.

Parameters
F_star
Returns
template<typename T >
MPI_Datatype hpp::getEulerAnglesTypeMPI ( )
template<typename T >
EulerAngles<T> hpp::getEulerZXZAngles ( Tensor2< T >  R)

Get Euler angles from rotation matrix.

Parameters
Rthe rotation matrix
Returns
the angles
template<typename T , unsigned int M, unsigned int N>
__device__ EulerAngles<T> hpp::getEulerZXZAngles ( const Tensor2CUDA< T, M, N > &  R)

Coresponds to the ZXZ Proper Euler Angles.

Parameters
Rthe rotation matrix
Returns
the angles
template<typename T , unsigned int N>
__device__ void hpp::getExpVal ( unsigned int *  spatialCoord,
SpectralCoordCUDA< N > &  coord,
unsigned int  gridDim,
expArgFactor,
T *  expValRe,
T *  expValIm 
)
template<typename T >
H5::DataType hpp::getHDF5Type ( )

Get HDF5 equivalent type.

For C++ interface

Returns
The HDF5 type
template<typename T >
hid_t hpp::getHDF5TypeC ( )

Get HDF5 equivalent type.

For the C interface

Returns
The HDF5 type
template<typename U >
Tensor4<U> hpp::getHomogenizedStiffnessELSC ( const std::vector< Tensor4< U >> &  cVec,
const std::vector< U > &  vVec,
const Tensor4< U > &  S 
)

Get the homogenized stiffness tensor using the elastic self-consistent method.

Parameters
cVeca vector of stiffness tensors to homogenize
vVeca vector of volumes occupied by the stiffness tensors
Sthe Eshelby tensor
Returns
the homogenized stiffness tensor
template<typename U >
Tensor4<U> hpp::getHomogenizedStiffnessVolumeAverage ( const std::vector< Tensor4< U >> &  cVec,
const std::vector< U > &  vVec 
)

Get the homogenized stiffness tensor using a volume average.

Parameters
cVeca vector of stiffness tensors to homogenize
vVeca vector of volumes occupied by the stiffness tensors
Returns
the homogenized stiffness tensor
template<typename T >
HDFReadWriteParamsC hpp::getReadWriteParametersForMultipleHDF5Arrays ( hid_t  dset_id,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims,
std::vector< hsize_t >  arrayCount 
)

Gets the parameters and ensures that they're compatible with the actual dataset.

Parameters
dset_id
gridOffset
arrayDims
arrayCount
Returns
template<typename T >
HDFReadWriteParams hpp::getReadWriteParametersForSingleHDF5Array ( H5::DataSet &  dataset,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims 
)
template<typename T >
HDFReadWriteParamsC hpp::getReadWriteParametersForSingleHDF5Array ( hid_t  dset_id,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims 
)

Gets the parameters and ensures that they're compatible with the actual dataset.

Parameters
dset_id
gridOffset
arrayDims
Returns
template<typename T , unsigned int P>
__device__ void hpp::getSpectralCrystalDatabaseCoordinate ( SpectralCrystalCUDA< T > &  crystal,
SpectralDatabaseUnifiedCUDA< T, 4, P > *  db,
Tensor2CUDA< T, 3, 3 > &  RStretchingTensor,
theta,
unsigned int *  dbCoord 
)
template<typename T >
StretchingTensorDecomposition<T> hpp::getStretchingTensorDecomposition ( const Tensor2< T > &  L)
Parameters
Lthe velocity gradient tensor
Returns
void hpp::hdfAssert ( herr_t  code,
const char *  file,
int  line 
)
inline
void hpp::hdfRestoreErrors ( )
void hpp::hdfSuppressErrors ( )
template<typename T >
Tensor2<T> hpp::histogramPoleEqualArea ( const std::vector< EulerAngles< T >> &  angles,
const std::vector< T > &  planeNormal 
)
template<typename T >
Tensor2<T> hpp::identityTensor2 ( unsigned int  n)
template<typename T >
void hpp::identityTensor2InPlace ( unsigned int  n,
Tensor2< T > &  A 
)
template<typename T >
Tensor4<T> hpp::identityTensor4 ( unsigned int  n)
inline
template<typename T >
void hpp::identityTensor4InPlace ( unsigned int  n,
Tensor4< T > &  A 
)
template<typename U >
Tensor4<U> hpp::isotropicElasticityTensor ( mu,
kappa 
)

Construct a fourth order isotropic elasticity tensor.

Parameters
muthe elastic shear modulus \(\mu\)
kappathe elastic bulk modulus \(\kappa\)
Returns
The elasticity tensor
template<typename U >
void hpp::loadSpectralDatabase ( hpp::HDF5Handler infile,
std::string  dsetNameCoords,
std::string  dsetNameCoeffs,
std::vector< hsize_t >  componentIdx,
std::vector< std::vector< unsigned int >> &  coordsList,
std::vector< std::complex< U >> &  coeffList,
unsigned int  nTerms,
unsigned int  refineMult = 1 
)

Load a spectral database.

Parameters
infilethe file
dsetNameCoordsthe name of the dataset containing the spectral coordinates
dsetNameCoeffsthe name of the dataset containing the corresponding spectral coefficients
componentIdxthe vector or tensor index of the component of the quantity to load
coordsListthe list of coordinates to return
coeffListthe list of coefficients to return
nTermsthe number of Fourier terms to load
refineMultthe spectral refinement multiplier to use
std::mt19937 hpp::makeMt19937 ( bool  defaultSeed)
inline
template<typename T >
T hpp::max ( const std::vector< T > &  vec)
template<typename T >
T hpp::min ( const std::vector< T > &  vec)
bool hpp::MPIAllTrue ( bool  condition,
MPI_Comm  comm 
)

Determine if condition is true for all processes.

Parameters
condition
comm
Returns
True if true for all, false otherwise
template<typename T >
T hpp::MPIBroadcastFromRoot ( rootVal,
MPI_Comm  comm 
)
template<typename T >
std::vector<T> hpp::MPIBroadcastFromRoot ( std::vector< T >  rootVec,
MPI_Comm  comm 
)
template<typename T >
std::vector<T> hpp::MPIConcatOnRoot ( localVal,
MPI_Comm  comm 
)
template<typename T >
std::vector<T> hpp::MPIConcatOnRoot ( std::vector< T >  localVec,
MPI_Comm  comm 
)
template<typename T >
T hpp::MPIMax ( const T &  localVal,
MPI_Comm  comm 
)
template<typename T >
T hpp::MPIMin ( const T &  localVal,
MPI_Comm  comm 
)
template<typename T >
std::vector<T> hpp::MPISplitVectorEvenly ( const std::vector< T > &  rootVec,
MPI_Comm  comm 
)
template<typename T >
std::vector<T> hpp::MPISplitVectorEvenly ( const std::vector< T > &  rootVec,
MPI_Comm  comm,
MPI_Datatype  dtype 
)
template<typename T >
T hpp::MPISum ( const T &  localVal,
MPI_Comm  comm 
)
template<typename T >
Tensor2<T> hpp::MPISum ( Tensor2< T > &  local,
MPI_Comm  comm 
)
template<typename T >
MPI_Datatype hpp::MPIType ( )

Get the MPI datatype from the C++ type.

Returns
The MPI_Datatype
template<typename U >
hpp::Tensor2<U> hpp::newtonStressCorrection ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  A,
const std::vector< hpp::Tensor2< U >> &  C_alphas,
const std::vector< U > &  s_alphas,
const hpp::Tensor2< U > &  T_old,
const U  DT_max,
const U  dt 
)

See kalidindi1992.NewtonStressCorrection.

template<typename U >
void hpp::newtonStressCorrectionInPlace ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  A,
const std::vector< hpp::Tensor2< U >> &  C_alphas,
const std::vector< U > &  s_alphas,
const hpp::Tensor2< U > &  T_old,
const U  DT_max,
const U  dt,
std::vector< U > &  dumDgamma_alphas,
hpp::Tensor4< U > &  dumJ,
hpp::Tensor2< U > &  dumT_tr,
hpp::Tensor2< U > &  dumG,
hpp::Tensor2< U > &  dum2ndOrder,
hpp::Tensor4< U > &  dum4thOrder,
hpp::Tensor2< U > &  DT 
)

In-place version of hpp.newtonStressCorrection.

constexpr int hpp::nSlipSystems ( CrystalType  crystalType)
template<typename T >
std::vector<T> hpp::ones ( unsigned int  n)
template<typename T >
Tensor2<T> hpp::ones2 ( unsigned int  n)
template<typename T >
bool hpp::operator!= ( const EulerAngles< T > &  l,
const EulerAngles< T > &  r 
)
template<typename T >
bool hpp::operator!= ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T >
bool hpp::operator!= ( const Tensor4< T > &  A,
const Tensor4< T > &  B 
)
template<typename T >
std::vector<T> hpp::operator* ( const std::vector< T > &  vec,
const T  scalar 
)
template<typename T >
std::vector<T> hpp::operator* ( const T  scalar,
const std::vector< T > &  vec 
)
template<typename T >
std::vector<std::vector<T> > hpp::operator* ( scalar,
const std::vector< std::vector< T >> &  veclist 
)
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA<T,M,N> hpp::operator* ( const Tensor2CUDA< T, M, N > &  A,
scalar 
)
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA<T,M,N> hpp::operator* ( scalar,
const Tensor2CUDA< T, M, N > &  A 
)
template<typename T , unsigned int M, unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA<T,M,P> hpp::operator* ( const Tensor2CUDA< T, M, N > &  A,
const Tensor2CUDA< T, N, P > &  B 
)
template<typename T , unsigned int M, unsigned int N>
__host__ __device__ VecCUDA<T,M> hpp::operator* ( const Tensor2CUDA< T, M, N > &  A,
const VecCUDA< T, N > &  x 
)
template<typename T , unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA<T,N,P> hpp::operator* ( const Tensor2SymmCUDA< T, N > &  A,
const Tensor2CUDA< T, N, P > &  B 
)
template<typename T , unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA<T,N,P> hpp::operator* ( const Tensor2AsymmCUDA< T, N > &  A,
const Tensor2CUDA< T, N, P > &  B 
)
template<typename T >
Tensor2<T> hpp::operator* ( const Tensor2< T > &  A,
const T &  B 
)
template<typename T >
Tensor2<T> hpp::operator* ( const T &  A,
const Tensor2< T > &  B 
)
template<typename T , unsigned int N, unsigned int P>
__host__ __device__ Tensor2CUDA<T,N,P> hpp::operator* ( const Tensor2CUDA< T, N, P > &  A,
const Tensor2AsymmCUDA< T, P > &  B 
)
template<typename T >
Tensor2<T> hpp::operator* ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T >
Tensor4<T> hpp::operator* ( const Tensor4< T > &  A,
const T &  B 
)
template<typename T >
Tensor4<T> hpp::operator* ( const T &  A,
const Tensor4< T > &  B 
)
template<typename T >
std::vector<T> hpp::operator* ( const hpp::Tensor2< T > &  A,
const std::vector< T > &  x 
)
template<typename T >
std::valarray<T> hpp::operator* ( const hpp::Tensor2< T > &  A,
const std::valarray< T > &  x 
)
template<typename T >
hpp::Tensor2<T> hpp::operator* ( const std::vector< T > &  x,
const hpp::Tensor2< T > &  A 
)

Replicates the numpy behaviour of vec*mat.

Parameters
xthe vector
Athe matrix
Returns
x*A
template<typename T >
void hpp::operator*= ( std::vector< T > &  vec,
const T  scalar 
)
template<typename T >
void hpp::operator*= ( Tensor2< T > &  A,
const T &  B 
)
template<typename T >
void hpp::operator*= ( Tensor4< T > &  A,
const T &  B 
)
template<typename T >
GSHCoeffs<T> hpp::operator+ ( const GSHCoeffs< T > &  coeffs1,
const GSHCoeffs< T > &  coeffs2 
)
template<typename T >
__host__ __device__ GSHCoeffsCUDA<T> hpp::operator+ ( const GSHCoeffsCUDA< T > &  coeffs1,
const GSHCoeffsCUDA< T > &  coeffs2 
)
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA<T,M,N> hpp::operator+ ( const Tensor2CUDA< T, M, N > &  A,
const Tensor2CUDA< T, M, N > &  B 
)
template<typename T , unsigned int N>
__host__ __device__ Tensor2SymmCUDA<T,N> hpp::operator+ ( const Tensor2SymmCUDA< T, N > &  A,
const Tensor2SymmCUDA< T, N > &  B 
)
template<typename T >
Tensor2<T> hpp::operator+ ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)

Addition.

Parameters
Afirst tensor
Bsecond tensor
Template Parameters
Tthe scalar type
Returns
A+b
template<typename T >
Tensor2<T> hpp::operator+ ( const Tensor2< T > &  A,
const T &  B 
)
template<typename T >
Tensor2<T> hpp::operator+ ( const T &  A,
const Tensor2< T > &  B 
)
template<typename T , unsigned int N>
__host__ __device__ Tensor2AsymmCUDA<T,N> hpp::operator+ ( const Tensor2AsymmCUDA< T, N > &  A,
const Tensor2AsymmCUDA< T, N > &  B 
)
template<typename T >
Tensor4<T> hpp::operator+ ( const Tensor4< T > &  A,
const Tensor4< T > &  B 
)

Addition.

Parameters
Afirst tensor
Bsecond tensor
Template Parameters
Tthe scalar type
Returns
A+B
template<typename T >
Tensor4<T> hpp::operator+ ( const Tensor4< T > &  A,
const T &  B 
)
template<typename T >
Tensor4<T> hpp::operator+ ( const T &  A,
const Tensor4< T > &  B 
)
template<typename T >
void hpp::operator+= ( GSHCoeffs< T > &  A,
const GSHCoeffs< T > &  B 
)
template<typename T >
__host__ __device__ void hpp::operator+= ( GSHCoeffsCUDA< T > &  A,
const GSHCoeffsCUDA< T > &  B 
)
template<typename T , unsigned int M, unsigned N>
__host__ __device__ void hpp::operator+= ( Tensor2CUDA< T, M, N > &  A,
const Tensor2CUDA< T, M, N > &  B 
)
template<typename T >
void hpp::operator+= ( Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T >
void hpp::operator+= ( Tensor2< T > &  A,
const T &  B 
)
template<typename T >
void hpp::operator+= ( Tensor4< T > &  A,
const Tensor4< T > &  B 
)
template<typename T >
void hpp::operator+= ( Tensor4< T > &  A,
const T &  B 
)
template<typename T >
std::vector<T> hpp::operator- ( const std::vector< T > &  vec1,
const std::vector< T > &  vec2 
)
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA<T,M,N> hpp::operator- ( const Tensor2CUDA< T, M, N > &  A,
const Tensor2CUDA< T, M, N > &  B 
)
template<typename T , unsigned int N>
__host__ __device__ Tensor2SymmCUDA<T,N> hpp::operator- ( const Tensor2SymmCUDA< T, N > &  A,
const Tensor2SymmCUDA< T, N > &  B 
)
template<typename T >
Tensor2<T> hpp::operator- ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T , unsigned int N>
__host__ __device__ Tensor2AsymmCUDA<T,N> hpp::operator- ( const Tensor2AsymmCUDA< T, N > &  A,
const Tensor2AsymmCUDA< T, N > &  B 
)
template<typename T , unsigned int N>
__host__ __device__ Tensor2CUDA<T,N,N> hpp::operator- ( const Tensor2AsymmCUDA< T, N > &  A,
const Tensor2CUDA< T, N, N > &  B 
)
template<typename T >
Tensor2<T> hpp::operator- ( const Tensor2< T > &  A,
const T &  B 
)
template<typename T >
Tensor2<T> hpp::operator- ( const T &  A,
const Tensor2< T > &  B 
)
template<typename T >
Tensor4<T> hpp::operator- ( const Tensor4< T > &  A)
template<typename T >
Tensor4<T> hpp::operator- ( const Tensor4< T > &  A,
const Tensor4< T > &  B 
)
template<typename T >
Tensor4<T> hpp::operator- ( const Tensor4< T > &  A,
const T &  B 
)
template<typename T >
Tensor4<T> hpp::operator- ( const T &  A,
const Tensor4< T > &  B 
)
template<typename T >
void hpp::operator-= ( Tensor2< T > &  A,
const Tensor2< T > &  B 
)
Parameters
A
B
Template Parameters
Tthe scalar type
template<typename T >
void hpp::operator-= ( Tensor2< T > &  A,
const T &  B 
)
template<typename T >
void hpp::operator-= ( Tensor4< T > &  A,
const Tensor4< T > &  B 
)
Parameters
A
B
Template Parameters
Tthe scalar type
template<typename T >
void hpp::operator-= ( Tensor4< T > &  A,
const T &  B 
)
template<typename T >
std::vector<T> hpp::operator/ ( const std::vector< T > &  vec,
scalar 
)

Elementwise division.

Parameters
vec
scalar
Returns
template<typename T >
std::vector<std::vector<T> > hpp::operator/ ( const std::vector< std::vector< T >> &  veclist,
scalar 
)
template<typename T , unsigned int N>
__device__ VecCUDA<T,N> hpp::operator/ ( const VecCUDA< T, N > &  inVec,
scalar 
)
template<typename T >
GSHCoeffs<T> hpp::operator/ ( const GSHCoeffs< T > &  coeffs,
const T  val 
)
template<typename T >
__host__ __device__ GSHCoeffsCUDA<T> hpp::operator/ ( const GSHCoeffsCUDA< T > &  coeffs,
val 
)
template<typename T , unsigned int M, unsigned N>
__host__ __device__ Tensor2CUDA<T,M,N> hpp::operator/ ( const Tensor2CUDA< T, M, N > &  A,
scalar 
)
template<typename T >
Tensor2<T> hpp::operator/ ( const Tensor2< T > &  A,
const T &  B 
)
template<typename T >
Tensor4<T> hpp::operator/ ( const Tensor4< T > &  A,
const T &  B 
)
template<typename T >
void hpp::operator/= ( std::vector< T > &  vec,
const T  scalar 
)
template<typename T >
void hpp::operator/= ( GSHCoeffs< T > &  A,
const T  B 
)
template<typename T , unsigned int M, unsigned N>
__host__ __device__ void hpp::operator/= ( Tensor2CUDA< T, M, N > &  A,
scalar 
)
template<typename T >
void hpp::operator/= ( Tensor2< T > &  A,
const T &  B 
)
template<typename T >
void hpp::operator/= ( Tensor4< T > &  A,
const T &  B 
)
bool hpp::operator< ( const SpectralDatasetID l,
const SpectralDatasetID r 
)
template<typename T >
std::ostream& hpp::operator<< ( std::ostream &  out,
const EulerAngles< T > &  angles 
)
template<typename T >
std::ostream& hpp::operator<< ( std::ostream &  out,
const std::vector< T > &  vec 
)
template<typename T >
std::ostream& hpp::operator<< ( std::ostream &  out,
const std::valarray< T > &  vec 
)
template<typename T , unsigned int M, unsigned int N>
std::ostream& hpp::operator<< ( std::ostream &  out,
const Tensor2CUDA< T, M, N > &  A 
)
template<typename T >
std::ostream& hpp::operator<< ( std::ostream &  out,
const Tensor2< T > &  A 
)
template<typename T >
std::ostream& hpp::operator<< ( std::ostream &  out,
Tensor4< T > &  A 
)
template<typename T >
bool hpp::operator== ( const SpectralCrystalCUDA< T > &  l,
const SpectralCrystalCUDA< T > &  r 
)
template<typename T >
bool hpp::operator== ( const EulerAngles< T > &  l,
const EulerAngles< T > &  r 
)
bool hpp::operator== ( const SpectralDatasetID l,
const SpectralDatasetID r 
)
template<typename U >
bool hpp::operator== ( const Crystal< U > &  l,
const Crystal< U > &  r 
)
template<typename T >
bool hpp::operator== ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)

Test equality of two tensors.

Uses the arbitrary \(10 \epsilon_{\mathrm{machine}}\) to determine if floats are equal. Test coverage in testTensor2Basics().

Parameters
Afirst tensor
Bsecond tensor
Template Parameters
Tthe scalar type
Returns
true if they're equal, false if they're not
template<typename T >
bool hpp::operator== ( const Tensor4< T > &  A,
const Tensor4< T > &  B 
)

Test equality of two tensors.

Uses the arbitrary \(10 \epsilon_{\mathrm{machine}}\) to determine if floats are equal. Test coverage in testTensor4Basics().

Parameters
Afirst tensor
Bsecond tensor
Template Parameters
Tthe scalar type
Returns
true if they're equal, false if they're not
template<typename T >
Tensor2<T> hpp::outer ( const std::vector< T > &  A,
const std::vector< T > &  B 
)
template<typename T >
Tensor2<T> hpp::outer ( const std::valarray< T > &  A,
const std::valarray< T > &  B 
)
template<typename T >
Tensor4<T> hpp::outer ( const Tensor2< T > &  A,
const Tensor2< T > &  B 
)
template<typename T >
void hpp::outerInPlace ( const std::vector< T > &  A,
const std::vector< T > &  B,
Tensor2< T > &  C 
)
template<typename T >
void hpp::outerInPlace ( const Tensor2< T > &  A,
const Tensor2< T > &  B,
Tensor4< T > &  C 
)
template<typename U >
hpp::Tensor2<U> hpp::partialDGammaPartialT ( m,
gammadot_0,
const hpp::Tensor2< U > &  T,
const hpp::Tensor2< U > &  S_0_alpha,
s_alpha,
dt 
)

See kalidindi1992.PartialDGammaPartialT.

Parameters
m\(m\)
gammadot_0\(\dot{\gamma}_0\)
T\(\mathbf{T}^*(t_{i+1})\)
S_0_alpha\(\mathbf{S}_0^\alpha\)
s_alpha\(s^\alpha(t_{i+1})\)
dt\(\Delta t\)
Returns
template<typename U >
void hpp::partialDGammaPartialTInPlace ( m,
gammadot_0,
const hpp::Tensor2< U > &  T,
const hpp::Tensor2< U > &  S_0_alpha,
s_alpha,
dt,
hpp::Tensor2< U > &  pdgpt 
)

In-place version of hpp.partialDGammaPartialT.

Parameters
m\(m\)
gammadot_0\(\dot{\gamma}_0\)
T\(\mathbf{T}^*(t_{i+1})\)
S_0_alpha\(\mathbf{S}_0^\alpha\)
s_alpha\(s^\alpha(t_{i+1})\)
dt\(\Delta t\)
pdgpt
template<typename U >
Tensor2<U> hpp::planeStrainCompressionDeformationGradient ( t,
comp_rate 
)

The deformation gradient from a plane strain compression.

The deformation gradient is zero with the addition of \(\mathbf{F}(0,0) = \exp(1.0\dot{\varepsilon} t)\), \(\mathbf{F}(1,1) = 1.0\), \(\mathbf{F}(2,2) = \exp(-1.0\dot{\varepsilon} t)\).

Parameters
tThe current time \(t\)
comp_rate
Template Parameters
Uthe scalar type
Returns
\( \mathbf{F} \)
template<typename U >
Tensor2<U> hpp::planeStrainCompressionVelocityGradient ( t,
comp_rate 
)

The velocity gradient from a plane strain compression.

The velocity gradient is zero with the addition of \(\mathbf{L}(0,0) = 1.0\dot{\varepsilon}\), \(\mathbf{L}(2,2) = -1.0\dot{\varepsilon}\).

Parameters
tThe current time \(t\)
comp_ratethe compression rate \( \dot{\varepsilon}\)
Template Parameters
Uthe scalar type
Returns
\( \mathbf{L} \)
template<typename U >
hpp::Tensor2<U> hpp::plasticDeformationGradientUpdate ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  F_p_prev_time,
const std::vector< U > &  Dgamma_alphas 
)
template<typename U >
U hpp::plasticShearingRate ( tau_alpha,
s_alpha,
gammadot_0,
m 
)

See kalidindi1992.PlasticShearingRate.

Parameters
tau_alpha\(\tau^\alpha\)
s_alpha\(s^\alpha\)
gammadot_0\(\dot{\gamma}_0\)
m\(m\)
Returns
\(\dot{\gamma}^\alpha\)
template<typename U >
Tensor2<U> hpp::plasticSpinTensor ( CrystalProperties< U >  mprops,
std::vector< U >  gammadot_alphas,
std::vector< std::vector< U >>  m_alphas,
std::vector< std::vector< U >>  n_alphas 
)

Mihaila 2014 equation 14.

Parameters
mprops
gammadot_alphas
m_alphasan array of \(\mathbf{m}^\alpha\) for each slip system \(\alpha\)
n_alphasan array of \(\mathbf{n}^\alpha\) for each slip system \(\alpha\)
Returns
template<typename T >
EulerAngles<T> hpp::polarToEuler ( theta,
phi 
)

Convert polar angles to Euler angles.

Parameters
thetathe azimuthal angle
phithe zenithal angle
Returns
The Euler angles
FFTWConfigRealND hpp::prepareFFTWConfigRealND ( const std::vector< ptrdiff_t > &  realDims,
MPI_Comm  comm 
)
template<typename T >
T hpp::prod ( const std::vector< T > &  vec)

Product of all the entries of a vector.

Parameters
vecthe vector
Returns
the product
template<typename T >
EulerAngles<T> hpp::quaternionToEulerAngles ( isoi::Quaternion &  q)
template<typename T >
Tensor2<T> hpp::randomRotationTensor ( unsigned int  dim,
bool  defaultSeed = false 
)
template<typename T >
void hpp::randomRotationTensorArvo1992 ( std::mt19937 &  gen,
Tensor2< T > &  A 
)
template<typename T >
void hpp::randomRotationTensorArvo1992 ( Tensor2< T > &  A,
bool  defaultSeed = false 
)
template<typename T >
void hpp::randomRotationTensorInPlace ( unsigned int  dim,
Tensor2< T > &  A,
bool  defaultSeed = false 
)

Generate a random rotation tensor.

Parameters
dim
Athe output tensor to populate
defaultSeedif true, use a default random seed, otherwise generate a truly random one
template<typename T >
Tensor2<T> hpp::randomRotationTensorStewart1980 ( unsigned int  dim,
bool  defaultSeed = false 
)

Generate a random rotation tensor.

This is translated from scipy/linalg/tests/test_decomp.py The description there, verbatim is: Return a random rotation matrix, drawn from the Haar distribution (the only uniform distribution on SO(n)). The algorithm is described in the paper Stewart, G.W., 'The efficient generation of random orthogonal matrices with an application to condition estimators', SIAM Journal on Numerical Analysis, 17(3), pp. 403-409, 1980. For more information see http://en.wikipedia.org/wiki/Orthogonal_matrix#Randomization

Parameters
dim
defaultSeedif true, use a default random seed, otherwise generate a truly random one
Returns
The rotation tensor
template<typename T >
void hpp::readHDF5SimpleArray ( H5::DataSet &  dataset,
HDFReadWriteParams  parms,
T *  output 
)
template<typename T >
void hpp::readHDF5SimpleArray ( hid_t  dset_id,
hid_t  plist_id,
HDFReadWriteParamsC  parms,
T *  output 
)
template<typename T >
void hpp::readMultipleHDF5Arrays ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims,
std::vector< hsize_t >  arrayCount,
T *  output 
)
template<typename T >
void hpp::readSingleHDF5Array ( H5::DataSet &  dataset,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims,
T *  output 
)
template<typename T >
void hpp::readSingleHDF5Array ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims,
T *  output 
)
template<typename T >
void hpp::readSingleHDF5Array ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  arrayDims,
T *  output 
)
template<typename T >
void hpp::readSingleHDF5Value ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  gridOffset,
T *  output 
)
template<typename T >
void hpp::readWriteHDF5SimpleArray ( H5::DataSet &  dataset,
HDFReadWriteParams  parms,
T *  output,
HDFReadWrite  mode 
)
template<typename T >
void hpp::readWriteHDF5SimpleArray ( hid_t  dset_id,
hid_t  plist_id,
HDFReadWriteParamsC  parms,
T *  output,
HDFReadWrite  mode 
)
template<typename U >
CrystalProperties<U> hpp::rotate ( const CrystalProperties< U > &  propsOld,
hpp::Tensor2< U >  rotTensor 
)
template<typename T >
void hpp::rotationTensorFrom3UniformRandoms ( Tensor2< T > &  A,
x1,
x2,
x3 
)
template<typename T >
void hpp::rotationTensorFrom3UniformRandomsArvo1992 ( Tensor2< T > &  A,
x1,
x2,
x3 
)
unsigned int hpp::roundUpToMultiple ( unsigned int  val,
unsigned int  multiple 
)
template<typename U >
U hpp::setTimestepByShearRate ( dt_old,
const std::vector< U > &  Dgamma_alphas,
Dgamma_goal 
)

The exact method in kalidindi1992.

Equation 38 in kalidindi1992

Parameters
dt_oldthe old timestep
Dgamma_alphasthe strain increments
Dgamma_goalthe strain increment goal
Returns
The new timestep
template<typename U >
std::vector<U> hpp::shearStrainIncrements ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  T,
const std::vector< U > &  s_alphas,
const U  dt 
)
inline

See kalidindi1992.ShearStrainIncrements.

Parameters
propsthe properties of the material defined in hpp.CrystalProperties
T\(\mathbf{T}^*\)
s_alphasa list of \(s^\alpha(t_{i+1})\) for each slip system \(\alpha\)
dt\(\Delta t\)
Returns
A list of \(\Delta \gamma^\alpha\) for each \(\alpha\)
template<typename U >
void hpp::shearStrainIncrementsInPlace ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  T,
const std::vector< U > &  s_alphas,
const U  dt,
std::vector< U > &  Dgamma_alphas 
)
inline

In-place version of hpp.shearStrainIncrements.

Parameters
propsthe properties of the material defined in hpp.CrystalProperties
T\(\mathbf{T}^*\)
s_alphasa list of \(s^\alpha(t_{i+1})\) for each slip system \(\alpha\)
dt\(\Delta t\)
Dgamma_alphaslist of \(\Delta \gamma^\alpha\) for each \(\alpha\) to return
template<typename U >
std::vector<U> hpp::shearStrainRates ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  T,
const std::vector< U > &  s_alphas 
)

See kalidindi1992.ShearStrainRates.

Parameters
propsthe properties of the material defined in hpp.CrystalProperties
T\(\mathbf{T}^*\)
s_alphasa list of \(s^\alpha(t_{i+1})\) for each slip system \(\alpha\)
Returns
An array of \(\dot{\gamma}^\alpha\) for each slip system \(\alpha\)
template<typename U >
void hpp::shearStrainRatesInPlace ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  T,
const std::vector< U > &  s_alphas,
std::vector< U > &  gammadot_alphas 
)
inline

In-place version of hpp.shearStrainRates.

Parameters
propsthe properties of the material defined in hpp.CrystalProperties
T\(\mathbf{T}^*\)
s_alphasa list of \(s^\alpha(t_{i+1})\) for each slip system \(\alpha\)
gammadot_alphasAn array of \(\dot{\gamma}^\alpha\) for each slip system \(\alpha\) that is returned
template<typename U >
Tensor2<U> hpp::simpleCompressionDeformationGradient ( t,
comp_rate 
)

The deformation gradient from a simple compression.

The deformation gradient is zero with the addition of \(\mathbf{F}(0,0) = \exp(-0.5\dot{\varepsilon} t)\), \(\mathbf{F}(1,1) = \exp(-0.5\dot{\varepsilon} t)\), \(\mathbf{F}(2,2) = \exp(1.0\dot{\varepsilon} t)\).

Parameters
tThe current time \(t\)
comp_ratethe compression rate \( \dot{\varepsilon}\)
Template Parameters
Uthe scalar type
Returns
\( \mathbf{F} \)
template<typename U >
Tensor2<U> hpp::simpleCompressionVelocityGradient ( t,
comp_rate 
)

The velocity gradient from a simple compression.

The velocity gradient is zero with the addition of \(\mathbf{L}(0,0) = -0.5\dot{\varepsilon}\), \(\mathbf{L}(1,1) = -0.5\dot{\varepsilon}\), \(\mathbf{L}(2,2) = 1.0\dot{\varepsilon}\).

Parameters
tThe current time \(t\)
comp_ratethe compression rate \( \dot{\varepsilon}\)
Template Parameters
Uthe scalar type
Returns
\( \mathbf{L} \)
template<typename U >
Tensor2<U> hpp::simpleShearDeformationGradient ( t,
shear_rate 
)

The deformation gradient from a simple shear.

The deformation gradient is the identity with the addition of \( \mathbf{F}(0,1) = \dot{\varepsilon} t \).

Parameters
tThe current time \(t\)
shear_rateThe shear rate \( \dot{\varepsilon}\)
Template Parameters
Uthe scalar type
Returns
\( \mathbf{F} \)
template<typename U >
Tensor2<U> hpp::simpleShearVelocityGradient ( t,
shear_rate 
)

The velocity gradient from a simple shear.

The velocity gradient is zero with the addition of \( \mathbf{L}(0,1) = \dot{\varepsilon} \).

Parameters
tThe current time \(t\)
shear_rateThe shear rate \( \dot{\varepsilon}\)
Template Parameters
Uthe scalar type
Returns
\( \mathbf{L} \)
template<typename T , unsigned int N>
__device__ T hpp::slipDeformationResistanceStepSpectralSolver ( const CrystalPropertiesCUDA< T, N > *  props,
const T  s_alpha,
const T  gammaSum,
const T  dt 
)
template<typename U >
U hpp::slipDeformationResistanceStepSpectralSolver ( const CrystalProperties< U > &  props,
const U  s_alpha,
const U  gammaSum,
const U  dt 
)

From kalidindi2006 Equation 5.

Parameters
props
s_alpha
gammaSum
dt
Returns
template<typename U >
std::vector<U> hpp::slipDeformationResistanceUpdate ( const CrystalProperties< U > &  props,
const hpp::Tensor2< U > &  T_next,
const std::vector< U > &  s_alphas_current_time,
const std::vector< U > &  s_alphas_prev_iter,
const U  dt 
)
template<typename U >
U hpp::slipHardeningRate ( h_0,
s_s,
a,
s_beta 
)
template<typename T , unsigned int N>
__global__ void hpp::SPECTRAL_POLYCRYSTAL_STEP ( unsigned int  nCrystals,
SpectralCrystalCUDA< T > *  crystals,
CrystalPropertiesCUDA< T, N > *  props,
Tensor2CUDA< T, 3, 3 >  RStretchingTensor,
Tensor2AsymmCUDA< T, 3 >  WNext,
theta,
eDot,
dt,
SpectralDatabaseCUDA< T, 4 > *  db,
Tensor2CUDA< T, 3, 3 > *  TCauchyPerBlockSums 
)
template<typename T , unsigned int N, unsigned int P>
__global__ void hpp::SPECTRAL_POLYCRYSTAL_STEP_UNIFIED ( unsigned int  nCrystals,
SpectralCrystalCUDA< T > *  crystals,
CrystalPropertiesCUDA< T, N > *  props,
Tensor2CUDA< T, 3, 3 >  RStretchingTensor,
Tensor2AsymmCUDA< T, 3 >  WNext,
theta,
eDot,
dt,
SpectralDatabaseUnifiedCUDA< T, 4, P > *  db,
Tensor2CUDA< T, 3, 3 > *  TCauchyPerBlockSums 
)
template<typename U >
Tensor2<U> hpp::staticDeformationGradient ( t)

The deformation gradient for no deformation.

The deformation gradient is the identity.

Parameters
tThe current time \(t\)
Template Parameters
Uthe scalar type
Returns
\( \mathbf{F} \)
template<typename U >
Tensor2<U> hpp::staticVelocityGradient ( t)

The velocity gradient for no deformation.

The velocity gradient is zero.

Parameters
tThe current time \(t\)
Template Parameters
Uthe scalar type
Returns
\( \mathbf{L} \)
template<typename U >
hpp::Tensor2<U> hpp::strainHardeningRates ( const CrystalProperties< U > &  props,
const std::vector< U > &  s_alphas 
)
template<typename U >
U hpp::strainHardeningRateVoce ( const CrystalProperties< U > &  props,
const std::vector< U > &  s_alphas 
)

Voce hardening law.

identical for each slip system

Parameters
props
s_alphas
Returns
template<typename U >
bool hpp::strainRateLowEnough ( const std::vector< U > &  Dgamma_alphas,
Dgamma_goal 
)
template<typename U >
Tensor2<U> hpp::stretchingVelocityGradient ( theta,
eDot 
)
template<typename U >
hpp::Tensor2<U> hpp::tensorA ( hpp::Tensor2< U >  F_p,
hpp::Tensor2< U >  F_next 
)

See kalidindi1992.TensorA.

Parameters
F_p\(\mathbf{F}^{p}(t_{i})\)
F_next\(\mathbf{F}(t_{i+1})\)
Template Parameters
Uthe scalar type
Returns
\(\mathbf{A}\)
template<typename U >
std::vector<hpp::Tensor2<U> > hpp::tensorC_alphas ( const hpp::Tensor4< U > &  L,
const hpp::Tensor2< U > &  A,
const std::vector< hpp::Tensor2< U >> &  S_0,
const unsigned int  n_alpha 
)

See kalidindi1992.TensorC_alphas.

Parameters
L\(\mathcal{L}\)
A\(\mathbf{A}\)
S_0List of \(\mathbf{S}_0^\alpha\)
n_alphathe number of slip systems
Template Parameters
Uthe scalar type
Returns
\(\mathbf{C}^\alpha\)
template<typename U >
void hpp::tensorC_alphasInPlace ( const hpp::Tensor4< U > &  L,
const hpp::Tensor2< U > &  A,
const std::vector< hpp::Tensor2< U >> &  S_0,
const unsigned int  n_alpha,
hpp::Tensor2< U > &  dumB_alpha,
std::vector< hpp::Tensor2< U >> &  C_alphas 
)

In-place version of hpp.tensorC_alphas.

The input C_alphas do not need to be zeroed

Parameters
L\(\mathcal{L}\)
A\(\mathbf{A}\)
S_0List of \(\mathbf{S}_0^\alpha\) for each \(\alpha\)
n_alphathe number of slip systems
dumB_alphaa dummy tensor of size 3x3 for holding intermediate values of \(\mathbf{B}^\alpha\)
C_alphaslist of \(\mathbf{C}^\alpha\) for each \(\alpha\) to be returned
Template Parameters
Uthe scalar type
template<typename U >
hpp::Tensor2<U> hpp::tensorG ( const hpp::Tensor4< U > &  L,
const hpp::Tensor2< U > &  A,
const hpp::Tensor2< U > &  T_prev_iter,
const std::vector< U > &  Dgamma_alphas,
const std::vector< hpp::Tensor2< U >> &  C_alphas 
)
inline

See kalidindi1992.TensorG.

template<typename U >
void hpp::tensorGInPlace ( const hpp::Tensor4< U > &  L,
const hpp::Tensor2< U > &  A,
const hpp::Tensor2< U > &  T_prev_iter,
const std::vector< U > &  Dgamma_alphas,
const std::vector< hpp::Tensor2< U >> &  C_alphas,
Tensor2< U > &  dum2ndOrder,
hpp::Tensor2< U > &  dumT_tr,
hpp::Tensor2< U > &  G 
)
inline

In-place version of hpp.TensorG.

template<typename U >
hpp::Tensor4<U> hpp::tensorJ ( const CrystalProperties< U > &  props,
const std::vector< hpp::Tensor2< U >> &  C_alphas,
const std::vector< U > &  s_alphas,
const hpp::Tensor2< U > &  T,
const U  dt 
)
inline

See kalidindi1992.TensorJ.

Parameters
propsthe properties of the material defined in hpp.CrystalProperties
C_alphaslist of \(\mathbf{C}^\alpha\) for each \(\alpha\)
s_alphasa list of \(s^\alpha(t_{i+1})\) for each slip system \(\alpha\)
T\(\mathbf{T}^*(t_{i+1})\)
dt\(\Delta t\)
Returns
\(\mathbf{J}\)
template<typename U >
void hpp::tensorJInPlace ( const CrystalProperties< U > &  props,
const std::vector< hpp::Tensor2< U >> &  C_alphas,
const std::vector< U > &  s_alphas,
const hpp::Tensor2< U > &  T,
const U  dt,
hpp::Tensor2< U > &  dum2ndOrder,
hpp::Tensor4< U > &  dum4thOrder,
hpp::Tensor4< U > &  J 
)
inline

In-place version of hpp.tensorJ.

See kalidindi1992.TensorJ.

Parameters
propsthe properties of the material defined in hpp.CrystalProperties
C_alphaslist of \(\mathbf{C}^\alpha\) for each \(\alpha\)
s_alphasa list of \(s^\alpha(t_{i+1})\) for each slip system \(\alpha\)
T\(\mathbf{T}^*(t_{i+1})\)
dt\(\Delta t\)
Returns
\(\mathbf{J}\)
Parameters
dum2ndOrdera 3x3 dummy tensor
dum4thOrdera 3x3x3x3 dummy tensor
Ja 3x3x3x3 tensor for the output of J
template<class T >
boost::python::list hpp::toPythonList ( const std::vector< T > &  vec)
template<class T >
std::vector<T> hpp::toStdVector ( const boost::python::list &  list)
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA<T,M,M> hpp::transformIntoFrame ( const Tensor2CUDA< T, M, M > &  A,
const Tensor2CUDA< T, M, M > &  Q 
)

Transform tensor \( \mathbf{A} \) into the frame given by the columns of \( \mathbf{Q} \).

Parameters
A\( \mathbf{A} \)
Q\( \mathbf{Q} \)
Returns
\( \mathbf{A}^* \)
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA<T,M,M> hpp::transformIntoFrame ( const Tensor2AsymmCUDA< T, M > &  A,
const Tensor2CUDA< T, M, M > &  Q 
)
template<typename T >
Tensor2<T> hpp::transformIntoFrame ( const Tensor2< T > &  A,
const Tensor2< T > &  Q 
)

Transform tensor \( \mathbf{A} \) into the frame given by the columns of \( \mathbf{Q} \).

Parameters
A\( \mathbf{A} \)
Q\( \mathbf{Q} \)
Returns
\( \mathbf{A}^* \)
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA<T,M,M> hpp::transformOutOfFrame ( const Tensor2CUDA< T, M, M > &  A_star,
const Tensor2CUDA< T, M, M > &  Q 
)

Transform tensor \( \mathbf{A}^* \) out of the frame given by the columns of \( \mathbf{Q} \).

Parameters
A_star\( \mathbf{A}^* \)
Q\( \mathbf{Q} \)
Returns
\( \mathbf{A} \)
template<typename T , unsigned int M>
__host__ __device__ Tensor2CUDA<T,M,M> hpp::transformOutOfFrame ( const Tensor2AsymmCUDA< T, M > &  A_star,
const Tensor2CUDA< T, M, M > &  Q 
)
template<typename T >
Tensor2<T> hpp::transformOutOfFrame ( const Tensor2< T > &  A_star,
const Tensor2< T > &  Q 
)

Transform tensor \( \mathbf{A}^* \) out of the frame given by the columns of \( \mathbf{Q} \).

Parameters
A_star\( \mathbf{A}^* \)
Q\( \mathbf{Q} \)
Returns
\( \mathbf{A} \)
template<typename T >
Tensor4<T> hpp::transformOutOfFrame ( const Tensor4< T > &  E_star,
const Tensor2< T > &  Q 
)

Transform tensor \( \mathbf{E}^* \) out of the frame given by the columns of \( \mathbf{Q} \).

Parameters
E_star\( \mathbf{E}^* \)
Q\( \mathbf{Q} \)
Returns
\( \mathbf{E} \)
Todo:
Use symmetries to reduce the computation
template<typename U >
U hpp::TruncatedIFFT4DSingleSquareReal ( const unsigned int  spatialDim,
const U  expTable[],
const int  spectralCoordsList[],
const U  coeffs[],
const std::vector< unsigned int > &  spatialCoord,
unsigned int  nTerms 
)

AVX implementation of hpp.TruncatedIFFT4DSingleSquareReal.

AVX2 implementation of hpp.TruncatedIFFT4DSingleSquareReal

Parameters
spatialDim
expTable16/32-byte (float/double) aligned array of real-im pairs from the exp table
spectralCoordsList16/32-byte (float/double) aligned array of spectral indices
coeffs16/32-byte (float/double) aligned array of real-im pairs from the coeffs list in the same order as the spectral indices
spatialCoordthe spatial coordinate
nTermsthe number of Fourier terms to use
Returns
template<typename U >
U hpp::TruncatedIFFT4DSingleSquareReal ( const unsigned int  spatialDim,
const int  spectralCoordsList[],
const U  coeffs[],
const std::vector< unsigned int > &  spatialCoord,
unsigned int  nTerms 
)
Todo:
CPU SIMD implementations for arbitrary dimension
template<typename U >
std::complex<U> hpp::TruncatedIFFTNDSingleSquare ( const unsigned int  spatialDim,
const std::vector< std::complex< U >> &  expTable,
const std::vector< std::vector< unsigned int >> &  spectralCoordsList,
const std::vector< std::complex< U >> &  coeffs,
std::vector< unsigned int >  spatialCoord,
unsigned int  nTerms 
)

TruncatedIFFTND where all dimensions have the same number of points.

AVX implementation of hpp.TruncatedIFFT4DSingleSquareReal. AVX2 implementation of hpp.TruncatedIFFT4DSingleSquareReal

Parameters
spatialDim
expTable16/32-byte (float/double) aligned array of real-im pairs from the exp table
spectralCoordsList16/32-byte (float/double) aligned array of spectral indices
coeffs16/32-byte (float/double) aligned array of real-im pairs from the coeffs list in the same order as the spectral indices
spatialCoordthe spatial coordinate
nTermsthe number of Fourier terms to use
Returns
template<typename U >
U hpp::TruncatedIFFTNDSingleSquareReal ( const unsigned int  spatialDim,
const std::vector< std::complex< U >> &  expTable,
const std::vector< std::vector< unsigned int >> &  spectralCoordsList,
const std::vector< std::complex< U >> &  coeffs,
std::vector< unsigned int >  spatialCoord,
unsigned int  nTerms 
)
template<typename U >
U hpp::TruncatedIFFTNDSingleSquareReal ( const unsigned int  nDims,
const unsigned int  spatialDim,
const U  expTable[],
const int  spectralCoordsList[],
const U  coeffs[],
const std::vector< unsigned int > &  spatialCoord,
unsigned int  nTerms 
)
Todo:
CPU SIMD implementations for arbitrary dimension
template<typename U >
U hpp::TruncatedIFFTNDSingleSquareReal ( const unsigned int  nDims,
const unsigned int  spatialDim,
const int  spectralCoordsList[],
const U  coeffs[],
const std::vector< unsigned int > &  spatialCoord,
unsigned int  nTerms 
)
Todo:
CPU SIMD implementations for arbitrary dimension
template<typename U >
U hpp::TruncatedIFFTNDSingleSquareRealScalar ( const unsigned int  nDims,
const unsigned int  spatialDim,
const U  expTable[],
const int  spectralCoordsList[],
const U  coeffs[],
std::vector< unsigned int >  spatialCoord,
unsigned int  nTerms 
)

Scalar implementation of hpp.TruncatedIFFTNDSingleSquareReal.

No guarantees that the compiler won't use vector instructions.

template<typename U >
U hpp::TruncatedIFFTNDSingleSquareRealScalar ( const unsigned int  nDims,
const unsigned int  spatialDim,
const int  spectralCoordsList[],
const U  coeffs[],
std::vector< unsigned int >  spatialCoord,
unsigned int  nTerms 
)
idx2d hpp::unflat ( unsigned int  flatIdx,
unsigned int  n1,
unsigned int  n2 
)
inline
idx4d hpp::unflat ( unsigned int  flatIdx,
unsigned int  n1,
unsigned int  n2,
unsigned int  n3,
unsigned int  n4 
)
inline
template<typename T , typename U >
std::vector<T> hpp::unflatC ( flatIdx,
std::vector< U >  dims 
)
template<typename U >
Tensor4<U> hpp::volumeAverage ( const std::vector< Tensor4< U >> &  tVec,
const std::vector< U > &  vVec 
)

Volume average fourth order tensors.

Parameters
tVeca vector of the tensors
vVeca vector of the volumes in which each tensor holds
Returns
the volume average
template<typename T >
__device__ GSHCoeffsCUDA<T> hpp::warpReduceSumGSHCoeffs ( GSHCoeffsCUDA< T >  coeffs)
inline

Note that this is valid only on CC 3.0 and above. Adapted from https://devblogs.nvidia.com/parallelforall/faster-parallel-reductions-kepler/

Parameters
val
Returns
template<typename T , unsigned int M, unsigned N>
__device__ Tensor2CUDA<T,M,N> hpp::warpReduceSumTensor2 ( Tensor2CUDA< T, M, N >  A)
inline

Note that this is valid only on CC 3.0 and above. Adapted from https://devblogs.nvidia.com/parallelforall/faster-parallel-reductions-kepler/

Parameters
val
Returns
template<typename T >
void hpp::writeHDF5SimpleArray ( H5::DataSet &  dataset,
HDFReadWriteParams  parms,
T *  output 
)
template<typename T >
void hpp::writeHDF5SimpleArray ( hid_t  dset_id,
hid_t  plist_id,
HDFReadWriteParamsC  parms,
T *  output 
)
template<typename T >
void hpp::writeMultipleHDF5Arrays ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims,
std::vector< hsize_t >  arrayCount,
T *  output 
)
template<typename T >
void hpp::writePoleHistogramHistoryHDF5 ( H5::H5File &  outfile,
const std::string &  dsetBaseName,
std::vector< Tensor2< T >> &  history,
const std::vector< T > &  pole 
)
Parameters
outfile
dsetBaseName
history
pole
template<typename T >
void hpp::writeSingleHDF5Array ( H5::DataSet &  dataset,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims,
T *  output 
)
template<typename T >
void hpp::writeSingleHDF5Array ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  gridOffset,
std::vector< hsize_t >  arrayDims,
T *  output 
)
template<typename T >
void hpp::writeSingleHDF5Array ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  arrayDims,
T *  output 
)
template<typename T >
void hpp::writeSingleHDF5Value ( hid_t  dset_id,
hid_t  plist_id,
std::vector< hsize_t >  gridOffset,
T *  output 
)
template<typename T >
void hpp::writeVectorToHDF5Array ( H5::H5File &  file,
const std::string &  dsetName,
std::vector< T > &  vec 
)

Variable Documentation

void * hpp::dummyHDFClientData
H5E_auto2_t hpp::dummyHDFErrorHandler
hid_t hpp::hdf_complex_id = getComplexHDFType()