PSCF v1.2

Element-wise vector operations performed on the GPU for FTS classes. More...

Functions

void mcftsScale (DeviceArray< Prdc::Cuda::cudaReal > &a, Prdc::Cuda::cudaReal const b)
 Rescale array a from [0,1] to [-b, b], GPU kernel wrapper.
 
void fourierMove (DeviceArray< Prdc::Cuda::cudaComplex > &a, DeviceArray< Prdc::Cuda::cudaReal > const &b, DeviceArray< Prdc::Cuda::cudaReal > const &c)
 Add array b to real part of a and array c to imaginary part of a.
 
void computeDField (DeviceArray< Prdc::Cuda::cudaReal > &d, DeviceArray< Prdc::Cuda::cudaReal > const &Wc, DeviceArray< Prdc::Cuda::cudaReal > const &Cc, Prdc::Cuda::cudaReal const a, Prdc::Cuda::cudaReal const b, Prdc::Cuda::cudaReal const s)
 Compute d field (functional derivative of H[w])
 
void computeForceBias (DeviceArray< Prdc::Cuda::cudaReal > &result, DeviceArray< Prdc::Cuda::cudaReal > const &di, DeviceArray< Prdc::Cuda::cudaReal > const &df, DeviceArray< Prdc::Cuda::cudaReal > const &dwc, Prdc::Cuda::cudaReal mobility)
 Compute force bias.
 

Detailed Description

Element-wise vector operations performed on the GPU for FTS classes.

CUDA kernels that perform the operations are defined in VecOpFts.cu in an anonymous namespace, so they are not directly accessible. Kernel wrapper functions to be called by the host CPU, which call the kernel internally, are public.

The output (the LHS of the vector operation) will always be the first parameter passed to the function.

Function Documentation

◆ mcftsScale()

void Pscf::Rpg::VecOpFts::mcftsScale ( DeviceArray< Prdc::Cuda::cudaReal > & a,
Prdc::Cuda::cudaReal const b )

Rescale array a from [0,1] to [-b, b], GPU kernel wrapper.

Parameters
aarray containing data between 0 and 1 (LHS)
bscalar bound for rescaled array (RHS)

Definition at line 77 of file VecOpFts.cu.

References Pscf::DeviceArray< Data >::capacity(), Pscf::DeviceArray< Data >::cArray(), and Pscf::ThreadArray::setThreadsLogical().

Referenced by Pscf::Rpg::FourierMove< D >::attemptMove(), and Pscf::Rpg::RealMove< D >::attemptMove().

◆ fourierMove()

void Pscf::Rpg::VecOpFts::fourierMove ( DeviceArray< Prdc::Cuda::cudaComplex > & a,
DeviceArray< Prdc::Cuda::cudaReal > const & b,
DeviceArray< Prdc::Cuda::cudaReal > const & c )

Add array b to real part of a and array c to imaginary part of a.

Parameters
aoutput array of cudaComplex values
binput array of cudaReal values
cinput array of cudaReal values

Definition at line 92 of file VecOpFts.cu.

References Pscf::DeviceArray< Data >::capacity(), Pscf::DeviceArray< Data >::cArray(), Pscf::ThreadArray::setThreadsLogical(), and UTIL_CHECK.

Referenced by Pscf::Rpg::FourierMove< D >::attemptMove().

◆ computeDField()

void Pscf::Rpg::VecOpFts::computeDField ( DeviceArray< cudaReal > & d,
DeviceArray< cudaReal > const & Wc,
DeviceArray< cudaReal > const & Cc,
cudaReal const a,
cudaReal const b,
cudaReal const s )

Compute d field (functional derivative of H[w])

Definition at line 112 of file VecOpFts.cu.

References Pscf::DeviceArray< Data >::capacity(), Pscf::DeviceArray< Data >::cArray(), Pscf::ThreadArray::setThreadsLogical(), and UTIL_CHECK.

Referenced by Pscf::Rpg::Simulator< D >::computeDc().

◆ computeForceBias()

void Pscf::Rpg::VecOpFts::computeForceBias ( DeviceArray< cudaReal > & result,
DeviceArray< cudaReal > const & di,
DeviceArray< cudaReal > const & df,
DeviceArray< cudaReal > const & dwc,
cudaReal mobility )