Autodocs
-
template<class T = double>
class Vector The Vector class creates a simple vector that is used for shared and distributed memory applications.
If MPI is available the vector is distributed across the different processes. If not, is is simply stored in aligned memory. The syntax is equivalent in both cases.
Template parameters
class T: Primitive datatype.
Default: double
Things to know
The real datatype of T can be accessed via:
Vector<T>::T_real
If T is
std::complex<double>
T_real refers todouble
. If T isdouble
T_real refers todouble
.The Vector is usually created by the ShellMatrix class.
The debugging level can be adjusted from 0 (no debugging) to 10 (maximal debugging):
#define dbug_level 10; // 10 maximal debugging, 0 is no debugging
Public Types
-
using T_real = decltype(std::real(std::declval<T>()))
Real part of T.
The real T can be accessed via:
Vector<T>::T_real
If T is
std::complex<double>
T_real refers todouble
. If T isdouble
T_real refers todouble
.
Public Functions
-
Vector(uint64_t dim_)
Creates vector of size dim_ with zeros.
- Parameters:
dim_ – Dimension
-
Vector(const Mat &mat)
Creates vector from Mat object of Petsc.
This constructor is only available if Petsc is included. The Vector uses the same memory layout as the input matrix.
- Parameters:
mat – Input Mat from Petsc
-
Vector(const Vec &vec)
Creates vector from Vec of Petsc.
This constructor is only available if Petsc is included. The Vector uses the same memory layout as the input vector.
- Parameters:
vec – Input Vec from Petsc
-
Vector(void)
Creates vector of size zero.
-
Vector(const std::vector<uint64_t> ownership_per_rank_)
Constructor with given memory layout.
Initialized a vector with zeros in the given memory layout. Only available if MPI is activated.
- Parameters:
ownership_per_rank_ – Sets ownership_per_rank
-
bool operator==(const Vector<T> &input) const
Overloaded == operator for comparison.
- Parameters:
input – Input vector
- Returns:
True
if input equals this,False
otherwise
-
Vector<T> &operator=(const Vector<T> &input)
Overloaded = operator.
- Parameters:
input – Input vector
- Returns:
New instance of Input vector
-
Vector<T> &operator=(const T &input)
Overloaded = operator.
- Parameters:
input – Scalar
- Returns:
New instance where all entries equal input
-
std::vector<std::vector<T>> get_density_matrix(const int32_t root = 0) const
Returns the density matrix (if possible).
We represent the density matrix \(\rho\) of dimension \(d\) as a large vector \(v\) with dimension \(d^2\) when we work with the Lindbladian OperatorType. This function returns the density matrix using a two-dimensional vector if the total dimension dim is the quadrat of an integer.
\[\rho[i,j] = v[i\cdot d+j]\]The density matrix will be only returned to the root rank if MPI is enabled.
- Parameters:
root – Designated MPI rank (Default: 0)
- Returns:
density matrix
-
EigenMatrixType get_density_matrix_in_eigen(const int32_t root = 0) const
Returns the density matrix (if possible).
This is a equivalent to get_density_matrix(root). It returns the density matrix using Eigen which has to be available.
- Parameters:
root – Designated MPI rank (Default: 0)
- Returns:
density matrix
-
int32_t print(void) const
Prints the vector.
- Returns:
error_code
-
T *get_data_ptr(const uint64_t offset = 0UL)
Pointer to the local data with an offset.
Points to the local data with an offset.
Note
If MPI is enabled, each rank points to its part with the same offset.
- Returns:
Pointer to data
-
const T *get_const_data_ptr(const uint64_t offset = 0UL) const
const Pointer to the local data with an offset.
Points to the local data with an offset. This allows only read access.
Note
If MPI is enabled, each rank points to its part with the same offset.
- Returns:
const Pointer to data
-
int32_t info(void) const
Prints information about the vector.
- Returns:
error_code
-
const T operator[](const uint64_t index) const
Overloaded [] operator for readout.
Returns value at index that can range from 0 to dim. The function is compatible with MPI and each processes obtains the same value stored at index.
Note
It is not possible to set an entry with this operator! Use set_index(index,value) to do this.
Warning
The performance is bad if multiple values are requested, use get_const_data_ptr(…) instead.
- Parameters:
index – Index
- Returns:
Vales
-
EigenVectorType create_EigenVector(const int32_t root = 0) const
Returns the Vector as an instance in Eigen.
This is only available if Eigen is included.
using MyVector=danceq::internal::Vector<ScalarType>; // Do somthing to V MyVector V(100); V.make_random(); // Retrieve the vector as EigenVectorType auto EigenV = create_EigenVector(void);
The vector matrix will be only returned to the root rank if MPI is enabled.
- Parameters:
root – Designated MPI rank (Default: 0)
- Returns:
Copy in form of EigenVectorType
-
Vec create_PetscVec(void) const
Returns Vec from Petsc.
This is only available if Petsc is included. It uses the same memory layout.
using MyVector=danceq::internal::Vector<ScalarType>; // Create a Mat object from Petsc (H is a Hamiltonian) Mat H_mat = H.create_PetscSparseMatrix(); // Create a Vec object from Petsc using H_mat MyVector V(H_mat); // Do somthing to V V.make_random(); // Retrieve the Petsc Vec Vec Psi = V.create_PetscVec();
- Returns:
Copy in form of Vec
-
int32_t set_index(uint64_t index, T value)
Sets a single number.
Warning
The performance is bad if multiple values should be modified, use get_data_ptr(…) instead.
- Parameters:
index – Index to be set
value – Value to be set
- Returns:
error_code
-
int32_t set_all(T value)
Sets all values to value.
- Returns:
error_code
-
int32_t make_random(uint32_t seed = -1)
Generates a random and normalized vector.
A random seed is generated if it is
-1
. Values are Gaussian random numbers with mean 0 and standard deviation 1. The vector is normalized.- Parameters:
seed – Seed for random generator
- Returns:
error_code
-
int32_t conj(void)
Conjugates the vector.
- Returns:
error_code
-
int32_t scale(const T scalar)
Scales the vector by scalar.
- Parameters:
scalar – Scalar
- Returns:
error_code
-
int32_t add(const Vector<T> &input)
Adds another vector input.
- Parameters:
input – Input vector
- Returns:
error_code
-
int32_t add_and_scale(const Vector<T> &input, const T scalar)
Adds another vector input that is scaled by scalar.
- Parameters:
input – Input vector
scalar – Scalar to scale the input Input vector
- Returns:
error_code
-
T_real get_norm(void) const
Returns the norm.
The norm is returned by the real part of the template parameter T. If T is real T_real and T are identical.
- Returns:
norm
-
Vector<T> operator+(const Vector<T> &input) const
Overloaded + operator.
Creates a New vector when calling \(C = A + B\).
- Parameters:
input – Input vector
- Returns:
New vector
-
Vector<T> operator-(const Vector<T> &input) const
Overloaded + operator.
Creates a New vector when calling \(C = A - B\).
- Parameters:
input – Input vector
- Returns:
New vector
-
Vector<T> operator*(const T scalar) const
Overloaded * operator to scale the vector.
Creates a New vector when calling \(B = s*A\). Here, \(s\) is a scalar.
- Parameters:
scalar – Scalar
- Returns:
New vector
-
T operator*(const Vector<T> &input) const
Overloaded * operator for the scalar product.
Returns
\[\sum_{i=0}^{\textbf{dim}-1} \texttt{std::conj}(a[i])*b[i]\]where \(A\) and \(B\) are vectors with entries \(a\) and \(b\), respectively.
- Parameters:
input – Input vector
- Returns:
inner product
-
std::vector<T> get_data(void) const
Returns data.
This returns only the local data.
- Returns:
data
-
std::vector<T> gather_data(const int32_t root = 0) const
Gathers the data into one vector.
This is equivalent get_data() without MPI. The input root is meaningless in this case.
- Parameters:
root – Designated MPI rank (Default: 0)
- Returns:
data
-
uint64_t get_dim(void) const
Returns dim.
- Returns:
dim
-
std::vector<uint64_t> get_ownership_per_rank(void) const
Returns ownership_per_rank.
- Returns:
ownership_per_rank
-
int32_t get_myrank(void) const
Returns ownership_per_rank.
- Returns:
ownership_per_rank
-
uint64_t get_start(void) const
Returns ownership_per_rank.
- Returns:
ownership_per_rank
-
uint64_t get_end(void) const
Returns start.
- Returns:
start
-
int32_t get_world_size(void) const
Returns end.
- Returns:
end
-
uint64_t get_mydim(void) const
Returns ownership_per_rank.
- Returns:
ownership_per_rank
-
MPI_Datatype get_MPI_SCALAR(void) const
Returns ownership_per_rank.
- Returns:
ownership_per_rank
-
MPI_Datatype get_MPI_SCALAR_REAL(void) const
Returns ownership_per_rank.
- Returns:
ownership_per_rank
Private Members
-
std::vector<T> data
Local data.
If MPI is used this refers only to the local data that is distributed across processes.
-
const uint64_t dim
Total dimension.
-
const MPI_Datatype MPI_SCALAR
Datatype used by MPI.
-
const MPI_Datatype MPI_SCALAR_REAL
Real datatype used by MPI.
-
int32_t world_size
Number of MPI processes.
-
int32_t myrank
MPI rank of this process.
-
uint64_t mydim
Local dimension of rank.
-
uint64_t start
Start of local rows enumerated by the Basis in ptr.
-
uint64_t end
End of local rows enumerated by the Basis in ptr.
-
const std::vector<uint64_t> ownership_per_rank
Rows that mark the memory distribution per rank.
Rank i is in charge of rows ownership_per_rank[i-1] to (not including the last) ownership_per_rank[i]. Note that ownership_per_rank[-1] is zero.