HILA
Loading...
Searching...
No Matches
Field< T > Class Template Reference

The field class implements the standard methods for accessing Fields. Hilapp replaces the parity access patterns, Field[par] with a loop over the appropriate sites. Since the Field class is one of the more important functionalities of HILA, extensive general instructions on the Field class can be read at HILA Functionality documentation. More...

#include <field.h>

Public Member Functions

void check_alloc ()
 Allocate Field if it is not already allocated.
 
void clear ()
 Destroys field data.
 
Field< T > conj () const
 Returns field with all elements conjugated depending how conjugate is defined for type.
 
template<typename A >
void copy_boundary_condition (const Field< A > &rhs)
 Copy the boundary condition from another field.
 
template<typename R = T, typename A = decltype(::dagger(std::declval<R>()))>
Field< A > dagger () const
 Returns dagger or Hermitian conjugate of Field depending on how it is defined for Field type T.
 
Field< T > FFT (const CoordinateVector &dirs, fft_direction fdir=fft_direction::forward) const
 Field method for performing FFT.
 
Field< Complex< hila::arithmetic_type< T > > > FFT_real_to_complex (fft_direction fdir=fft_direction::forward) const
 
 Field ()
 Field constructor.
 
hila::bc get_boundary_condition (Direction dir) const
 Get the boundary condition of the Field.
 
const T get_element (const CoordinateVector &coord) const
 Get singular element which will be broadcast to all nodes.
 
std::vector< T > get_elements (const std::vector< CoordinateVector > &coord_list, bool broadcast=false) const
 Get a list of elements and store them into a vector.
 
std::vector< T > get_slice (const CoordinateVector &c, bool broadcast=false) const
 Get a slice (subvolume)
 
std::vector< T > get_subvolume (const CoordinateVector &cmin, const CoordinateVector &cmax, bool broadcast=false) const
 Get a subvolume of the field elements to all nodes.
 
auto get_value_at (const unsigned i) const
 Get an individual element outside a loop. This is also used as a getter in the vanilla code.
 
gpu_minmax (bool min_or_max, Parity par, CoordinateVector &loc) const
 Declare gpu_reduce here, defined only for GPU targets.
 
template<typename R = T, typename A = decltype(::imag(std::declval<R>()))>
Field< A > imag () const
 Returns imaginary part of Field.
 
bool is_allocated () const
 Returns true if the Field data has been allocated.
 
bool is_initialized (Parity p) const
 Returns true if the Field has been assigned to.
 
minmax (bool is_min, Parity par, CoordinateVector &loc) const
 Function to perform min or max operations.
 
Lattice mylattice () const
 Return the lattice to which this field instance belongs to.
 
double norm ()
 Norm.
 
template<typename A , typename B >
auto operator* (const Field< A > &lhs, const Field< B > &rhs) -> Field< hila::type_mul< A, B > >
 Multiplication operator.
 
template<typename A , std::enable_if_t< std::is_convertible< hila::type_mul< T, A >, T >::value, int > = 0>
Field< T > & operator*= (const Field< A > &rhs)
 Product assignment operator.
 
Field< T > operator+ () const
 Summation operator.
 
template<typename A , std::enable_if_t< std::is_convertible< hila::type_plus< T, A >, T >::value, int > = 0>
Field< T > & operator+= (const Field< A > &rhs)
 Addition assignment operator.
 
Field< T > operator- () const
 Subtraction operator.
 
template<typename A , std::enable_if_t< std::is_convertible< hila::type_minus< T, A >, T >::value, int > = 0>
Field< T > & operator-= (const Field< A > &rhs)
 Subtraction assignment operator.
 
template<typename A , typename B >
auto operator/ (const Field< A > &l, const Field< B > &r) -> Field< hila::type_div< A, B > >
 Division operator.
 
template<typename A , std::enable_if_t< std::is_convertible< hila::type_div< T, A >, T >::value, int > = 0>
Field< T > & operator/= (const Field< A > &rhs)
 Division assignment operator.
 
Field< T > & operator= (const Field< T > &rhs)
 Assignment operator.
 
Field< T > & operator= (Field< T > &&rhs)
 Move Assignment.
 
template<typename S >
bool operator== (const Field< S > &rhs) const
 Field comparison operator.
 
product (Parity par=Parity::all, bool allreduce=true) const
 Product reduction of Field.
 
void read (std::ifstream &inputfile)
 Read the Field from a stream.
 
void read (std::ifstream &inputfile, const CoordinateVector &insize)
 Read the Field from a stream.
 
template<typename R = T, typename A = decltype(::real(std::declval<R>()))>
Field< A > real () const
 Returns real part of Field.
 
Field< T > reflect (const CoordinateVector &dirs) const
 Reflect the Field around the desired axis.
 
void set_boundary_condition (Direction dir, hila::bc bc)
 Set the boundary condition in a given Direction (periodic or antiperiodic)
 
template<typename A , std::enable_if_t< std::is_assignable< T &, A >::value, int > = 0>
const T set_element (const CoordinateVector &coord, const A &value)
 Get singular element which will be broadcast to all nodes.
 
template<typename A >
void set_value_at (const A &value, unsigned i)
 Set an individual element outside a loop. This is also used as a getter in the vanilla code.
 
Field< T > & shift (const CoordinateVector &v, Field< T > &r, Parity par) const
 Create a periodically shifted copy of the field.
 
double squarenorm () const
 Squarenorm.
 
dir_mask_t start_gather (Direction d, Parity p=ALL) const
 
sum (Parity par=Parity::all, bool allreduce=true) const
 Sum reduction of Field.
 
void unblock_to (Field< T > &target) const
 Copy content to the argument Field on blocked (sparse) sites. a.unblock_to(b) is the inverse of a.block_from(b)
 
void write (const std::string &filename, bool binary=true, int precision=8) const
 Write the Field to a named file replacing the file.
 
void write (std::ofstream &outputfile, bool binary=true, int precision=8) const
 Write the field to a file stream.
 
void write_subvolume (std::ofstream &outputfile, const CoordinateVector &cmin, const CoordinateVector &cmax, int precision=6) const
 
Min functions

Minimum value of Field. If CoordinateVector is passed to function then location of minimum value will be stored in said CoordinateVector

Parameters
parParity
Returns
T Minimum value of Field
min (Parity par=ALL) const
 Find minimum value from Field.
 
min (CoordinateVector &loc) const
 Find minimum value and location from Field.
 
min (Parity par, CoordinateVector &loc) const
 Find minimum value and location from Field.
 
Max functions

Maximum value of Field. If CoordinateVector is passed to function then location of maximum value will be stored in said CoordinateVector

Parameters
parParity
Returns
T Minimum value of Field
max (Parity par=ALL) const
 Find maximum value from Field.
 
max (CoordinateVector &loc) const
 Find maximum value and location from Field.
 
max (Parity par, CoordinateVector &loc) const
 Find maximum value and location from Field.
 

Public Attributes

field_struct *__restrict__ fs
 Field::fs holds all of the field content in Field::field_struct.
 

Related Symbols

(Note that these are not member symbols.)

Mathematical methods.
Template Parameters
TField element type of arg
RField return type
Parameters
arginput Field
Returns
Field<R>

See field documentation

template<typename T , typename R = decltype(exp(std::declval<T>()))>
Field< R > exp (const Field< T > &arg)
 Exponential.
 
template<typename T , typename R = decltype(log(std::declval<T>()))>
Field< R > log (const Field< T > &arg)
 Logarithm.
 
template<typename T , typename R = decltype(sin(std::declval<T>()))>
Field< R > sin (const Field< T > &arg)
 Sine.
 
template<typename T , typename R = decltype(cos(std::declval<T>()))>
Field< R > cos (const Field< T > &arg)
 Cosine.
 
template<typename T , typename R = decltype(tan(std::declval<T>()))>
Field< R > tan (const Field< T > &arg)
 Tangent.
 
template<typename T , typename R = decltype(asin(std::declval<T>()))>
Field< R > asin (const Field< T > &arg)
 Arcsine.
 
template<typename T , typename R = decltype(acos(std::declval<T>()))>
Field< R > acos (const Field< T > &arg)
 Arccosine.
 
template<typename T , typename R = decltype(atan(std::declval<T>()))>
Field< R > atan (const Field< T > &arg)
 Arctangent.
 
template<typename T , typename R = decltype(abs(std::declval<T>()))>
Field< R > abs (const Field< T > &arg)
 Absolute value.
 
template<typename T , typename P , typename R = decltype(pow(std::declval<T>()), std::declval<P>())>
Field< R > pow (const Field< T > &arg, const P p)
 Power.
 
template<typename T >
double squarenorm (const Field< T > &arg)
 Squared norm \(|f|^2\).
 
template<typename A , typename B , typename R = decltype(std::declval<A>() - std::declval<B>())>
double squarenorm_relative (const Field< A > &a, const Field< B > &b)
 Squarenorm relative \(|a-b|^2\).
 

Detailed Description

template<typename T>
class Field< T >

The field class implements the standard methods for accessing Fields. Hilapp replaces the parity access patterns, Field[par] with a loop over the appropriate sites. Since the Field class is one of the more important functionalities of HILA, extensive general instructions on the Field class can be read at HILA Functionality documentation.

The Field class contains member functions used by hilapp, which are marked internal, as well as members that are useful for application developers.

The Field mainly implements the interface to the Field and not the content.

The Field contains a pointer to Field::field_struct, which implements MPI communication of the Field boundaries. Though generally the application developer does not need to worry about the Field::field_struct class, hence it is marked as internal and it's documentation cannot be viewed in the Web documentation.

The Field::field_struct points to a field_storage, which is defined by each backend. It implements storing and accessing the Field data, including buffers for storing haloes returned from MPI communication.

Template Parameters
TField content type
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h, and lattice.block.

Definition at line 62 of file field.h.

Constructor & Destructor Documentation

◆ Field()

template<typename T >
Field< T >::Field ( )
inline

Field constructor.

The following ways of constructing a Field object are:

Default constructor:

Assigns Field::fs to nullptr.

The field class implements the standard methods for accessing Fields. Hilapp replaces the parity acce...
Definition field.h:62

Copy constructor:

Copy from already existing field of similar type if conversion exists. For example for float to double.

Assignment from constant constructor:

Assign constant to field if conversion exists.

.
. Let a be a constant of type MyType1
.
Field<MyType2> f = a;
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 278 of file field.h.

Member Function Documentation

◆ check_alloc()

template<typename T >
void Field< T >::check_alloc ( )
inline

Allocate Field if it is not already allocated.

Check that Field is allocated, and if not do it (if not const). Must be called before the var is actually used. "hilapp" will generate these calls as needed!

If Field is const assert that the Field is allocated.

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 458 of file field.h.

◆ clear()

template<typename T >
void Field< T >::clear ( )
inline

Destroys field data.

don't call destructors when exiting - either MPI or cuda can already be off.

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h, and lattice.block.

Definition at line 406 of file field.h.

◆ conj()

template<typename T >
Field< T > Field< T >::conj ( ) const
inline

Returns field with all elements conjugated depending how conjugate is defined for type.

Returns
Field<T>
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1175 of file field.h.

◆ copy_boundary_condition()

template<typename T >
template<typename A >
void Field< T >::copy_boundary_condition ( const Field< A > &  rhs)
inline

Copy the boundary condition from another field.

Template Parameters
Athe type of the field which we are copying from
Parameters
rhsthe ohter Field
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 651 of file field.h.

◆ dagger()

template<typename T >
template<typename R = T, typename A = decltype(::dagger(std::declval<R>()))>
Field< A > Field< T >::dagger ( ) const
inline

Returns dagger or Hermitian conjugate of Field depending on how it is defined for Field type T.

Returns
Field<T>
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1188 of file field.h.

◆ FFT()

template<typename T >
Field< T > Field< T >::FFT ( const CoordinateVector dirs,
fft_direction  fftdir = fft_direction::forward 
) const

Field method for performing FFT.

By default calling without arguments will execute FFT in all directions.

.
. // Field f is defined
.
auto res = f.FFT() //Forward transform
auto res_2 = res.FFT(fft_direction::back) // res_2 is same as f

One can also specify the direction of the FFT with a coordinate vector:

.
. // Field f is defined
.
auto res = f.FFT(e_x) //Forward transform in x-direction
auto res_2 = res.FFT(e_X,fft_direction::back) // res_2 is same as f

With this in mind f.FFT(e_x+e_y+e_z) = f.FFT()

Template Parameters
T
Parameters
dirsDirection to perform FFT in, default is all directions
fftdirfft_direction::forward (default) or fft_direction::back
Returns
Field<T> Transformed field

Definition at line 449 of file fft.h.

◆ FFT_real_to_complex()

template<typename T >
Field< Complex< hila::arithmetic_type< T > > > Field< T >::FFT_real_to_complex ( fft_direction  fftdir = fft_direction::forward) const

FFT_real_to_complex: Field must be a real-valued field, result is a complex-valued field of the same type Implemented just by doing a FFT with a complex field with im=0; fft_direction::back gives a complex conjugate of the forward transform Result is f(-x) = f(L - x) = f(x)^*

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 474 of file fft.h.

◆ get_boundary_condition()

template<typename T >
hila::bc Field< T >::get_boundary_condition ( Direction  dir) const
inline

Get the boundary condition of the Field.

Parameters
dirBoundary condition in certain direction
Returns
hila::bc The boundary condition of the Field
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 627 of file field.h.

◆ get_element()

template<typename T >
const T Field< T >::get_element ( const CoordinateVector coord) const
inline

Get singular element which will be broadcast to all nodes.

Parameters
coordcoordinate of fetched element
Returns
const T
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1285 of file field.h.

◆ get_value_at()

template<typename T >
auto Field< T >::get_value_at ( const unsigned  i) const
inline

Get an individual element outside a loop. This is also used as a getter in the vanilla code.

Parameters
i
Returns
auto
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 685 of file field.h.

◆ imag()

template<typename T >
template<typename R = T, typename A = decltype(::imag(std::declval<R>()))>
Field< A > Field< T >::imag ( ) const
inline

Returns imaginary part of Field.

Template Parameters
RField current type
AField imaginary part type
Returns
Field
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1216 of file field.h.

◆ is_allocated()

template<typename T >
bool Field< T >::is_allocated ( ) const
inline

Returns true if the Field data has been allocated.

Returns
bool
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 421 of file field.h.

◆ is_initialized()

template<typename T >
bool Field< T >::is_initialized ( Parity  p) const
inline

Returns true if the Field has been assigned to.

Parameters
pField parity
Returns
bool
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 430 of file field.h.

◆ mylattice()

template<typename T >
Lattice Field< T >::mylattice ( ) const
inline

Return the lattice to which this field instance belongs to.

Useful for switching lattices,for example lattice.switch_to(a.mylattice()); changes active lattice to the one to which Field variable a belongs to.

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1599 of file field.h.

◆ norm()

template<typename T >
double Field< T >::norm ( )
inline

Norm.

Norm of Field depending on how it is defined for Field type T

\(|f| = \sqrt{\sum_{\forall x \in f} |x|^2}\)

If norm is defined for type T then the definition can be seen on the types documentation page.

Returns
double
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1166 of file field.h.

◆ operator*()

template<typename A , typename B >
auto operator* ( const Field< A > &  lhs,
const Field< B > &  rhs 
) -> Field<hila::type_mul<A, B>>

Multiplication operator.

Multiplication operator can be called in the following ways as long as types are convertible.

Multiplication with field:

. . .
f = f * g;

Multiplication with scalar:

. . .
f = f * a;
f = a * f;
Complex definition.
Definition cmplx.h:58
Template Parameters
AType of lhs Field
BType of rhs field
Parameters
lhsFirst (left) Field
rhsSecond (right) Field
Returns
Field<hila::type_plus<A, B>> Multiplied Field

Definition at line 1705 of file field.h.

◆ operator*=()

template<typename T >
template<typename A , std::enable_if_t< std::is_convertible< hila::type_mul< T, A >, T >::value, int > = 0>
Field< T > & Field< T >::operator*= ( const Field< A > &  rhs)
inline

Product assignment operator.

Product assignment operator can be called in the following ways as long as types are convertible.

Product assignment with field:

. . .
f *= g;

Product assignment with scalar:

. . .
f *= a;
Template Parameters
AType of r.h.s element
Parameters
rhsField to compute product with
Returns
Field<T>&
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 928 of file field.h.

◆ operator+()

template<typename T >
Field< T > operator+ ( ) const
inline

Summation operator.

Summation operator can be called in the following ways as long as types are convertible.

unary operator:

Acts as identity of field

. . .
f = +f;

Summation with field:

. . .
f = f + g;

Summation with scalar:

. . .
f = f + a;
Template Parameters
AType of lhs Field
BType of rhs field
Parameters
lhsFirst (left) Field
rhsSecond (right) Field
Returns
Field<hila::type_plus<A, B>> Summed Field
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1067 of file field.h.

◆ operator+=()

template<typename T >
template<typename A , std::enable_if_t< std::is_convertible< hila::type_plus< T, A >, T >::value, int > = 0>
Field< T > & Field< T >::operator+= ( const Field< A > &  rhs)
inline

Addition assignment operator.

Addition assignmen operator can be called in the following ways as long as types are convertible.

Addition assignment with field:

. . .
f += g;

Addition assignment with scalar:

. . .
f += a;
Template Parameters
AType of r.h.s element
Parameters
rhsField to sum
Returns
Field<T>&
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 862 of file field.h.

◆ operator-()

template<typename T >
Field< T > operator- ( ) const
inline

Subtraction operator.

Subtraction operator can be called in the following ways as long as types are convertible.

unary operator:

Acts as negation of field

Subtraction with field:

. . .
f = f - g;

Subtraction with scalar:

. . .
f = f - a;
f = a - f;
Template Parameters
AType of lhs Field
BType of rhs field
Parameters
lhsFirst (left) Field
rhsSecond (right) Field
Returns
Field<hila::type_plus<A, B>> Subtracted Field
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1111 of file field.h.

◆ operator-=()

template<typename T >
template<typename A , std::enable_if_t< std::is_convertible< hila::type_minus< T, A >, T >::value, int > = 0>
Field< T > & Field< T >::operator-= ( const Field< A > &  rhs)
inline

Subtraction assignment operator.

Subtraction assignment operator can be called in the following ways as long as types are convertible.

Subtraction assignment with field:

. . .
f -= g;

Subtraction assignment with scalar:

. . .
f -= a;
Template Parameters
AType of r.h.s element
Parameters
rhsField to subtract
Returns
Field<T>&
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 895 of file field.h.

◆ operator/()

template<typename A , typename B >
auto operator/ ( const Field< A > &  l,
const Field< B > &  r 
) -> Field<hila::type_div<A, B>>

Division operator.

Division operator can be called in the following ways as long as types are convertible.

Division with field:

. . .
f = f / g;

Division with scalar:

. . .
f = f / a;
f = a / f;
Template Parameters
AType of lhs Field
BType of rhs field
Parameters
lhsFirst (left) Field
rhsSecond (right) Field
Returns
Field<hila::type_plus<A, B>> Subtracted Field

Definition at line 1765 of file field.h.

◆ operator/=()

template<typename T >
template<typename A , std::enable_if_t< std::is_convertible< hila::type_div< T, A >, T >::value, int > = 0>
Field< T > & Field< T >::operator/= ( const Field< A > &  rhs)
inline

Division assignment operator.

Division assignment operator can be called in the following ways as long as types are convertible.

Division assignment with field:

. . .
f /= g;

Division assignment with scalar:

. . .
f /= a;
Template Parameters
AType of r.h.s element
Parameters
rhsField to divide with
Returns
Field<T>&
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 961 of file field.h.

◆ operator=() [1/2]

template<typename T >
Field< T > & Field< T >::operator= ( const Field< T > &  rhs)
inline

Assignment operator.

Assignment can be performed in the following ways:

Assignment from field:

Assignment from field is possible if types are the same or convertible

f = g;

Assignment from scalar:

Assignment from scalar is possible if scalar is convertible to Field type

f = a; // a is a scalar
Parameters
rhs
Returns
Field<T>& Assigned field
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 750 of file field.h.

◆ operator=() [2/2]

template<typename T >
Field< T > & Field< T >::operator= ( Field< T > &&  rhs)
inline

Move Assignment.

Parameters
rhs
Returns
Field<T>&

Definition at line 823 of file field.h.

◆ operator==()

template<typename T >
template<typename S >
bool Field< T >::operator== ( const Field< S > &  rhs) const
inline

Field comparison operator.

Fields are equal if their content is element-by-element equal

Parameters
rhsField to compare Field with
Returns
true
false
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1126 of file field.h.

◆ product()

template<typename T >
T Field< T >::product ( Parity  par = Parity::all,
bool  allreduce = true 
) const

Product reduction of Field.

The product in the reduction is defined by the Field type T.

Parameters
parParity
allreduceSwitch to turn on or off MPI allreduce
Returns
T Field element type reduced to one element
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 296 of file reduction.h.

◆ real()

template<typename T >
template<typename R = T, typename A = decltype(::real(std::declval<R>()))>
Field< A > Field< T >::real ( ) const
inline

Returns real part of Field.

Template Parameters
RField current type
AField real part type
Returns
Field
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1202 of file field.h.

◆ reflect()

template<typename T >
Field< T > Field< T >::reflect ( const CoordinateVector dirs) const

Reflect the Field around the desired axis.

Field<T>::reflect() reflects the field around the desired axis This is here because it uses similar communications as fft TODO: refactorise so that there is separate "make columns" class!

Can be called in the following ways:

Reflect on all axes:

.
.
.
f.reflect()
Template Parameters
T
Parameters
dirs
Returns
Field<T>

Definition at line 594 of file fft.h.

◆ set_boundary_condition()

template<typename T >
void Field< T >::set_boundary_condition ( Direction  dir,
hila::bc  bc 
)
inline

Set the boundary condition in a given Direction (periodic or antiperiodic)

Parameters
dirDirection of boundary condition
bcField boundary condition
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 579 of file field.h.

◆ set_element()

template<typename T >
template<typename A , std::enable_if_t< std::is_assignable< T &, A >::value, int > = 0>
const T Field< T >::set_element ( const CoordinateVector coord,
const A &  value 
)
inline

Get singular element which will be broadcast to all nodes.

Set a single element. Assuming that each node calls this with the same value, it is sufficient to set the element locally

Works as long as value type A and field type T match

Template Parameters
Atype
Parameters
coordCoordinate of element to be set
valueValue to be set
Returns
const T
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1268 of file field.h.

◆ set_value_at()

template<typename T >
template<typename A >
void Field< T >::set_value_at ( const A &  value,
unsigned  i 
)
inline

Set an individual element outside a loop. This is also used as a getter in the vanilla code.

Parameters
i
Returns
auto
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 710 of file field.h.

◆ shift()

template<typename T >
Field< T > & Field< T >::shift ( const CoordinateVector v,
Field< T > &  r,
Parity  par 
) const

Create a periodically shifted copy of the field.

this is currently OK only for short moves, very inefficient for longer moves

.
. //Field<MyType> f is defined
.
CoordinateVector v = {1,1,0}
f.shift(v) // Shift all elements of field once in the x direction and once in the y direction
f.shift(v,EVEN) // Shift all EVEN elements of field once in the x direction and once in the y
direction
constexpr Parity EVEN
bit pattern: 001
Parameters
v
par
r
Returns
Field<T>& returns a reference to res
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 162 of file field_comm.h.

◆ squarenorm()

template<typename T >
double Field< T >::squarenorm ( ) const
inline

Squarenorm.

Squarenorm of Field \(f\) is dependent on how it is defined for Field type T.

\(|f|^2 = \sum_{\forall x \in f} |x|^2\)

If squarenorm is defined for type T then the definition can be seen on the types documentation page.

Returns
double
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 1148 of file field.h.

◆ start_gather()

template<typename T >
dir_mask_t Field< T >::start_gather ( Direction  d,
Parity  p = ALL 
) const

start_gather(): Communicate the field at Parity par from Direction d. Uses accessors to prevent dependency on the layout. return the Direction mask bits where something is happening

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 262 of file field_comm.h.

◆ sum()

template<typename T >
T Field< T >::sum ( Parity  par = Parity::all,
bool  allreduce = true 
) const

Sum reduction of Field.

The sum in the reduction is defined by the Field type T

Parameters
parParity
allreduceSwitch to turn on or off MPI allreduce
Returns
T Field element type reduced to one element
Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 287 of file reduction.h.

◆ unblock_to()

template<typename T >
void Field< T >::unblock_to ( Field< T > &  target) const

Copy content to the argument Field on blocked (sparse) sites. a.unblock_to(b) is the inverse of a.block_from(b)

Leaves other sites of the argument Field unmodified.

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h, and lattice.block.

Definition at line 971 of file field_comm.h.

◆ write_subvolume()

template<typename T >
void Field< T >::write_subvolume ( std::ofstream &  outputfile,
const CoordinateVector cmin,
const CoordinateVector cmax,
int  precision = 6 
) const

Write a "subspace" of the original lattice Each element is written on a single line TODO: more formatting?

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 280 of file field_io.h.

Friends And Related Symbol Documentation

◆ pow()

template<typename T , typename P , typename R = decltype(pow(std::declval<T>()), std::declval<P>())>
Field< R > pow ( const Field< T > &  arg,
const P  p 
)
related

Power.

Parameters
pexponent to which Field element is raised to.

Definition at line 1931 of file field.h.

◆ squarenorm_relative()

template<typename A , typename B , typename R = decltype(std::declval<A>() - std::declval<B>())>
double squarenorm_relative ( const Field< A > &  a,
const Field< B > &  b 
)
related

Squarenorm relative \(|a-b|^2\).

Template Parameters
AElement type of Field a
BElement type of Field b
Parameters
aField a
bField b
Returns
double

Definition at line 1987 of file field.h.

Member Data Documentation

◆ fs

template<typename T >
field_struct* __restrict__ Field< T >::fs

Field::fs holds all of the field content in Field::field_struct.

The field_struct class is nonetheless marked as internal so the user need not worry about it.

Examples
/home/runner/work/HILA/HILA/libraries/plumbing/field.h.

Definition at line 244 of file field.h.


The documentation for this class was generated from the following files: