C++ API Reference¶
-
class DerivativeBase¶
- #include <derivative.hh>
A base class that represents a derivative.
This class provides the basic functionalities for a derivative.
Subclassed by muFFT::DiscreteDerivative, muFFT::FourierDerivative
Public Types
-
using Vector = Eigen::Matrix<Real, Eigen::Dynamic, 1>¶
Alias for Eigen::Matrix<Real, Eigen::Dynamic, 1>
Public Functions
-
DerivativeBase() = delete¶
Deleted default constructor.
This constructor is deleted because a DerivativeBase object requires a spatial dimension to be properly initialized.
-
explicit DerivativeBase(Index_t spatial_dimension)¶
Constructor that takes the spatial dimension as an argument.
- Parameters:
spatial_dimension – The spatial dimension of the derivative.
-
DerivativeBase(const DerivativeBase &other) = default¶
Default copy constructor.
-
DerivativeBase(DerivativeBase &&other) = default¶
Default move constructor.
-
virtual ~DerivativeBase() = default¶
Default destructor.
-
DerivativeBase &operator=(const DerivativeBase &other) = delete¶
Deleted copy assignment operator.
-
DerivativeBase &operator=(DerivativeBase &&other) = delete¶
Deleted move assignment operator.
Protected Attributes
-
Index_t spatial_dimension¶
The spatial dimension of the problem.
-
using Vector = Eigen::Matrix<Real, Eigen::Dynamic, 1>¶
-
class DerivativeError : public RuntimeError¶
- #include <derivative.hh>
A class that represents exceptions related to derivatives.
This class is derived from the RuntimeError class and is used to handle exceptions related to derivatives.
Public Functions
-
inline explicit DerivativeError(const std::string &what)¶
A constructor that takes a string as an argument.
- Parameters:
what – A string that describes the error.
-
inline explicit DerivativeError(const char *what)¶
A constructor that takes a character array as an argument.
- Parameters:
what – A character array that describes the error.
-
inline explicit DerivativeError(const std::string &what)¶
-
class DiscreteDerivative : public muFFT::DerivativeBase¶
- #include <derivative.hh>
A class that represents a discrete derivative.
This class is derived from the DerivativeBase class and provides functionalities for discrete derivatives.
Public Types
-
using Parent = DerivativeBase¶
Base class alias.
Public Functions
-
DiscreteDerivative() = delete¶
Default constructor is deleted as a DiscreteDerivative object requires raw stencil information for proper initialization
-
DiscreteDerivative(DynCcoord_t nb_pts, DynCcoord_t lbounds, const std::vector<Real> &stencil)¶
Constructor with raw stencil information.
This constructor initializes a DiscreteDerivative object with the provided stencil size, relative starting point of stencil, and stencil coefficients.
- Parameters:
nb_pts – The stencil size.
lbounds – The relative starting point of stencil. For example, a value of (-2,) means that the stencil starts two pixels to the left of where the derivative should be computed.
stencil – The stencil coefficients.
-
DiscreteDerivative(DynCcoord_t nb_pts, DynCcoord_t lbounds, const Eigen::ArrayXd &stencil)¶
Constructor with raw stencil information.
This constructor initializes a DiscreteDerivative object with the provided stencil size, relative starting point of stencil, and stencil coefficients.
- Parameters:
nb_pts – The stencil size.
lbounds – The relative starting point of stencil. For example, a value of (-2,) means that the stencil starts two pixels to the left of where the derivative should be computed.
stencil – The stencil coefficients.
-
DiscreteDerivative(const DiscreteDerivative &other) = default¶
Default copy constructor.
-
DiscreteDerivative(DiscreteDerivative &&other) = default¶
Default move constructor.
-
virtual ~DiscreteDerivative() = default¶
Default destructor.
-
DiscreteDerivative &operator=(const DiscreteDerivative &other) = delete¶
Copy assignment operator is deleted as copying is not allowed for DiscreteDerivative objects
-
DiscreteDerivative &operator=(DiscreteDerivative &&other) = delete¶
Move assignment operator is deleted as moving is not allowed for DiscreteDerivative objects
-
inline Real operator()(const DynCcoord_t &dcoord) const¶
Returns the stencil value at a given coordinate.
This function returns the stencil value at the provided coordinate.
- Parameters:
dcoord – The coordinate at which the stencil value is to be returned.
- Returns:
The stencil value at the provided coordinate.
-
inline const Dim_t &get_dim() const¶
Returns the dimension of the stencil.
This function returns the dimension of the stencil, which is the number of spatial dimensions in which the stencil operates. The stencil is a finite-differences stencil used for computing derivatives.
- Returns:
A constant reference to the Dim_t object that contains the dimension of the stencil.
-
inline const DynCcoord_t &get_nb_pts() const¶
Returns the number of grid points in the stencil.
This function returns the number of grid points in the stencil, which is the size of the stencil. The stencil is a finite-differences stencil used for computing derivatives.
- Returns:
A constant reference to the DynamicCoordinate object that contains the number of grid points in the stencil.
-
inline const DynCcoord_t &get_lbounds() const¶
Returns the lower bounds of the stencil.
This function returns the lower bounds of the stencil, which represent the relative starting point of the stencil. For example, a value of (-2,) means that the stencil starts two pixels to the left of where the derivative should be computed.
- Returns:
A constant reference to the DynamicCoordinate object that contains the lower bounds of the stencil.
-
inline const muGrid::CcoordOps::DynamicPixels &get_pixels() const¶
Returns the pixels class that allows to iterate over pixels.
This function returns the DynamicPixels object from the muGrid::CcoordOps namespace. This object is used to iterate over the pixels of the stencil.
- Returns:
A constant reference to the DynamicPixels object that allows to iterate over the pixels of the stencil.
-
template<typename T>
inline void apply(const muGrid::TypedFieldBase<T> &in_field, Index_t in_dof, muGrid::TypedFieldBase<T> &out_field, Index_t out_dof, Real fac = 1.0) const¶ Apply the “stencil” to a component (degree-of-freedom) of a field and store the result to a select component of a second field.
This function applies the stencil to a component of an input field and stores the result in a selected component of an output field. It performs various checks to ensure the fields are global and the specified degrees of freedom are within range. It then loops over the field pixel iterator and the stencil to compute the derivative and store it in the output field. Note that this function is designed to be inlined by the compiler to optimize loops over degrees of freedom.
Note
This function currently only works without MPI parallelization. If parallelization is needed, apply the stencil in Fourier space using the
fourier
method. Currently, this method is only used in the serial tests.- Template Parameters:
T – The type of the field elements.
- Parameters:
in_field – The input field to which the stencil is applied.
in_dof – The degree of freedom in the input field to which the stencil is applied.
out_field – The output field where the result is stored.
out_dof – The degree of freedom in the output field where the result is stored.
fac – A factor that is multiplied with the derivative. The default value is 1.0.
- Throws:
DerivativeError – If the input or output field is not global, or if the specified degree of freedom is out of range, or if the input and output fields live on incompatible grids.
-
inline virtual Complex fourier(const Vector &phase) const¶
Returns the Fourier representation of the stencil.
Any translationally invariant linear combination of grid values (as expressed through the “stencil”) becomes a multiplication with a number in Fourier space. This method returns the Fourier representation of this stencil.
- Parameters:
phase – The phase is the wavevector times cell dimension, but lacking a factor of 2 π.
- Returns:
The Fourier representation of the stencil.
-
DiscreteDerivative rollaxes(int distance = 1) const¶
Returns a new stencil with rolled axes.
Given a stencil on a three-dimensional grid with axes (x, y, z), the stencil that has been “rolled” by distance one has axes (z, x, y). This is a simple implementation of a rotation operation. For example, given a stencil that described the derivative in the x-direction, rollaxes(1) gives the derivative in the y-direction and rollaxes(2) gives the derivative in the z-direction.
- Parameters:
distance – The distance to roll the axes. Default value is 1.
- Returns:
A new DiscreteDerivative object with rolled axes.
-
inline const Eigen::ArrayXd &get_stencil() const¶
Returns the stencil data.
This function returns the finite-differences stencil data which is used for computing derivatives.
- Returns:
A constant reference to the Eigen::ArrayXd object that contains the stencil data.
-
using Parent = DerivativeBase¶
-
template<Dim_t dim>
class FFT_freqs¶ - #include <fft_utils.hh>
A class that encapsulates the creation and retrieval of wave vectors.
This class is templated on the dimensionality of the wave vectors. It provides methods to get unnormalized, normalized, and complex wave vectors. It also provides methods to get the number of grid points in a specific dimension.
- Template Parameters:
dim – The dimensionality of the wave vectors.
Public Types
-
using CcoordVector = Eigen::Matrix<Dim_t, dim, 1>¶
A type alias for an Eigen matrix with dimensions equivalent to Ccoord_t.
Public Functions
-
FFT_freqs() = delete¶
Default constructor is deleted to prevent creating an object without specifying the grid points.
-
inline explicit FFT_freqs(Ccoord_t<dim> nb_grid_pts)¶
Constructor that initializes the object with the number of grid points.
- Parameters:
nb_grid_pts – The number of grid points in each dimension.
-
inline FFT_freqs(Ccoord_t<dim> nb_grid_pts, std::array<Real, dim> lengths)¶
Constructor that initializes the object with the number of grid points and the lengths of the domain.
- Parameters:
nb_grid_pts – The number of grid points in each dimension.
lengths – The lengths of the domain in each dimension.
-
FFT_freqs(const FFT_freqs &other) = delete¶
Copy constructor is deleted to prevent copying of the object.
-
virtual ~FFT_freqs() = default¶
Destructor.
-
FFT_freqs &operator=(const FFT_freqs &other) = delete¶
Copy assignment operator is deleted to prevent copying of the object.
-
inline Vector get_xi(const Ccoord_t<dim> ccoord) const¶
Get the unnormalized wave vector in sampling units.
Get the xi value for the FFT frequencies.
This function is a member of the FFT_freqs class template. It takes a constant Ccoord_t object as an argument, which represents the coordinates in the frequency domain. The function iterates over each dimension and assigns the corresponding frequency to the return vector.
- Parameters:
ccoord – The coordinates in the frequency domain.
ccoord – A constant reference to a Ccoord_t object representing the coordinates in the frequency domain.
- Returns:
The unnormalized wave vector.
- Template Parameters:
dim – The dimensionality of the FFT operation.
- Returns:
A Vector object where each element is the frequency corresponding to the coordinate in the same dimension.
-
inline VectorComplex get_complex_xi(const Ccoord_t<dim> ccoord) const¶
Get the unnormalized complex wave vector in sampling units.
- Parameters:
ccoord – The coordinates in the frequency domain.
- Returns:
The unnormalized complex wave vector.
-
inline Vector get_unit_xi(const Ccoord_t<dim> ccoord) const¶
Get the normalized wave vector.
- Parameters:
ccoord – The coordinates in the frequency domain.
- Returns:
The normalized wave vector.
-
inline Index_t get_nb_grid_pts(Index_t i) const¶
Get the number of grid points in a specific dimension.
- Parameters:
i – The index of the dimension.
- Returns:
The number of grid points in the specified dimension.
-
class FFTEngineBase¶
- #include <fft_engine_base.hh>
Virtual base class for FFT engines. To be implemented by all FFT_engine implementations.
Subclassed by muFFT::FFTWEngine, muFFT::FFTWMPIEngine, muFFT::PFFTEngine, muFFT::PocketFFTEngine
Public Types
-
using GFieldCollection_t = muGrid::GlobalFieldCollection¶
global FieldCollection
-
using Pixels_t = typename GFieldCollection_t::DynamicPixels¶
pixel iterator
-
using RealField_t = muGrid::TypedFieldBase<Real>¶
Field type on which to apply the projection. This is a TypedFieldBase because it need to be able to hold either TypedField or a WrappedField.
-
using FourierField_t = muGrid::TypedFieldBase<Complex>¶
Field type holding a Fourier-space representation of a real-valued second-order tensor field
-
using iterator = typename GFieldCollection_t::DynamicPixels::iterator¶
iterator over Fourier-space discretisation point
Public Functions
-
FFTEngineBase() = delete¶
Default constructor.
-
FFTEngineBase(const DynCcoord_t &nb_grid_pts, Communicator comm = Communicator(), const FFT_PlanFlags &plan_flags = FFT_PlanFlags::estimate, bool allow_temporary_buffer = true, bool allow_destroy_input = false, bool engine_has_rigid_memory_layout = true)¶
Constructs an FFTEngineBase object with the specified parameters.
This constructor initializes an FFTEngineBase object with the given number of grid points, communicator, FFT planner flags, and buffer options. The constructor does not perform any FFT computations; it merely sets up the object for future FFT operations.
- Parameters:
nb_grid_pts – A DynCcoord_t object representing the number of grid points in each dimension of the global grid.
comm – An optional Communicator object for MPI communication. Defaults to an empty Communicator.
plan_flags – An optional FFT_PlanFlags object representing the FFT planner flags. Defaults to FFT_PlanFlags::estimate.
allow_temporary_buffer – An optional boolean flag indicating whether the creation of temporary buffers is allowed if the input buffer has the wrong memory layout. Defaults to true.
allow_destroy_input – An optional boolean flag indicating whether the input buffers can be invalidated during the FFT. Defaults to false.
engine_has_rigid_memory_layout – An optional boolean flag indicating whether the underlying FFT engine requires a fixed memory layout. Defaults to true.
-
FFTEngineBase(const FFTEngineBase &other) = delete¶
Copy constructor.
-
FFTEngineBase(FFTEngineBase &&other) = delete¶
Move constructor.
-
virtual ~FFTEngineBase() = default¶
Destructor.
-
FFTEngineBase &operator=(const FFTEngineBase &other) = delete¶
Copy assignment operator.
-
FFTEngineBase &operator=(FFTEngineBase &&other) = delete¶
Move assignment operator.
-
virtual void create_plan(const Index_t &nb_dof_per_pixel) = 0¶
prepare a plan for a transform with nb_dof_per_pixel entries per pixel. Needs to be called for every different sized transform
-
void create_plan(const Shape_t &shape = Shape_t{})¶
prepare a plan for a transform with shape entries per pixel. Needs to be called for every different sized transform
-
void fft(const RealField_t &input_field, FourierField_t &output_field)¶
forward transform, performs copy of buffer if required
-
void ifft(const FourierField_t &input_field, RealField_t &output_field)¶
inverse transform, performs copy of buffer if required
-
void hcfft(const RealField_t &input_field, RealField_t &output_field)¶
forward transform using half-complex data storage,
forward transform, performs copy of buffer if required
-
void ihcfft(const RealField_t &input_field, RealField_t &output_field)¶
inverse transform using half-complex data storage,
inverse transform, performs copy of buffer if required
-
virtual muGrid::ComplexField ®ister_fourier_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel)¶
Create a Fourier-space field with the ideal strides and dimensions for this engine. Fields created this way are meant to be reused again and again, and they will stay in the memory of the
muFFT::FFTEngineBase
’s field collection until the engine is destroyed.
-
virtual muGrid::ComplexField ®ister_fourier_space_field(const std::string &unique_name, const Shape_t &shape = Shape_t{})¶
Create a Fourier-space field with the ideal strides and dimensions for this engine. Fields created this way are meant to be reused again and again, and they will stay in the memory of the
muFFT::FFTEngineBase
’s field collection until the engine is destroyed.
-
FourierField_t &fourier_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel)¶
Fetches a Fourier-space field with the ideal strides and dimensions for this engine. If the field does not exist, it is created using
register_fourier_space_field
.
-
FourierField_t &fourier_space_field(const std::string &unique_name, const Shape_t &shape = Shape_t{})¶
Fetches a Fourier-space field with the ideal strides and dimensions for this engine. If the field does not exist, it is created using
register_fourier_space_field
.
-
virtual RealField_t ®ister_halfcomplex_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel)¶
Create a Fourier-space field with the ideal strides and dimensions for this engine. Fields created this way are meant to be reused again and again, and they will stay in the memory of the
muFFT::FFTEngineBase
’s field collection until the engine is destroyed.
-
virtual RealField_t ®ister_halfcomplex_field(const std::string &unique_name, const Shape_t &shape = Shape_t{})¶
Create a Fourier-space field with the ideal strides and dimensions for this engine. Fields created this way are meant to be reused again and again, and they will stay in the memory of the
muFFT::FFTEngineBase
’s field collection until the engine is destroyed.
-
RealField_t &halfcomplex_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel)¶
Fetches a Fourier-space field with the ideal strides and dimensions for this engine. If the field does not exist, it is created using
register_fourier_space_field
.
-
RealField_t &halfcomplex_field(const std::string &unique_name, const Shape_t &shape = Shape_t{})¶
Fetches a Fourier-space field with the ideal strides and dimensions for this engine. If the field does not exist, it is created using
register_fourier_space_field
.
-
virtual RealField_t ®ister_real_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel)¶
Create a real-space field with the ideal strides and dimensions for this engine. Fields created this way are meant to be reused again and again, and they will stay in the memory of the
muFFT::FFTEngineBase
’s field collection until the engine is destroyed.
-
virtual RealField_t ®ister_real_space_field(const std::string &unique_name, const Shape_t &shape = Shape_t{})¶
Create a real-space field with the ideal strides and dimensions for this engine. Fields created this way are meant to be reused again and again, and they will stay in the memory of the
muFFT::FFTEngineBase
’s field collection until the engine is destroyed.
-
RealField_t &real_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel)¶
Fetches a real-space field with the ideal strides and dimensions for this engine. If the field does not exist, it is created using
register_real_space_field
.
-
RealField_t &real_space_field(const std::string &unique_name, const Shape_t &shape = Shape_t{})¶
Fetches a real-space field with the ideal strides and dimensions for this engine. If the field does not exist, it is created using
register_real_space_field
.
-
inline virtual bool has_grid_pts() const¶
return whether this engine is active
-
const Pixels_t &get_fourier_pixels() const¶
iterators over only those pixels that exist in frequency space (i.e. about half of all pixels, see rfft)
-
size_t size() const¶
nb of pixels (mostly for debugging)
-
size_t fourier_size() const¶
nb of pixels in Fourier space
-
size_t workspace_size() const¶
nb of pixels in the work space (may contain a padding region)
-
inline const Communicator &get_communicator() const¶
return the communicator object
-
inline const DynCcoord_t &get_nb_subdomain_grid_pts() const¶
returns the process-local number of grid points in each direction of the cell
-
inline const DynCcoord_t &get_nb_domain_grid_pts() const¶
returns the global number of grid points in each direction of the cell
-
inline const DynCcoord_t &get_subdomain_locations() const¶
returns the process-local locations of the cell
-
inline const DynCcoord_t &get_subdomain_strides() const¶
returns the data layout of the process-local grid
-
inline const DynCcoord_t &get_nb_fourier_grid_pts() const¶
returns the process-local number of grid points in each direction of the cell in Fourier space
-
inline const DynCcoord_t &get_fourier_locations() const¶
returns the process-local locations of the cell in Fourier space
-
inline const DynCcoord_t &get_fourier_strides() const¶
returns the data layout of the cell in Fourier space
-
inline GFieldCollection_t &get_real_field_collection()¶
returns the field collection handling fields in real space
-
inline GFieldCollection_t &get_halfcomplex_field_collection()¶
returns the field collection handling fields confirming with
-
inline GFieldCollection_t &get_fourier_field_collection()¶
returns the field collection handling fields in Fourier space
-
inline Real normalisation() const¶
factor by which to multiply projection before inverse transform (this is typically 1/nb_pixels for so-called unnormalized transforms (see, e.g. http://www.fftw.org/fftw3_doc/Multi_002dDimensional-DFTs-of-Real-Data.html#Multi_002dDimensional-DFTs-of-Real-Data or https://docs.scipy.org/doc/numpy-1.13.0/reference/routines.fft.html . Rather than scaling the inverse transform (which would cost one more loop), FFT engines provide this value so it can be used in the projection operator (where no additional loop is required)
-
const Index_t &get_spatial_dim() const¶
return the number of spatial dimensions
-
inline bool get_allow_temporary_buffer() const¶
return whether temporary buffers are allowed
-
virtual std::unique_ptr<FFTEngineBase> clone() const = 0¶
perform a deep copy of the engine (this should never be necessary in c++)
-
bool has_plan_for(const Index_t &nb_dof_per_pixel) const¶
check whether a plan for nb_dof_per_pixel exists
Protected Functions
-
void initialise_field_collections()¶
calls initialize of the real, hc and fourier field collections
-
virtual void compute_fft(const RealField_t &input_field, FourierField_t &output_field) = 0¶
forward transform, assumes that the buffer has the correct memory layout
-
virtual void compute_ifft(const FourierField_t &input_field, RealField_t &output_field) = 0¶
inverse transform, assumes that the buffer has the correct memory layout
-
virtual void compute_hcfft(const RealField_t &input_field, RealField_t &output_field)¶
forward half complex transform
-
virtual void compute_ihcfft(const RealField_t &input_field, RealField_t &output_field)¶
inverse half complex transform
-
virtual bool check_real_space_field(const RealField_t &field, FFTDirection direction) const¶
check whether real-space buffer has the correct memory layout
-
virtual bool check_fourier_space_field(const FourierField_t &field, FFTDirection direction) const¶
check whether Fourier-space buffer has the correct memory layout
-
virtual bool check_halfcomplex_field(const RealField_t &field, FFTDirection direction) const¶
check whether the half-complex buffer has the correct memory layout
Protected Attributes
-
Index_t spatial_dimension¶
spatial dimension of the grid
-
Communicator comm¶
Field collection in which to store fields associated with Fourier-space points communicator
-
GFieldCollection_t real_field_collection¶
Field collection for real-space fields.
-
GFieldCollection_t fourier_field_collection¶
Field collection for Fourier-space fields.
-
GFieldCollection_t halfcomplex_field_collection¶
Field collection for half-complex-space fields in the r2hc transform real fields and fourier fields are identical In serial the hc_field is identical to the real field, But in parallel, the hc_field has no padding region.
-
const DynCcoord_t nb_domain_grid_pts¶
nb_grid_pts of the full domain of the cell
-
DynCcoord_t nb_subdomain_grid_pts¶
nb_grid_pts of the process-local (subdomain) portion of the cell
-
DynCcoord_t subdomain_locations¶
location of the process-local (subdomain) portion of the cell
-
DynCcoord_t subdomain_strides¶
data layout of the porcess-local portion of the cell
-
DynCcoord_t nb_fourier_grid_pts¶
nb_grid_pts of the process-local (subdomain) portion of the Fourier transformed data
-
DynCcoord_t fourier_locations¶
location of the process-local (subdomain) portion of the Fourier transformed data
-
DynCcoord_t fourier_strides¶
data layout of the process-local (subdomain) portion of the Fourier transformed data
-
bool allow_temporary_buffer¶
allow the FFTEngine to create temporary copies (if it cannot work with a specific memory layout)
-
bool allow_destroy_input¶
allow the FFTEngine to destroy input buffers
-
bool engine_has_rigid_memory_layout¶
the underlying FFT engine requires a fixed memory layout
-
const Real norm_factor¶
normalisation coefficient of fourier transform
-
const FFT_PlanFlags plan_flags¶
FFT planner flags.
-
std::set<Index_t> planned_nb_dofs = {}¶
number of degrees of freedom per pixel for which this field collection has been primed. Can be queried. Corresponds to the number of sub-points per pixel multiplied by the number of components per sub-point
-
using GFieldCollection_t = muGrid::GlobalFieldCollection¶
-
class FFTEngineError : public RuntimeError¶
- #include <fft_engine_base.hh>
base class for FFTEngine-related exceptions
-
class FFTWEngine : public muFFT::FFTEngineBase¶
- #include <fftw_engine.hh>
implements the
muFFT::FftEngine_Base
interface using the FFTW libraryPublic Types
-
using Parent = FFTEngineBase¶
base class
Public Functions
-
FFTWEngine() = delete¶
Default constructor.
-
FFTWEngine(const DynCcoord_t &nb_grid_pts, Communicator comm = Communicator(), const FFT_PlanFlags &plan_flags = FFT_PlanFlags::estimate, bool allow_temporary_buffer = true, bool allow_destroy_input = false)¶
Constructs a FFTWEngine object with the specified parameters.
This constructor initializes a FFTWEngine object with the given number of grid points, communicator, FFT planner flags, and buffer options. The constructor does not perform any FFT computations; it merely sets up the object for future FFT operations.
- Parameters:
nb_grid_pts – A DynCcoord_t object representing the number of grid points in each dimension of the global grid.
comm – An optional Communicator object for MPI communication. Defaults to an empty Communicator.
plan_flags – An optional FFT_PlanFlags object representing the FFT planner flags. Defaults to FFT_PlanFlags::estimate.
allow_temporary_buffer – An optional boolean flag indicating whether the creation of temporary buffers is allowed if the input buffer has the wrong memory layout. Defaults to true.
allow_destroy_input – An optional boolean flag indicating whether the input buffers can be invalidated during the FFT. Defaults to false.
-
FFTWEngine(const DynCcoord_t &nb_grid_pts, const FFT_PlanFlags &plan_flags, bool allow_temporary_buffer = true, bool allow_destroy_input = false)¶
Constructs a FFTWEngine object with the specified parameters.
This constructor initializes a FFTWEngine object with the given number of grid points and FFT planner flags. The constructor does not perform any FFT computations; it merely sets up the object for future FFT operations.
- Parameters:
nb_grid_pts – A DynCcoord_t object representing the number of grid points in each dimension of the global grid.
plan_flags – An FFT_PlanFlags object representing the FFT planner flags.
allow_temporary_buffer – A boolean flag indicating whether the creation of temporary buffers is allowed if the input buffer has the wrong memory layout. Defaults to true.
allow_destroy_input – A boolean flag indicating whether the input buffers can be invalidated during the FFT. Defaults to false.
-
FFTWEngine(const FFTWEngine &other) = delete¶
Copy constructor.
-
FFTWEngine(FFTWEngine &&other) = delete¶
Move constructor.
-
virtual ~FFTWEngine() noexcept¶
Destructor.
-
FFTWEngine &operator=(const FFTWEngine &other) = delete¶
Copy assignment operator.
-
FFTWEngine &operator=(FFTWEngine &&other) = delete¶
Move assignment operator.
-
virtual void create_plan(const Index_t &nb_dof_per_pixel) override¶
prepare a plan for a transform with nb_dof_per_pixel entries per pixel. Needs to be called for every different sized transform
-
virtual std::unique_ptr<FFTEngineBase> clone() const final¶
perform a deep copy of the engine (this should never be necessary in c++)
Protected Functions
-
virtual void compute_fft(const RealField_t &input_field, FourierField_t &output_field) override¶
forward transform
-
virtual void compute_ifft(const FourierField_t &input_field, RealField_t &output_field) override¶
inverse transform
-
virtual void compute_hcfft(const RealField_t &input_field, RealField_t &output_field) override¶
forward half complex transform
-
virtual void compute_ihcfft(const RealField_t &input_field, RealField_t &output_field) override¶
inverse half complex transform
-
virtual bool check_fourier_space_field(const FourierField_t &field, FFTDirection direction) const final¶
check whether Fourier-space buffer has the correct memory layout
Protected Attributes
-
std::map<Index_t, fftw_plan> fft_plans = {}¶
holds the plans for forward fourier transforms
-
std::map<Index_t, fftw_plan> ifft_plans = {}¶
holds the plans for inversefourier transforms
-
std::map<Index_t, fftw_plan> hcfft_plans = {}¶
holds the plans for forward half-complex fourier transforms
-
std::map<Index_t, fftw_plan> ihcfft_plans = {}¶
holds the plans for inverse half-complex fourier transforms
Protected Static Attributes
-
static Index_t nb_engines = {0}¶
global counter for total number of initialized engines
-
using Parent = FFTEngineBase¶
-
class FFTWMPIEngine : public muFFT::FFTEngineBase¶
- #include <fftwmpi_engine.hh>
implements the
muFFT::FFTEngineBase
interface using the FFTW libraryPublic Types
-
using Parent = FFTEngineBase¶
base class
Public Functions
-
FFTWMPIEngine() = delete¶
Default constructor.
-
FFTWMPIEngine(const DynCcoord_t &nb_grid_pts, Communicator comm = Communicator(), const FFT_PlanFlags &plan_flags = FFT_PlanFlags::estimate, bool allow_temporary_buffer = true, bool allow_destroy_input = false)¶
Constructs a FFTWMPIEngine object with the specified parameters.
This constructor initializes a FFTWMPIEngine object with the given number of grid points, communicator, FFT planner flags, and buffer options. The constructor does not perform any FFT computations; it merely sets up the object for future FFT operations.
- Parameters:
nb_grid_pts – A DynCcoord_t object representing the number of grid points in each dimension of the global grid.
comm – An optional Communicator object for MPI communication. Defaults to an empty Communicator.
plan_flags – An optional FFT_PlanFlags object representing the FFT planner flags. Defaults to FFT_PlanFlags::estimate.
allow_temporary_buffer – An optional boolean flag indicating whether the creation of temporary buffers is allowed if the input buffer has the wrong memory layout. Defaults to true.
allow_destroy_input – An optional boolean flag indicating whether the input buffers can be invalidated during the FFT. Defaults to false.
-
FFTWMPIEngine(const FFTWMPIEngine &other) = delete¶
Copy constructor.
-
FFTWMPIEngine(FFTWMPIEngine &&other) = delete¶
Move constructor.
-
virtual ~FFTWMPIEngine() noexcept¶
Destructor.
-
FFTWMPIEngine &operator=(const FFTWMPIEngine &other) = delete¶
Copy assignment operator.
-
FFTWMPIEngine &operator=(FFTWMPIEngine &&other) = delete¶
Move assignment operator.
-
virtual void create_plan(const Index_t &nb_dof_per_pixel) override¶
prepare a plan for a transform with nb_dof_per_pixel entries per pixel. Needs to be called for every different sized transform
-
inline virtual bool has_grid_pts() const override¶
return whether this engine is active (an engine is active if it has more than zero grid points. FFTWMPI sometimes assigns zero grid points)
-
virtual std::unique_ptr<FFTEngineBase> clone() const final¶
perform a deep copy of the engine (this should never be necessary in c++)
-
virtual RealField_t ®ister_real_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel) final¶
need to override this method here, since FFTWMPI requires field padding
-
virtual RealField_t ®ister_real_space_field(const std::string &unique_name, const Shape_t &shape) final¶
need to override this method here, since FFTWMPI requires field padding
-
virtual muGrid::ComplexField ®ister_fourier_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel) final¶
need to override this method here, since FFTWMPI requires field padding
-
virtual muGrid::ComplexField ®ister_fourier_space_field(const std::string &unique_name, const Shape_t &shape) final¶
need to override this method here, since FFTWMPI requires field padding
Protected Functions
-
virtual void compute_fft(const RealField_t &field, FourierField_t &output_field) override¶
forward transform
-
virtual void compute_ifft(const FourierField_t &input_field, RealField_t &output_field) override¶
inverse transform
-
virtual bool check_real_space_field(const RealField_t &field, FFTDirection direction) const final¶
check whether real-space buffer has the correct memory layout
-
virtual bool check_fourier_space_field(const FourierField_t &field, FFTDirection direction) const final¶
check whether Fourier-space buffer has the correct memory layout
Protected Attributes
-
std::map<Index_t, fftw_plan> fft_plans = {}¶
holds the plans for forward fourier transforms
-
std::map<Index_t, fftw_plan> ifft_plans = {}¶
holds the plans for inversefourier transforms
-
std::map<Index_t, Index_t> required_workspace_sizes = {}¶
holds the fourier field sizes including padding for different transforms
-
bool active = {true}¶
FFTWMPI sometimes assigns zero grid points
-
std::vector<ptrdiff_t> nb_fourier_non_transposed = {}¶
Input to local_size_many_transposed.
Protected Static Attributes
-
static int nb_engines = {0}¶
number of times this engine has been instantiated
-
using Parent = FFTEngineBase¶
-
class FourierDerivative : public muFFT::DerivativeBase¶
- #include <derivative.hh>
A class that represents a derivative computed by Fourier interpolation.
This class is derived from the DerivativeBase class and provides functionalities for Fourier interpolated derivatives.
Public Types
-
using Parent = DerivativeBase¶
Alias for the base class.
Public Functions
-
FourierDerivative() = delete¶
Deleted default constructor.
This constructor is deleted because a FourierDerivative object requires a spatial dimension and direction to be properly initialized.
-
explicit FourierDerivative(Index_t spatial_dimension, Index_t direction)¶
Constructor that takes the spatial dimension and direction as arguments.
- Parameters:
spatial_dimension – The spatial dimension of the derivative.
direction – The direction of the derivative.
-
explicit FourierDerivative(Index_t spatial_dimension, Index_t direction, const Eigen::ArrayXd &shift)¶
Constructor that takes the spatial dimension, direction, and shift info as arguments.
- Parameters:
spatial_dimension – The spatial dimension of the derivative.
direction – The direction of the derivative.
shift – The shift information for the derivative.
-
FourierDerivative(const FourierDerivative &other) = default¶
Default copy constructor.
-
FourierDerivative(FourierDerivative &&other) = default¶
Default move constructor.
-
virtual ~FourierDerivative() = default¶
Default destructor.
-
FourierDerivative &operator=(const FourierDerivative &other) = delete¶
Deleted copy assignment operator.
-
FourierDerivative &operator=(FourierDerivative &&other) = delete¶
Deleted move assignment operator.
-
inline virtual Complex fourier(const Vector &phase) const¶
Returns the Fourier derivative.
Returns the derivative of the Fourier interpolated field.
- Parameters:
phase – The phase is the wavevector times cell dimension, but lacking a factor of 2 π.
- Returns:
The Fourier representation of the derivative.
-
using Parent = DerivativeBase¶
-
class PFFTEngine : public muFFT::FFTEngineBase¶
- #include <pfft_engine.hh>
implements the
muFFT::FFTEngineBase
interface using the FFTW libraryPublic Types
-
using Parent = FFTEngineBase¶
base class
Public Functions
-
PFFTEngine() = delete¶
Default constructor.
-
PFFTEngine(const DynCcoord_t &nb_grid_pts, Communicator comm = Communicator(), const FFT_PlanFlags &plan_flags = FFT_PlanFlags::estimate, bool allow_temporary_buffer = true, bool allow_destroy_input = false)¶
Constructor with the domain’s number of grid points in each direction and the communicator
- Parameters:
nb_grid_pts – number of grid points of the global grid
comm – MPI communicator object
plan_flags – MPI planner flags
allow_temporary_buffer – allow the creation of temporary buffers if the input buffer has the wrong memory layout
allow_destroy_input – allow that the input buffers are invalidated during the FFT
-
PFFTEngine(const PFFTEngine &other) = delete¶
Copy constructor.
-
PFFTEngine(PFFTEngine &&other) = delete¶
Move constructor.
-
virtual ~PFFTEngine() noexcept¶
Destructor.
-
PFFTEngine &operator=(const PFFTEngine &other) = delete¶
Copy assignment operator.
-
PFFTEngine &operator=(PFFTEngine &&other) = delete¶
Move assignment operator.
-
virtual void create_plan(const Index_t &nb_dof_per_pixel) override¶
prepare a plan for a transform with nb_dof_per_pixel entries per pixel. Needs to be called for every different sized transform
-
virtual std::unique_ptr<FFTEngineBase> clone() const final¶
perform a deep copy of the engine (this should never be necessary in c++)
-
virtual RealField_t ®ister_real_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel) final¶
need to override this method here, since FFTWMPI requires field padding
-
virtual RealField_t ®ister_real_space_field(const std::string &unique_name, const Shape_t &shape) final¶
need to override this method here, since FFTWMPI requires field padding
-
virtual muGrid::ComplexField ®ister_fourier_space_field(const std::string &unique_name, const Index_t &nb_dof_per_pixel) final¶
need to override this method here, since FFTWMPI requires field padding
-
virtual muGrid::ComplexField ®ister_fourier_space_field(const std::string &unique_name, const Shape_t &shape) final¶
need to override this method here, since FFTWMPI requires field padding
Protected Functions
-
virtual void compute_fft(const RealField_t &field, FourierField_t &output_field) override¶
forward transform
-
virtual void compute_ifft(const FourierField_t &input_field, RealField_t &output_field) override¶
inverse transform
-
virtual bool check_real_space_field(const RealField_t &field, FFTDirection direction) const final¶
check whether real-space buffer has the correct memory layout
-
virtual bool check_fourier_space_field(const FourierField_t &field, FFTDirection direction) const final¶
check whether Fourier-space buffer has the correct memory layout
Protected Attributes
-
MPI_Comm mpi_comm¶
MPI communicator
-
std::map<Index_t, pfft_plan> fft_plans = {}¶
holds the plans for forward fourier transforms
-
std::map<Index_t, pfft_plan> ifft_plans = {}¶
holds the plans for inversefourier transforms
-
std::map<Index_t, Index_t> required_workspace_sizes = {}¶
holds the fourier field sizes including padding for different transforms
-
bool active = {true}¶
PFFT sometimes assigns zero grid points.
Protected Static Attributes
-
static int nb_engines = {0}¶
number of times this engine has been instantiated
-
using Parent = FFTEngineBase¶
-
class PocketFFTEngine : public muFFT::FFTEngineBase¶
- #include <pocketfft_engine.hh>
implements the
muFFT::FftEngine_Base
interface using PocketFFT (that is shipped as part of this code). See pocketfft/LICENSE.md for PocketFFT’s license.Public Types
-
using Parent = FFTEngineBase¶
base class
Public Functions
-
PocketFFTEngine() = delete¶
Default constructor.
-
PocketFFTEngine(const DynCcoord_t &nb_grid_pts, Communicator comm = Communicator(), const FFT_PlanFlags &plan_flags = FFT_PlanFlags::estimate, bool allow_temporary_buffer = true, bool allow_destroy_input = false)¶
Constructs a PocketFFTEngine object with the specified parameters.
This constructor initializes a PocketFFTEngine object with the given number of grid points, communicator, FFT planner flags, and buffer options. The constructor does not perform any FFT computations; it merely sets up the object for future FFT operations.
- Parameters:
nb_grid_pts – A DynCcoord_t object representing the number of grid points in each dimension of the global grid.
comm – An optional Communicator object for MPI communication. Defaults to an empty Communicator.
plan_flags – An optional FFT_PlanFlags object representing the FFT planner flags. Defaults to FFT_PlanFlags::estimate.
allow_temporary_buffer – An optional boolean flag indicating whether the creation of temporary buffers is allowed if the input buffer has the wrong memory layout. Defaults to true.
allow_destroy_input – An optional boolean flag indicating whether the input buffers can be invalidated during the FFT. Defaults to false.
-
PocketFFTEngine(const DynCcoord_t &nb_grid_pts, const FFT_PlanFlags &plan_flags, bool allow_temporary_buffer = true, bool allow_destroy_input = false)¶
Constructs a PocketFFTEngine object with the specified parameters.
This constructor initializes a PocketFFTEngine object with the given number of grid points and FFT planner flags. The constructor does not perform any FFT computations; it merely sets up the object for future FFT operations.
- Parameters:
nb_grid_pts – A DynCcoord_t object representing the number of grid points in each dimension of the global grid.
plan_flags – An FFT_PlanFlags object representing the FFT planner flags.
allow_temporary_buffer – A boolean flag indicating whether the creation of temporary buffers is allowed if the input buffer has the wrong memory layout. Defaults to true.
allow_destroy_input – A boolean flag indicating whether the input buffers can be invalidated during the FFT. Defaults to false.
-
PocketFFTEngine(const PocketFFTEngine &other) = delete¶
Copy constructor.
-
PocketFFTEngine(PocketFFTEngine &&other) = delete¶
Move constructor.
-
virtual ~PocketFFTEngine() noexcept¶
Destructor.
-
PocketFFTEngine &operator=(const PocketFFTEngine &other) = delete¶
Copy assignment operator.
-
PocketFFTEngine &operator=(PocketFFTEngine &&other) = delete¶
Move assignment operator.
-
virtual void create_plan(const Index_t &nb_dof_per_pixel) override¶
prepare a plan for a transform with nb_dof_per_pixel entries per pixel. Needs to be called for every different sized transform
-
virtual std::unique_ptr<FFTEngineBase> clone() const final¶
perform a deep copy of the engine (this should never be necessary in c++)
Protected Functions
-
virtual void compute_fft(const RealField_t &input_field, FourierField_t &output_field) override¶
forward transform
-
virtual void compute_ifft(const FourierField_t &input_field, RealField_t &output_field) override¶
inverse transform
-
using Parent = FFTEngineBase¶
-
namespace muFFT¶
Typedefs
-
using Gradient_t = std::vector<std::shared_ptr<DerivativeBase>>¶
convenience alias
-
using FFTEngine_ptr = std::shared_ptr<FFTEngineBase>¶
reference to fft engine is safely managed through a
std::shared_ptr
Enums
-
enum class FFT_PlanFlags¶
Planner flags for FFT.
This enumeration follows the FFTW library’s convention for planning flags. The hope is that this choice will be compatible with alternative FFT implementations.
Values:
-
enumerator estimate¶
Represents the cheapest plan for the slowest execution.
This flag is used when the priority is to minimize the planning cost, even if it results in slower FFT execution.
-
enumerator measure¶
Represents a more expensive plan for faster execution.
This flag is used when the priority is to balance the planning cost and FFT execution speed.
-
enumerator patient¶
Represents the most expensive plan for the fastest execution.
This flag is used when the priority is to maximize the FFT execution speed, even if it results in higher planning cost.
-
enumerator estimate¶
-
enum class FFTDirection¶
Represents the direction of FFT transformation.
This enum class is used to define the direction of the Fast Fourier Transform (FFT) operation. It defines two possible directions: Forward and Reverse.
Values:
-
enumerator forward¶
Represents the forward direction of FFT transformation.
This value is used when the FFT operation is to be performed in the forward direction.
-
enumerator reverse¶
Represents the reverse direction of FFT transformation.
This value is used when the FFT operation is to be performed in the reverse direction.
-
enumerator forward¶
Functions
-
std::ostream &operator<<(std::ostream &os, const DiscreteDerivative &derivative)¶
Allows inserting
muFFT::DiscreteDerivative
s intostd::ostream
s
-
Gradient_t make_fourier_gradient(const Index_t &spatial_dimension)¶
convenience function to build a spatial_dimension-al gradient operator using exact Fourier differentiation
- Parameters:
spatial_dimension – number of spatial dimensions
-
std::valarray<Real> fft_freqs(size_t nb_samples)¶
Compute nondimensional FFT frequencies.
This function computes the FFT frequencies for a given index and number of samples. The frequency is normalized to the number of samples.
- Parameters:
nb_samples – The number of samples.
- Returns:
A valarray containing the FFT frequencies for the given number of samples.
-
std::valarray<Real> fft_freqs(size_t nb_samples, Real length)¶
Compute FFT frequencies in correct length or time units.
This function computes the FFT frequencies for a given number of samples, taking into account the correct length or time units. The length refers to the total size of the domain over which the FFT is taken (for instance the length of an edge of an RVE).
- Parameters:
nb_samples – The number of samples.
length – The length of the domain over which the FFT is taken.
- Returns:
A valarray containing the FFT frequencies for the given number of samples, in the correct length or time units.
-
template<size_t dim>
constexpr Ccoord_t<dim> get_nb_hermitian_grid_pts(Ccoord_t<dim> full_nb_grid_pts)¶ Returns the hermitian grid corresponding to a full grid.
This function template returns the hermitian grid that corresponds to a real-valued field on a full grid, assuming that the last dimension is not fully represented in reciprocal space because of symmetry.
- Template Parameters:
dim – The dimensionality of the grid.
- Parameters:
full_nb_grid_pts – A Ccoord_t object representing the number of grid points in the full grid.
- Returns:
A Ccoord_t object representing the hermitian grid.
-
template<size_t MaxDim>
inline muGrid::DynCcoord<MaxDim> get_nb_hermitian_grid_pts(muGrid::DynCcoord<MaxDim> full_nb_grid_pts)¶ Returns the hermitian grid corresponding to a full grid.
This function template returns the hermitian grid that corresponds to a real-valued field on a full grid, assuming that the last dimension is not fully represented in reciprocal space because of symmetry. It supports dynamic dimensionality.
- Template Parameters:
MaxDim – The maximum dimensionality of the grid.
- Parameters:
full_nb_grid_pts – A muGrid::DynCcoord object representing the number of grid points in the full grid.
- Throws:
RuntimeError – if the dimensionality of the grid is not 1, 2, or 3.
- Returns:
A muGrid::DynCcoord object representing the hermitian grid.
-
inline Int fft_freq(Int i, size_t nb_samples)¶
Compute nondimensional FFT frequency.
This function computes the FFT frequency for a given index and number of samples. The frequency is normalized to the number of samples.
- Parameters:
i – The index for which to compute the FFT frequency.
nb_samples – The number of samples.
- Returns:
The FFT frequency for the given index and number of samples.
-
inline Real fft_freq(Int i, size_t nb_samples, Real length)¶
Compute FFT frequency in correct length or time units.
This function computes the FFT frequency for a given index and number of samples, taking into account the correct length or time units. The length refers to the total size of the domain over which the FFT is taken (for instance the length of an edge of an RVE).
- Parameters:
i – The index for which to compute the FFT frequency.
nb_samples – The number of samples.
length – The length of the domain over which the FFT is taken.
- Returns:
The FFT frequency for the given index and number of samples, in the correct length or time units.
-
template<size_t dim>
inline std::array<std::valarray<Real>, dim> fft_freqs(Ccoord_t<dim> nb_grid_pts)¶ Compute multidimensional nondimensional FFT frequencies.
This function computes the FFT frequencies for a given index and number of samples for multidimensional fields. The frequency is normalized to the number of samples.
- Template Parameters:
dim – The dimensionality of the grid.
- Parameters:
nb_grid_pts – A Ccoord_t object representing the number of grid points in each dimension.
- Returns:
An array of valarrays where each valarray contains the FFT frequencies for a specific dimension.
-
template<size_t dim>
inline std::array<std::valarray<Real>, dim> fft_freqs(Ccoord_t<dim> nb_grid_pts, std::array<Real, dim> lengths)¶ Compute multidimensional FFT frequencies for a grid in correct length or time units.
This function template computes the FFT frequencies for a grid, taking into account the correct length or time units. It iterates over each dimension of the grid, computes the FFT frequencies for that dimension, and stores them in a return array.
- Template Parameters:
dim – The dimensionality of the grid.
- Parameters:
nb_grid_pts – A Ccoord_t object representing the number of grid points in each dimension.
lengths – An array representing the lengths of the domain in each dimension.
- Returns:
An array of valarrays where each valarray contains the FFT frequencies for a specific dimension.
-
Index_t _get_offset(Index_t index, Shape_t shape, Shape_t strides)¶
-
using Gradient_t = std::vector<std::shared_ptr<DerivativeBase>>¶
-
namespace internal¶
-
namespace version¶
Functions
-
std::string info()¶
Returns a formatted text that can be printed to stdout or to output files.
This function returns a formatted text that contains the git commit hash and repository url used to compile µFFT and whether the current state was dirty or not.
- Returns:
A string that contains the formatted text.
-
const char *hash()¶
Returns the git commit hash.
This function returns the git commit hash used to compile µFFT.
- Returns:
A constant character pointer that points to the git commit hash.
-
const char *description()¶
Returns the version string.
This function returns the formatted string of the version of µFFT.
- Returns:
A constant character pointer that points to the repository description.
-
bool is_dirty()¶
Checks if the current state was dirty or not.
This function checks if the current state of the git repository was dirty or not when compiling µFFT.
- Returns:
A boolean value that indicates whether the current state was dirty or not.
-
std::string info()¶
-
namespace version¶
A namespace that contains functions related to version information.
- file derivative.cc
- #include <iostream>#include “derivative.hh”
Representation of finite-differences stencils.
- Author
Richard Leute richard.leute@imtek.uni-freiburg.de Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
05 June 2019
Copyright © 2019 Lars Pastewka
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file derivative.hh
- #include <memory>#include “libmugrid/ccoord_operations.hh”#include “libmugrid/exception.hh”#include “libmugrid/field_map.hh”#include “libmugrid/field_typed.hh”#include “libmugrid/field_collection_global.hh”#include “mufft_common.hh”
Representation of finite-differences stencils.
- Author
Richard Leute richard.leute@imtek.uni-freiburg.de Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
05 June 2019
Copyright © 2019 Lars Pastewka
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fft_engine_base.cc
- #include “fft_engine_base.hh”#include “fft_utils.hh”#include “libmugrid/ccoord_operations.hh”
implementation for FFT engine base class
- Author
Till Junge till.junge@altermail.ch
- Date
03 Dec 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fft_engine_base.hh
- #include “libmugrid/ccoord_operations.hh”#include “libmugrid/field_collection_global.hh”#include “libmugrid/field_typed.hh”#include “libmugrid/communicator.hh”#include “mufft_common.hh”#include <set>
Interface for FFT engines.
- Author
Till Junge till.junge@epfl.ch
- Date
01 Dec 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fft_utils.cc
- #include “fft_utils.hh”
implementation of fft utilities
- Author
Till Junge till.junge@altermail.ch
- Date
11 Dec 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fft_utils.hh
- #include “libmugrid/exception.hh”#include “mufft_common.hh”#include “Eigen/Dense”#include <array>#include <valarray>
collection of functions used in the context of spectral operations
- Author
Till Junge till.junge@epfl.ch
- Date
06 Dec 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fftw_engine.cc
- #include <sstream>#include “libmugrid/ccoord_operations.hh”#include “libmugrid/exception.hh”#include “fftw_engine.hh”
implements the fftw engine
- Author
Till Junge till.junge@altermail.ch
- Date
03 Dec 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fftw_engine.hh
- #include “fft_engine_base.hh”#include <fftw3.h>
FFT engine using FFTW.
- Author
Till Junge till.junge@altermail.ch
- Date
03 Dec 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fftwmpi_engine.cc
- #include “libmugrid/ccoord_operations.hh”#include “libmugrid/exception.hh”#include “fftwmpi_engine.hh”#include “fft_utils.hh”
implements the MPI-parallel fftw engine
- Author
Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
06 Mar 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file fftwmpi_engine.hh
- #include “fft_engine_base.hh”#include <fftw3-mpi.h>
FFT engine using MPI-parallel FFTW.
- Author
Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
06 Mar 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file mufft_common.hh
- #include “libmugrid/grid_common.hh”
Small definitions of commonly used types throughout µFFT.
- Author
Till Junge till.junge@epfl.ch
- Date
24 Jan 2019
Copyright © 2019 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file pfft_engine.cc
- #include “libmugrid/ccoord_operations.hh”#include “libmugrid/exception.hh”#include “pfft_engine.hh”
implements the MPI-parallel pfft engine
- Author
Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
06 Mar 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file pfft_engine.hh
- #include “fft_engine_base.hh”#include <pfft.h>
FFT engine using MPI-parallel PFFT.
- Author
Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
06 Mar 2017
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file pocketfft_engine.cc
- #include <iostream>#include <sstream>#include “pocketfft_hdronly.h”#include “libmugrid/ccoord_operations.hh”#include “pocketfft_engine.hh”
implements the PocketFFT engine
- Author
Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
20 Nov 2022
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
Defines
-
POCKETFFT_NO_MULTITHREADING¶
- file pocketfft_engine.hh
- #include “fft_engine_base.hh”
FFT engine using PocketFFT.
- Author
Lars Pastewka lars.pastewka@imtek.uni-freiburg.de
- Date
20 Nov 2022
Copyright © 2017 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330,
Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- file version.cc
- #include <string>#include <sstream>#include “libmugrid/grid_common.hh”
File written at compile time containing git state.
- Author
Till Junge till.junge@epfl.ch
- Date
04 Feb 2020
Copyright © 2020 Till Junge
µFFT is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3, or (at your option) any later version.
µFFT is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with µFFT; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
Additional permission under GNU GPL version 3 section 7
If you modify this Program, or any covered work, by linking or combining it with proprietary FFT implementations or numerical libraries, containing parts covered by the terms of those libraries’ licenses, the licensors of this Program grant you additional permission to convey the resulting work.
- dir /home/runner/work/muFFT/muFFT/src/libmufft
- dir /home/runner/work/muFFT/muFFT/src