SPI
SLEPc PETSc Interface is used for rapid development and intuitive matrix operations similar to MatLab or NumPy style of operations
Data Structures | Typedefs | Enumerations | Functions
SPI Namespace Reference

Data Structures

struct  SPIbaseflow
 
struct  SPIgrid1D
 Class to contain various grid parameters. More...
 
struct  SPIgrid2D
 
struct  SPIMat
 
struct  SPIparams
 
struct  SPIVec
 

Typedefs

template<class T >
using Block2D = std::vector< std::vector< T > >
 

Enumerations

enum  gridtype {
  FD, FDperiodic, Fourier, Chebyshev,
  UltraS
}
 enumeration of grid types More...
 

Functions

SPIbaseflow blasius (SPIparams &params, SPIgrid1D &grid)
 
int _bblf (const PetscScalar input[3], PetscScalar output[3])
 
SPIbaseflow channel (SPIparams &params, SPIgrid1D &grid)
 
PetscInt factorial (PetscInt n)
 compute the factorial of n. This is needed for get_D_Coeffs function More...
 
SPIVec get_D_Coeffs (SPIVec &s, PetscInt d)
 get the coefficients of the given stencil. More...
 
SPIMat map_D (SPIMat D, SPIVec y, PetscInt d, PetscInt order)
 map the derivative operator to the proper y grid More...
 
SPIMat set_D (SPIVec &y, PetscInt d, PetscInt order, PetscBool uniform)
 set the derivative operator for the proper y grid if uniform=false. Uses map_D function More...
 
SPIMat set_D_periodic (SPIVec &y, PetscInt d, PetscInt order)
 set the derivative operator for the proper periodic grid assuming uniform discretization. More...
 
SPIVec set_FD_stretched_y (PetscScalar y_max, PetscInt ny, PetscScalar delta)
 set stretched grid from [0,y_max] using tanh stretching for use with finite difference operators More...
 
SPIMat set_D_Chebyshev (SPIVec &x, PetscInt d, PetscBool need_map)
 set a Chebyshev collocated operator acting with respect to the collocated grid More...
 
SPIMat map_D_Chebyshev (SPIVec &x, PetscInt d)
 map a Chebyshev collocated operator acting with respect to the stretched collocated grid More...
 
SPIVec set_Cheby_stretched_y (PetscScalar y_max, PetscInt ny, PetscScalar yi)
 create a stretched Chebyshev grid from [0,y_max] with yi being the midpoint location for the number of points More...
 
SPIVec set_Cheby_mapped_y (PetscScalar a, PetscScalar b, PetscInt ny)
 create a stretched Chebyshev grid from [a,b] using default Chebfun like mapping More...
 
SPIVec set_Cheby_y (PetscInt ny)
 create a Chebyshev grid from [-1,1] More...
 
SPIVec set_Fourier_t (PetscScalar T, PetscInt nt)
 create a Fourier grid from [0,T] More...
 
SPIMat set_D_Fourier (SPIVec t, PetscInt d)
 create a Fourier derivative operator acting on grid t More...
 
std::tuple< SPIMat, SPIMatset_D_UltraS (SPIVec &x, PetscInt d)
 set a UltraSpherical operator acting with respect to the collocated grid (keeps everything in UltraSpherical space), take in Chebyshev coefficients and outputs coefficients in C(d) coefficient space More...
 
std::tuple< SPIMat, SPIMatset_T_That (PetscInt n)
 set a T and That operators acting with respect to the Chebyshev collocated grid, That: physical -> Chebyshev coefficient, T: Chebyshev coefficient -> physical More...
 
SPIVec SPIVec1Dto2D (SPIgrid2D &grid2D, SPIVec &u)
 expand a 1D vector to a 2D vector copying data along time dimension More...
 
PetscScalar integrate (const SPIVec &a, SPIgrid1D &grid)
 integrate a vector of chebyshev Coefficients on a physical grid More...
 
PetscScalar integrate (const SPIVec &a, SPIgrid2D &grid)
 integrate a vector of chebyshev Coefficients on a physical grid More...
 
SPIVec proj (SPIVec &u, SPIVec &v, SPIgrid1D &grid)
 
SPIVec proj (SPIVec &u, SPIVec &v, SPIgrid2D &grid)
 
std::vector< SPIVecorthogonalize (std::vector< SPIVec > &x, SPIgrid1D &grid)
 
std::vector< SPIVecorthogonalize (std::vector< SPIVec > &x, SPIgrid2D &grid)
 
SPIMat interp1D_Mat (SPIgrid1D &grid1, SPIgrid1D &grid2)
 
SPIMat interp1D_Mat (SPIVec &y1, SPIVec &y2)
 
SPIMat dft (PetscInt nt)
 
std::tuple< SPIMat, SPIMat, SPIMat, SPIMatdft_dftinv_Ihalf_Ihalfn (PetscInt nt)
 
std::tuple< PetscScalar, SPIVecLST_temporal (SPIparams &params, SPIgrid1D &grid, SPIbaseflow &baseflow, SPIVec q)
 solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with alpha being pure real, and omega the eigenvalue More...
 
std::tuple< PetscScalar, SPIVecLST_spatial (SPIparams &params, SPIgrid1D &grid, SPIbaseflow &baseflow, SPIVec q)
 solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue More...
 
std::tuple< PetscScalar, PetscScalar, SPIVec, SPIVecLST_spatial_cg (SPIparams &params, SPIgrid1D &grid, SPIbaseflow &baseflow)
 solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue More...
 
std::tuple< PetscScalar, PetscScalar, SPIVec, SPIVecLSTNP_spatial (SPIparams &params, SPIgrid1D &grid, SPIbaseflow &baseflow, SPIVec ql, SPIVec qr)
 solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue More...
 
std::tuple< PetscScalar, SPIVecLSTNP_spatial_right (SPIparams &params, SPIgrid1D &grid, SPIbaseflow &baseflow, SPIVec qr)
 solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue More...
 
std::tuple< PetscScalar, SPIVecLSTNP_spatial_right2 (SPIparams &params, SPIgrid1D &grid, SPIbaseflow &baseflow, SPIVec qr)
 solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue More...
 
std::tuple< std::vector< PetscScalar >, std::vector< SPIVec > > LSTNP_spatials_right (SPIparams &params, SPIgrid1D &grid, SPIbaseflow &baseflow, std::vector< PetscScalar > &alphas, std::vector< SPIVec > &qrs)
 solve the local stability theory problem for the linearized Navier-Stokes equations using non-parallel baseflow with omega being pure real, and alpha the eigenvalue More...
 
SPIMat operator* (const PetscScalar a, const SPIMat A)
 a*A operation to be equivalent to A*a More...
 
SPIMat operator* (const SPIMat A, const PetscScalar a)
 A*a operation to be equivalent to A*a. More...
 
SPIVec operator/ (const SPIVec &b, const SPIMat &A)
 Solve linear system, Ax=b using b/A notation. More...
 
SPIMat operator^ (const PetscScalar a, const SPIMat &A)
 Y=a^A operation. More...
 
SPIVec solve (const SPIMat &A, const SPIVec &b)
 Solve linear system, Ax=b using solve(A,b) notation. More...
 
SPIMat eye (const PetscInt n)
 create, form, and return identity matrix of size n More...
 
SPIMat inv (const SPIMat &A)
 Solve linear system A*Ainv=B using MatMatSolve. More...
 
SPIMat zeros (const PetscInt m, const PetscInt n)
 create, form, and return zeros matrix of size mxn More...
 
SPIMat diag (const SPIVec &d, const PetscInt k)
 set diagonal of matrix More...
 
SPIMat kron (const SPIMat &A, const SPIMat &B)
 set kronecker inner product of two matrices More...
 
std::tuple< std::vector< PetscReal >, std::vector< SPIVec >, std::vector< SPIVec > > svd (const SPIMat &A)
 solve SVD problem of A=U*E*V^H for skinny A matrix More...
 
SPIVec lstsq (const SPIMat &A, SPIVec &y)
 solve least squares problem of A*x=y for skinny A matrix using SVD More...
 
SPIVec lstsq (const std::vector< SPIVec > &A, SPIVec &y)
 solve least squares problem of A*x=y for skinny A matrix (or vectors of columns vectors) using SVD More...
 
std::tuple< PetscScalar, SPIVec, SPIVeceig (const SPIMat &A, const SPIMat &B, const PetscScalar target, const PetscReal tol, const PetscInt max_iter)
 solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) More...
 
std::tuple< PetscScalar, SPIVeceig_right (const SPIMat &A, const SPIMat &B, const PetscScalar target, const PetscReal tol, const PetscInt max_iter)
 solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) More...
 
std::tuple< PetscScalar, SPIVec, SPIVeceig_init (const SPIMat &A, const SPIMat &B, const PetscScalar target, const SPIVec &ql, const SPIVec &qr, PetscReal tol, const PetscInt max_iter)
 solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) More...
 
std::tuple< PetscScalar, SPIVeceig_init_right (const SPIMat &A, const SPIMat &B, const PetscScalar target, const SPIVec &qr, PetscReal tol, const PetscInt max_iter)
 solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) More...
 
std::tuple< std::vector< PetscScalar >, std::vector< SPIVec > > eig_init_rights (const SPIMat &A, const SPIMat &B, const std::vector< PetscScalar > targets, const std::vector< SPIVec > &qrs, PetscReal tol, const PetscInt max_iter)
 solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) More...
 
std::tuple< PetscScalar, SPIVecpolyeig (const std::vector< SPIMat > &As, const PetscScalar target, const PetscReal tol, const PetscInt max_iter)
 solve general polynomial eigenvalue problem of (A0 + A1*alpha + A2*alpha^2 + ...)*x = 0 and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) More...
 
std::tuple< PetscScalar, SPIVecpolyeig_init (const std::vector< SPIMat > &As, const PetscScalar target, const SPIVec &qr, const PetscReal tol, const PetscInt max_iter)
 solve general polynomial eigenvalue problem of (A0 + A1*alpha + A2*alpha^2 + ...)*x = 0 and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) using initial subspace vector More...
 
SPIMat block ( const Block2D< SPIMat > Blocks)
 set block matrices using an input array of size rows*cols. Fills rows first More...
 
std::tuple< SPIMat, SPIMatmeshgrid (SPIVec &x, SPIVec &y)
 
PetscInt save (const SPIMat &A, const std::string filename)
 save matrix to filename in binary format (see Petsc documentation for format ) Format is (from Petsc documentation): int MAT_FILE_CLASSID int number of rows int number of columns int total number of nonzeros int *number nonzeros in each row int *column indices of all nonzeros (starting index is zero) PetscScalar *values of all nonzeros More...
 
PetscInt save (const std::vector< SPIMat > &As, const std::string filename)
 save matrices to filename in binary format (see Petsc documentation for format More...
 
PetscInt load (SPIMat &A, const std::string filename)
 load matrix from filename from binary format (works with save(SPIMat,std::string) function More...
 
PetscInt load (std::vector< SPIMat > &As, const std::string filename)
 load matrix from filename from binary format (works with save(SPIMat,std::string) function More...
 
PetscInt draw (const SPIMat &A)
 draw nonzero structure of matrix More...
 
template<class T >
SPIMat _Function_on_each_element (T(*f)(T const &), const SPIMat &A)
 take the function of each element in a matrix, e.g. (*f)(A(i,j)) for each i,j More...
 
SPIMat sin (const SPIMat &A)
 take the sin of each element in a matrix More...
 
SPIMat cos (const SPIMat &A)
 take the cos of each element in a matrix More...
 
SPIMat acos (const SPIMat &A)
 take the arccos of each element in a matrix More...
 
SPIMat tan (const SPIMat &A)
 take the tan of each element in a matrix More...
 
SPIMat abs (const SPIMat &A)
 take the abs of each element in a matrix More...
 
SPIMat operator% (const SPIMat &A, const SPIMat &B)
 take the elementwise multiplication of each element in a matrix More...
 
SPIMat orthogonalize (const std::vector< SPIVec > &x)
 
PetscInt printf (std::string msg,...)
 
PetscInt printfc (std::string msg, const PetscScalar val)
 
SPIVec operator/ (const PetscScalar a, const SPIVec &A)
 
SPIVec operator* (const PetscScalar a, const SPIVec &A)
 
SPIVec operator+ (const PetscScalar a, const SPIVec &A)
 
SPIVec operator- (const PetscScalar a, const SPIVec &A)
 
PetscInt save (const SPIVec &A, std::string filename)
 
PetscInt save (std::vector< PetscScalar > A, std::string variablename, std::string filename)
 
PetscInt save (std::vector< SPIVec > A, std::string variablename, std::string filename)
 
PetscInt load (SPIVec &A, const std::string filename)
 
SPIVec ones (const PetscInt rows)
 
SPIVec zeros (const PetscInt rows)
 
SPIVec conj (const SPIVec &A)
 
SPIVec real (const SPIVec &A)
 return the real part of the vector More...
 
SPIVec imag (const SPIVec &A)
 return the imaginary part of the vector More...
 
SPIVec linspace (const PetscScalar begin, const PetscScalar end, const PetscInt _rows)
 return linspace of number of rows equally spaced points between begin and end More...
 
SPIVec arange (const PetscScalar begin, const PetscScalar end, const PetscScalar stepsize)
 return a range of number of equally spaced points between begin and end by step size step More...
 
SPIVec arange (const PetscScalar end)
 
template<class T >
SPIVec _Function_on_each_element (T(*f)(T const &), const SPIVec &A)
 take the function of each element in a vector, e.g. (*f)(A(i)) for each i More...
 
SPIVec sin (const SPIVec &A)
 take the sin of each element in a vector More...
 
SPIVec cos (const SPIVec &A)
 take the cos of each element in a vector More...
 
SPIVec tan (const SPIVec &A)
 take the tan of each element in a vector More...
 
SPIVec exp (const SPIVec &A)
 take the exp of each element in a vector More...
 
SPIVec log (const SPIVec &A)
 take the log (natural log) of each element in a vector More...
 
SPIVec log10 (const SPIVec &A)
 take the log10 of each element in a vector More...
 
SPIVec sinh (const SPIVec &A)
 take the sinh of each element in a vector More...
 
SPIVec cosh (const SPIVec &A)
 take the cosh of each element in a vector More...
 
SPIVec tanh (const SPIVec &A)
 take the tanh of each element in a vector More...
 
SPIVec asin (const SPIVec &A)
 take the asin of each element in a vector More...
 
SPIVec acos (const SPIVec &A)
 take the acos of each element in a vector More...
 
SPIVec atan (const SPIVec &A)
 take the atan of each element in a vector More...
 
SPIVec asinh (const SPIVec &A)
 take the asinh of each element in a vector More...
 
SPIVec acosh (const SPIVec &A)
 take the acosh of each element in a vector More...
 
SPIVec atanh (const SPIVec &A)
 take the atanh of each element in a vector More...
 
SPIVec sqrt (const SPIVec &A)
 take the atanh of each element in a vector More...
 
SPIVec erf (const SPIVec &A)
 take the erf of each element in a vector More...
 
SPIVec erfc (const SPIVec &A)
 take the erfc(z) = 1-erf(z) of each element in a vector More...
 
template<class T >
SPIVec _Function_on_each_element (T(*f)(T const &, T const &), const SPIVec &A, SPIVec &B)
 function to take element by element of two vectors e.g. (*f)(A(i),B(i)) for all i More...
 
SPIVec pow (const SPIVec &A, SPIVec &B)
 take the pow of each element in the vectors More...
 
SPIVec pow (const SPIVec &A, PetscScalar b)
 take the pow of each element in the vector (A^b) More...
 
SPIVec abs (const SPIVec &A)
 take the absolute value of a vector More...
 
PetscScalar sum (SPIVec x1)
 take the sum of a vector More...
 
PetscScalar integrate_coeffs (const SPIVec &a)
 integrate a vector of chebyshev Coefficients on a grid More...
 
PetscScalar integrate_coeffs (const SPIVec &a, const SPIVec &y)
 integrate a vector of chebyshev Coefficients on a stretched grid More...
 
PetscScalar dot (SPIVec x, SPIVec y)
 take the inner product of the two vectors (i.e. y^H x) where ^H is the complex conjugate transpose More...
 
PetscReal L2 (SPIVec x1, const SPIVec x2, NormType type)
 calculate the \( L_2 \) norm of the difference between \(x_1\) and \(x_2\) vectors. More...
 
PetscReal L2 (const SPIVec x1, NormType type)
 calculate the \( L_2 \) norm of the vector More...
 
SPIVec diff (SPIVec x)
 diff of the vector (see numpy.diff) More...
 
PetscScalar trapz (SPIVec y)
 trapezoidal integration of y with unity coordinate spacing, \(\int y dx \) More...
 
PetscScalar trapz (SPIVec y, SPIVec x)
 trapezoidal integration of y with x coordinates, \(\int y dx \) More...
 
PetscInt draw (const SPIVec &x)
 draw nonzero structure and wait at command line input More...
 

Typedef Documentation

◆ Block2D

template<class T >
using SPI::Block2D = typedef std::vector<std::vector<T> >

Definition at line 15 of file SPIMat.hpp.

Enumeration Type Documentation

◆ gridtype

enumeration of grid types

Enumerator
FD 

finite difference grid

FDperiodic 

finite difference grid on periodic grid

Fourier 

Fourier transform collocated grid.

Chebyshev 

Chebyshev collocated grid.

UltraS 

UltraSpherical grid and derivatives.

Definition at line 25 of file SPIgrid.hpp.

Function Documentation

◆ _bblf()

int SPI::_bblf ( const PetscScalar  input[3],
PetscScalar  output[3] 
)

Definition at line 243 of file SPIbaseflow.cpp.

◆ _Function_on_each_element() [1/3]

template<class T >
SPIMat SPI::_Function_on_each_element ( T(*)(T const &)  f,
const SPIMat A 
)

take the function of each element in a matrix, e.g. (*f)(A(i,j)) for each i,j

Parameters
[in]f[in] function handle to pass in e.g. std::sin<PetscReal>
[in]A[in] matrix to perform function on each element

Definition at line 2011 of file SPIMat.cpp.

◆ _Function_on_each_element() [2/3]

template<class T >
SPIVec SPI::_Function_on_each_element ( T(*)(T const &)  f,
const SPIVec A 
)

take the function of each element in a vector, e.g. (*f)(A(i)) for each i

Parameters
[in]f[in] function handle to pass in e.g. std::sin<PetscReal>
[in]A[in] vector to perform function on each element

Definition at line 693 of file SPIVec.cpp.

◆ _Function_on_each_element() [3/3]

template<class T >
SPIVec SPI::_Function_on_each_element ( T(*)(T const &, T const &)  f,
const SPIVec A,
SPIVec B 
)

function to take element by element of two vectors e.g. (*f)(A(i),B(i)) for all i

Parameters
[in]f[in] function handle to pass in e.g. std::pow<PetscReal>
[in]A[in] first vector to perform function on each element
[in]B[in] second vector

Definition at line 762 of file SPIVec.cpp.

◆ abs() [1/2]

SPIMat SPI::abs ( const SPIMat A)

take the abs of each element in a matrix

Definition at line 2045 of file SPIMat.cpp.

◆ abs() [2/2]

SPIVec SPI::abs ( const SPIVec A)

take the absolute value of a vector

Definition at line 798 of file SPIVec.cpp.

◆ acos() [1/2]

SPIMat SPI::acos ( const SPIMat A)

take the arccos of each element in a matrix

Definition at line 2029 of file SPIMat.cpp.

◆ acos() [2/2]

SPIVec SPI::acos ( const SPIVec A)

take the acos of each element in a vector

Definition at line 731 of file SPIVec.cpp.

◆ acosh()

SPIVec SPI::acosh ( const SPIVec A)

take the acosh of each element in a vector

Definition at line 737 of file SPIVec.cpp.

◆ arange() [1/2]

SPIVec SPI::arange ( const PetscScalar  begin,
const PetscScalar  end,
const PetscScalar  stepsize = 1 
)

return a range of number of equally spaced points between begin and end by step size step

Parameters
[in]begin[in] beginning scalar of equally spaced points
[in]end[in] end scalar of equally spaced points
[in]stepsize[in] step size for equally spaced points

Definition at line 667 of file SPIVec.cpp.

◆ arange() [2/2]

SPIVec SPI::arange ( const PetscScalar  end)
Parameters
[in]end[in] end scalar of equally spaced points

Definition at line 684 of file SPIVec.cpp.

◆ asin()

SPIVec SPI::asin ( const SPIVec A)

take the asin of each element in a vector

Definition at line 729 of file SPIVec.cpp.

◆ asinh()

SPIVec SPI::asinh ( const SPIVec A)

take the asinh of each element in a vector

Definition at line 735 of file SPIVec.cpp.

◆ atan()

SPIVec SPI::atan ( const SPIVec A)

take the atan of each element in a vector

Definition at line 733 of file SPIVec.cpp.

◆ atanh()

SPIVec SPI::atanh ( const SPIVec A)

take the atanh of each element in a vector

Definition at line 739 of file SPIVec.cpp.

◆ blasius()

SPIbaseflow SPI::blasius ( SPIparams params,
SPIgrid1D grid 
)
Parameters
[in]params[in] parameters such as x and nu (freestream Uinf=1)
[in]grid[in] grid containing wall-normal points

Definition at line 90 of file SPIbaseflow.cpp.

◆ block()

SPIMat SPI::block ( const Block2D< SPIMat Blocks)

set block matrices using an input array of size rows*cols. Fills rows first

Returns
new matrix with inserted blocks
Parameters
[in]Blocks[in] array of matrices to set into larger matrix e.g. [ A, B, C, D ]

Definition at line 1857 of file SPIMat.cpp.

◆ channel()

SPIbaseflow SPI::channel ( SPIparams params,
SPIgrid1D grid 
)
Parameters
[in]params[in] parameters such as x and nu (freestream Uinf=1)
[in]grid[in] grid containing wall-normal points

Definition at line 253 of file SPIbaseflow.cpp.

◆ conj()

SPIVec SPI::conj ( const SPIVec A)

return the conjugate of the vector

Returns
conjugate of A
See also
SPIVec::conj()
Parameters
[in]A[in] vector to conjugate

Definition at line 622 of file SPIVec.cpp.

◆ cos() [1/2]

SPIMat SPI::cos ( const SPIMat A)

take the cos of each element in a matrix

Definition at line 2027 of file SPIMat.cpp.

◆ cos() [2/2]

SPIVec SPI::cos ( const SPIVec A)

take the cos of each element in a vector

Definition at line 708 of file SPIVec.cpp.

◆ cosh()

SPIVec SPI::cosh ( const SPIVec A)

take the cosh of each element in a vector

Definition at line 725 of file SPIVec.cpp.

◆ dft()

SPIMat SPI::dft ( PetscInt  nt)
Parameters
[in]nt[in] size of matrix to create

Definition at line 1573 of file SPIgrid.cpp.

◆ dft_dftinv_Ihalf_Ihalfn()

std::tuple< SPIMat, SPIMat, SPIMat, SPIMat > SPI::dft_dftinv_Ihalf_Ihalfn ( PetscInt  nt)
Parameters
[in]nt[in] size of matrix to create

Definition at line 1589 of file SPIgrid.cpp.

◆ diag()

SPIMat SPI::diag ( const SPIVec d,
const PetscInt  k 
)

set diagonal of matrix

Returns
new matrix with a diagonal vector set as the main diagonal
Parameters
[in]d[in] diagonal vector to set along main diagonal
[in]k[in] diagonal to set, k=0 is main diagonal, k=1 is offset one in positive direction

Definition at line 728 of file SPIMat.cpp.

◆ diff()

SPIVec SPI::diff ( SPIVec  x)

diff of the vector (see numpy.diff)

Returns
y[i] = x[i+1]-x[i] for i=0,1,...,x.rows-2
Parameters
[in]x[in] vector to diff (x[i+1]-x[i])

Definition at line 881 of file SPIVec.cpp.

◆ dot()

PetscScalar SPI::dot ( SPIVec  x,
SPIVec  y 
)

take the inner product of the two vectors (i.e. y^H x) where ^H is the complex conjugate transpose

Parameters
[in]x[in] first vector in inner product
[in]y[in] second vector in inner product (this one gets the complex conjugate transpose)

Definition at line 788 of file SPIVec.cpp.

◆ draw() [1/2]

PetscInt SPI::draw ( const SPIMat A)

draw nonzero structure of matrix

Returns
0 if successful
Parameters
[in]A[in] A to draw nonzero structure

Definition at line 1994 of file SPIMat.cpp.

◆ draw() [2/2]

PetscInt SPI::draw ( const SPIVec x)

draw nonzero structure and wait at command line input

Parameters
[in]x[in] vector to draw

Definition at line 932 of file SPIVec.cpp.

◆ eig()

std::tuple< PetscScalar, SPIVec, SPIVec > SPI::eig ( const SPIMat A,
const SPIMat B,
const PetscScalar  target,
const PetscReal  tol,
const PetscInt  max_iter 
)

solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector)

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha, eig_vector) = eig(A,B,0.1+0.4*PETSC_i)
Parameters
[in]A[in] A in Ax=kBx generalized eigenvalue problem
[in]B[in] B in Ax=kBx generalized eigenvalue problem
[in]target[in] target eigenvalue to solve for
[in]tol[in] tolerance of eigenvalue solver
[in]max_iter[in] maximum number of iterations

Definition at line 928 of file SPIMat.cpp.

◆ eig_init()

std::tuple< PetscScalar, SPIVec, SPIVec > SPI::eig_init ( const SPIMat A,
const SPIMat B,
const PetscScalar  target,
const SPIVec ql,
const SPIVec qr,
PetscReal  tol,
const PetscInt  max_iter 
)

solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector)

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha, eig_vector) = eig(A,B,0.1+0.4*PETSC_i) using initial subspace defined by a vector
Parameters
[in]A[in] A in Ax=kBx generalized eigenvalue problem
[in]B[in] B in Ax=kBx generalized eigenvalue problem
[in]target[in] target eigenvalue to solve for
[in]ql[in] initial subspace vector for EPS solver (left eigenvector)
[in]qr[in] initial subspace vector for EPS solver (right eigenvector)
[in]tol[in] tolerance of eigenvalue solver
[in]max_iter[in] maximum number of iterations

Definition at line 1217 of file SPIMat.cpp.

◆ eig_init_right()

std::tuple< PetscScalar, SPIVec > SPI::eig_init_right ( const SPIMat A,
const SPIMat B,
const PetscScalar  target,
const SPIVec qr,
PetscReal  tol,
const PetscInt  max_iter 
)

solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector)

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha, eig_vector) = eig(A,B,0.1+0.4*PETSC_i) using initial subspace defined by a vector
Parameters
[in]A[in] A in Ax=kBx generalized eigenvalue problem
[in]B[in] B in Ax=kBx generalized eigenvalue problem
[in]target[in] target eigenvalue to solve for
[in]qr[in] initial subspace vector for EPS solver (right eigenvector)
[in]tol[in] tolerance of eigenvalue solver
[in]max_iter[in] maximum number of iterations

Definition at line 1368 of file SPIMat.cpp.

◆ eig_init_rights()

std::tuple< std::vector< PetscScalar >, std::vector< SPIVec > > SPI::eig_init_rights ( const SPIMat A,
const SPIMat B,
const std::vector< PetscScalar >  targets,
const std::vector< SPIVec > &  qrs,
PetscReal  tol,
const PetscInt  max_iter 
)

solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector)

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha, eig_vector) = eig(A,B,0.1+0.4*PETSC_i) using initial subspace defined by a vector
Parameters
[in]A[in] A in Ax=kBx generalized eigenvalue problem
[in]B[in] B in Ax=kBx generalized eigenvalue problem
[in]targets[in] target eigenvalue to solve for
[in]qrs[in] initial subspace vector for EPS solver (right eigenvector)
[in]tol[in] tolerance of eigenvalue solver
[in]max_iter[in] maximum number of iterations

Definition at line 1521 of file SPIMat.cpp.

◆ eig_right()

std::tuple< PetscScalar, SPIVec > SPI::eig_right ( const SPIMat A,
const SPIMat B,
const PetscScalar  target,
const PetscReal  tol,
const PetscInt  max_iter 
)

solve general eigenvalue problem of Ax = kBx and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector)

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha, eig_vector) = eig(A,B,0.1+0.4*PETSC_i)
Parameters
[in]A[in] A in Ax=kBx generalized eigenvalue problem
[in]B[in] B in Ax=kBx generalized eigenvalue problem
[in]target[in] target eigenvalue to solve for
[in]tol[in] tolerance of eigenvalue solver
[in]max_iter[in] maximum number of iterations

Definition at line 1072 of file SPIMat.cpp.

◆ erf()

SPIVec SPI::erf ( const SPIVec A)

take the erf of each element in a vector

Definition at line 743 of file SPIVec.cpp.

◆ erfc()

SPIVec SPI::erfc ( const SPIVec A)

take the erfc(z) = 1-erf(z) of each element in a vector

Definition at line 752 of file SPIVec.cpp.

◆ exp()

SPIVec SPI::exp ( const SPIVec A)

take the exp of each element in a vector

Definition at line 712 of file SPIVec.cpp.

◆ eye()

SPIMat SPI::eye ( const PetscInt  n)

create, form, and return identity matrix of size n

Returns
identity matrix of size nxn
Parameters
[in]n[in] n size of square identity matrix

Definition at line 697 of file SPIMat.cpp.

◆ factorial()

PetscInt SPI::factorial ( PetscInt  n)

compute the factorial of n. This is needed for get_D_Coeffs function

Returns
factorial of n
Parameters
[in]n[in] integer to compute factorial of

Definition at line 6 of file SPIgrid.cpp.

◆ get_D_Coeffs()

SPIVec SPI::get_D_Coeffs ( SPIVec s,
PetscInt  d 
)

get the coefficients of the given stencil.

Returns
coefficients of stencil for derivative at 0 in s n
Parameters
[in]s[in] vector of stencil points e.g. [-3,-2,-1,0,1]
[in]d[in] order of desired derivative e.g. 1 or 2

Definition at line 15 of file SPIgrid.cpp.

◆ imag()

SPIVec SPI::imag ( const SPIVec A)

return the imaginary part of the vector

Parameters
[in]A[in] vector to take imaginary part of

Definition at line 638 of file SPIVec.cpp.

◆ integrate() [1/2]

PetscScalar SPI::integrate ( const SPIVec a,
SPIgrid1D grid 
)

integrate a vector of chebyshev Coefficients on a physical grid

Parameters
[in]a[in] vector to integrate over grid (Chebyshev coefficients or physical values)
[in]grid[in] respective grid

Definition at line 873 of file SPIgrid.cpp.

◆ integrate() [2/2]

PetscScalar SPI::integrate ( const SPIVec a,
SPIgrid2D grid 
)

integrate a vector of chebyshev Coefficients on a physical grid

Parameters
[in]a[in] vector to integrate over grid (Chebyshev coefficients or physical values)
[in]grid[in] respective grid

Definition at line 927 of file SPIgrid.cpp.

◆ integrate_coeffs() [1/2]

PetscScalar SPI::integrate_coeffs ( const SPIVec a)

integrate a vector of chebyshev Coefficients on a grid

Parameters
[in]a[in] chebyshev coefficients to integrate

Definition at line 813 of file SPIVec.cpp.

◆ integrate_coeffs() [2/2]

PetscScalar SPI::integrate_coeffs ( const SPIVec a,
const SPIVec y 
)

integrate a vector of chebyshev Coefficients on a stretched grid

Parameters
[in]a[in] chebyshev coefficients
[in]y[in] Cheby_mapped_y grid

Definition at line 846 of file SPIVec.cpp.

◆ interp1D_Mat() [1/2]

SPIMat SPI::interp1D_Mat ( SPIgrid1D grid1,
SPIgrid1D grid2 
)
Parameters
[in]grid1[in] grid to interpolate values from
[in]grid2[in] grid to interpolate values to

Definition at line 1506 of file SPIgrid.cpp.

◆ interp1D_Mat() [2/2]

SPIMat SPI::interp1D_Mat ( SPIVec y1,
SPIVec y2 
)
Parameters
[in]y1[in] grid to interpolate values from
[in]y2[in] grid to interpolate values to

Definition at line 1538 of file SPIgrid.cpp.

◆ inv()

SPIMat SPI::inv ( const SPIMat A)

Solve linear system A*Ainv=B using MatMatSolve.

Returns
Ainv matrix
Parameters
[in]A[in] A in A A^-1 = B

Definition at line 603 of file SPIMat.cpp.

◆ kron()

SPIMat SPI::kron ( const SPIMat A,
const SPIMat B 
)

set kronecker inner product of two matrices

Returns
kronecker inner product of the two matrices
Parameters
[in]A[in] A in A kron B operation
[in]B[in] B in A kron B operation

Definition at line 780 of file SPIMat.cpp.

◆ L2() [1/2]

PetscReal SPI::L2 ( const SPIVec  x1,
NormType  type 
)

calculate the \( L_2 \) norm of the vector

Returns
\(L_2\) norm of the vector
Parameters
[in]x1[in] \(x_1\) l
[in]type[in] type of norm (default NORM_2 \(\sqrt{\sum x_1^2}\))

Definition at line 871 of file SPIVec.cpp.

◆ L2() [2/2]

PetscReal SPI::L2 ( SPIVec  x1,
const SPIVec  x2,
NormType  type 
)

calculate the \( L_2 \) norm of the difference between \(x_1\) and \(x_2\) vectors.

Returns
\(L_2\) norm of the difference
Parameters
[in]x1[in] \(x_1\)
[in]x2[in] \(x_2\)
[in]type[in] type of norm (default NORM_2 \(\sqrt{\sum |x_1 - x_2|^2}\)) (NORM_1 denotes sum_i |x_i|), (NORM_2 denotes sqrt(sum_i |x_i|^2)), (NORM_INFINITY denotes max_i |x_i|)

Definition at line 859 of file SPIVec.cpp.

◆ linspace()

SPIVec SPI::linspace ( const PetscScalar  begin,
const PetscScalar  end,
const PetscInt  _rows 
)

return linspace of number of rows equally spaced points between begin and end

Parameters
[in]begin[in] beginning scalar of equally spaced points
[in]end[in] end scalar of equally spaced points
[in]_rows[in] how many points in array

Definition at line 645 of file SPIVec.cpp.

◆ load() [1/3]

PetscInt SPI::load ( SPIMat A,
const std::string  filename 
)

load matrix from filename from binary format (works with save(SPIMat,std::string) function

Returns
0 if successful
Parameters
[in,out]A[inout] A to load data into (must be initialized to the right size)
[in]filename[in] filename to read

Definition at line 1966 of file SPIMat.cpp.

◆ load() [2/3]

PetscInt SPI::load ( SPIVec A,
const std::string  filename 
)

load A from hdf5 filename using variable A.name, be sure it has the right size first before loading

Returns
0 if successful
Parameters
[in,out]A[inout] vector to load data into (must be initialized to the right size)
[in]filename[in] filename to read

Definition at line 591 of file SPIVec.cpp.

◆ load() [3/3]

PetscInt SPI::load ( std::vector< SPIMat > &  As,
const std::string  filename 
)

load matrix from filename from binary format (works with save(SPIMat,std::string) function

Returns
0 if successful
Parameters
[in,out]As[inout] matrices to load data into (must be initialized to the right size)
[in]filename[in] filename to read

Definition at line 1978 of file SPIMat.cpp.

◆ log()

SPIVec SPI::log ( const SPIVec A)

take the log (natural log) of each element in a vector

Definition at line 719 of file SPIVec.cpp.

◆ log10()

SPIVec SPI::log10 ( const SPIVec A)

take the log10 of each element in a vector

Definition at line 721 of file SPIVec.cpp.

◆ LST_spatial()

std::tuple< PetscScalar, SPIVec > SPI::LST_spatial ( SPIparams params,
SPIgrid1D grid,
SPIbaseflow baseflow,
SPIVec  q 
)

solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(omega,eig_vector) = LST_spatial(params,grid,baseflow). Will solve for closest eigenvalue to params.omega
Parameters
[in,out]params[inout] contains parameters including Re, alpha, and omega. Will overwrite omega with true omega value once solved
[in]grid[in] grid class containing the grid location and respective derivatives
[in]baseflow[in] baseflow for parallel flow
[in]q[in] initial guess for temporal problem

Definition at line 147 of file SPILST.cpp.

◆ LST_spatial_cg()

std::tuple< PetscScalar, PetscScalar, SPIVec, SPIVec > SPI::LST_spatial_cg ( SPIparams params,
SPIgrid1D grid,
SPIbaseflow baseflow 
)

solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(omega,eig_vector) = LST_spatial(params,grid,baseflow). Will solve for closest eigenvalue to params.omega
Parameters
[in,out]params[inout] contains parameters including Re, alpha, and omega. Will overwrite omega with true omega value once solved
[in]grid[in] grid class containing the grid location and respective derivatives
[in]baseflow[in] baseflow for parallel flow

Definition at line 349 of file SPILST.cpp.

◆ LST_temporal()

std::tuple< PetscScalar, SPIVec > SPI::LST_temporal ( SPIparams params,
SPIgrid1D grid,
SPIbaseflow baseflow,
SPIVec  q 
)

solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with alpha being pure real, and omega the eigenvalue

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(omega,eig_vector) = LST_temporal(params,grid,baseflow). Will solve for closest eigenvalue to params.omega
Parameters
[in,out]params[inout] contains parameters including Re, alpha, and omega. Will overwrite omega with true omega value once solved
[in]grid[in] grid class containing the grid location and respective derivatives
[in]baseflow[in] baseflow for parallel flow
[in]q[in] initial guess for temporal problem

Definition at line 5 of file SPILST.cpp.

◆ LSTNP_spatial()

std::tuple< PetscScalar, PetscScalar, SPIVec, SPIVec > SPI::LSTNP_spatial ( SPIparams params,
SPIgrid1D grid,
SPIbaseflow baseflow,
SPIVec  ql,
SPIVec  qr 
)

solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha,group_velocity,left_eig_vector,right_eig_vector) = LSTNP_spatial(params,grid,baseflow). Will solve for closest eigenvalue to params.alpha
Parameters
[in,out]params[inout] contains parameters including Re, alpha, and omega. Will overwrite omega with true omega value once solved
[in]grid[in] grid class containing the grid location and respective derivatives
[in]baseflow[in] baseflow for parallel flow
[in]ql[in] initial guess for spatial problem (adjoint) for left eigenfunction
[in]qr[in] initial guess for spatial problem

Definition at line 715 of file SPILST.cpp.

◆ LSTNP_spatial_right()

std::tuple< PetscScalar, SPIVec > SPI::LSTNP_spatial_right ( SPIparams params,
SPIgrid1D grid,
SPIbaseflow baseflow,
SPIVec  qr 
)

solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha,right_eig_vector) = LSTNP_spatial(params,grid,baseflow). Will solve for closest eigenvalue to params.alpha
Parameters
[in,out]params[inout] contains parameters including Re, alpha, and omega. Will overwrite omega with true omega value once solved
[in]grid[in] grid class containing the grid location and respective derivatives
[in]baseflow[in] baseflow for parallel flow
[in]qr[in] initial guess for spatial problem

Definition at line 1046 of file SPILST.cpp.

◆ LSTNP_spatial_right2()

std::tuple< PetscScalar, SPIVec > SPI::LSTNP_spatial_right2 ( SPIparams params,
SPIgrid1D grid,
SPIbaseflow baseflow,
SPIVec  qr 
)

solve the local stability theory problem for the linearized Navier-Stokes equations using parallel baseflow with omega being pure real, and alpha the eigenvalue

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha,right_eig_vector) = LSTNP_spatial(params,grid,baseflow). Will solve for closest eigenvalue to params.alpha
Parameters
[in,out]params[inout] contains parameters including Re, alpha, and omega. Will overwrite omega with true omega value once solved
[in]grid[in] grid class containing the grid location and respective derivatives
[in]baseflow[in] baseflow for parallel flow
[in]qr[in] initial guess for spatial problem

Definition at line 1276 of file SPILST.cpp.

◆ LSTNP_spatials_right()

std::tuple< std::vector< PetscScalar >, std::vector< SPIVec > > SPI::LSTNP_spatials_right ( SPIparams params,
SPIgrid1D grid,
SPIbaseflow baseflow,
std::vector< PetscScalar > &  alphas,
std::vector< SPIVec > &  qrs 
)

solve the local stability theory problem for the linearized Navier-Stokes equations using non-parallel baseflow with omega being pure real, and alpha the eigenvalue

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alphas,right_eig_vectors) = LSTNP_spatials_right(params,grid,baseflow). Will solve for closest eigenvalue to params.alpha
Parameters
[in,out]params[inout] contains parameters including Re and omega.
[in]grid[in] grid class containing the grid location and respective derivatives
[in]baseflow[in] baseflow for parallel flow
[in]alphas[in] vector of alpha guesses
[in]qrs[in] initial guesses for spatial problem

Definition at line 1482 of file SPILST.cpp.

◆ lstsq() [1/2]

SPIVec SPI::lstsq ( const SPIMat A,
SPIVec y 
)

solve least squares problem of A*x=y for skinny A matrix using SVD

Returns
x
Parameters
[in]A[in] A matrix in A*x=y least squares problem
[in]y[in] y in A*x=y least squares problem

Definition at line 905 of file SPIMat.cpp.

◆ lstsq() [2/2]

SPIVec SPI::lstsq ( const std::vector< SPIVec > &  A,
SPIVec y 
)

solve least squares problem of A*x=y for skinny A matrix (or vectors of columns vectors) using SVD

Returns
x
Parameters
[in]A[in] vectors of columns vectors making up A matrix
[in]y[in] y in A*x=y least squares problem

Definition at line 920 of file SPIMat.cpp.

◆ map_D()

SPIMat SPI::map_D ( SPIMat  D,
SPIVec  y,
PetscInt  d,
PetscInt  order 
)

map the derivative operator to the proper y grid

Returns
mapped derivative to streched grid
Parameters
[in]D[in] derivative operator on uniform grid from 0 to 1 (xi grid)
[in]y[in] potentially non-uniform grid from 0 to y_max
[in]d[in] order of derivative
[in]order[in] order of accuracy for finite difference stencil (default 4)

Definition at line 52 of file SPIgrid.cpp.

◆ map_D_Chebyshev()

SPIMat SPI::map_D_Chebyshev ( SPIVec x,
PetscInt  d 
)

map a Chebyshev collocated operator acting with respect to the stretched collocated grid

Returns
Chebyshev collocated derivative operator
Parameters
[in]x[in] grid (created from set_Cheby_stretched_y or set_Cheby_y or similar example)
[in]d[in] order of the derivative (default 1)

Definition at line 278 of file SPIgrid.cpp.

◆ meshgrid()

std::tuple< SPIMat, SPIMat > SPI::meshgrid ( SPIVec x,
SPIVec y 
)
Parameters
[in]x[in] x input array
[in]y[in] y input array

Definition at line 1894 of file SPIMat.cpp.

◆ ones()

SPIVec SPI::ones ( const PetscInt  rows)

create a vector of size rows full of ones

Returns
vector of size rows
Parameters
[in]rows[in] number of rows for vector size

Definition at line 605 of file SPIVec.cpp.

◆ operator%()

SPIMat SPI::operator% ( const SPIMat A,
const SPIMat B 
)

take the elementwise multiplication of each element in a matrix

Definition at line 2033 of file SPIMat.cpp.

◆ operator*() [1/3]

SPIMat SPI::operator* ( const PetscScalar  a,
const SPIMat  A 
)

a*A operation to be equivalent to A*a

Returns
new matrix of a*A
Parameters
[in]a[in] scalar a in a*A
[in]A[in] matrix A in a*A

Definition at line 517 of file SPIMat.cpp.

◆ operator*() [2/3]

SPIVec SPI::operator* ( const PetscScalar  a,
const SPIVec Y 
)

Z=a*Y operation to be equivalent to Y*a

Returns
Z in Z=a*Y operation
Parameters
[in]a[in] scalar a in a*Y operation
[in]Y[in] Y in Z=a*Y

Definition at line 510 of file SPIVec.cpp.

◆ operator*() [3/3]

SPIMat SPI::operator* ( const SPIMat  A,
const PetscScalar  a 
)

A*a operation to be equivalent to A*a.

Returns
new matrix of A*a operation

Definition at line 527 of file SPIMat.cpp.

◆ operator+()

SPIVec SPI::operator+ ( const PetscScalar  a,
const SPIVec Y 
)

Z=a+Y operation to be equivalent to Y+a

Returns
Z in Z=a+Y operation
Parameters
[in]a[in] scalar a in a+Y operation
[in]Y[in] Y in a+Y operation

Definition at line 521 of file SPIVec.cpp.

◆ operator-()

SPIVec SPI::operator- ( const PetscScalar  a,
const SPIVec Y 
)

Z=a-Y operation to be equivalent to Y-a

Returns
Z in Z=a-Y operation
Parameters
[in]a[in] scalar a in a-Y operation
[in]Y[in] Y in a-Y operation

Definition at line 532 of file SPIVec.cpp.

◆ operator/() [1/2]

SPIVec SPI::operator/ ( const PetscScalar  a,
const SPIVec Y 
)

Z=a/Y operation

Returns
Z in Z=a/Y operation
Parameters
[in]a[in] scalar a in a/Y operation
[in]Y[in] Y in Z=a/Y

Definition at line 497 of file SPIVec.cpp.

◆ operator/() [2/2]

SPIVec SPI::operator/ ( const SPIVec b,
const SPIMat A 
)

Solve linear system, Ax=b using b/A notation.

Returns
x vector in Ax=b solved linear system
Parameters
[in]b[in] b in Ax=b
[in]A[in] A in Ax=b

Definition at line 535 of file SPIMat.cpp.

◆ operator^()

SPIMat SPI::operator^ ( const PetscScalar  a,
const SPIMat A 
)

Y=a^A operation.

Returns
Y matrix in Y=a^A
Parameters
[in]a[in] a in Y=a^A
[in]A[in] A in Y=a^A

Definition at line 578 of file SPIMat.cpp.

◆ orthogonalize() [1/3]

SPIMat SPI::orthogonalize ( const std::vector< SPIVec > &  x)
Parameters
[in]x[in] array of vectors to orthogonalize

Definition at line 2057 of file SPIMat.cpp.

◆ orthogonalize() [2/3]

std::vector< SPIVec > SPI::orthogonalize ( std::vector< SPIVec > &  x,
SPIgrid1D grid 
)
Parameters
[in]x[in] array of vectors to orthogonalize
[in]grid[in] respective grid for integration using Gram-Schmidt

Definition at line 1323 of file SPIgrid.cpp.

◆ orthogonalize() [3/3]

std::vector< SPIVec > SPI::orthogonalize ( std::vector< SPIVec > &  x,
SPIgrid2D grid 
)
Parameters
[in]x[in] array of vectors to orthogonalize
[in]grid[in] respective grid for integration using Gram-Schmidt

Definition at line 1406 of file SPIgrid.cpp.

◆ polyeig()

std::tuple< PetscScalar, SPIVec > SPI::polyeig ( const std::vector< SPIMat > &  As,
const PetscScalar  target,
const PetscReal  tol,
const PetscInt  max_iter 
)

solve general polynomial eigenvalue problem of (A0 + A1*alpha + A2*alpha^2 + ...)*x = 0 and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector)

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha, eig_vector) = eig([A0,A1,A2...],0.1+0.4*PETSC_i)
Parameters
[in]As[in] [A0,A1,A2...] in generalized polynomial eigenvalue problem
[in]target[in] target eigenvalue to solve for
[in]tol[in] tolerance of eigenvalue solver
[in]max_iter[in] maximum number of iterations

Definition at line 1686 of file SPIMat.cpp.

◆ polyeig_init()

std::tuple< PetscScalar, SPIVec > SPI::polyeig_init ( const std::vector< SPIMat > &  As,
const PetscScalar  target,
const SPIVec qr,
const PetscReal  tol,
const PetscInt  max_iter 
)

solve general polynomial eigenvalue problem of (A0 + A1*alpha + A2*alpha^2 + ...)*x = 0 and return a tuple of tie(PetscScalar alpha, SPIVec eig_vector) using initial subspace vector

Returns
tuple of eigenvalue and eigenvector closest to the target value e.g. std::tie(alpha, eig_vector) = eig([A0,A1,A2...],0.1+0.4*PETSC_i)
Parameters
[in]As[in] [A0,A1,A2...] in generalized polynomial eigenvalue problem
[in]target[in] target eigenvalue to solve for
[in]qr[in] initial subspace vector for right eigenvalue problem
[in]tol[in] tolerance of eigenvalue solver
[in]max_iter[in] maximum number of iterations

Definition at line 1754 of file SPIMat.cpp.

◆ pow() [1/2]

SPIVec SPI::pow ( const SPIVec A,
PetscScalar  b 
)

take the pow of each element in the vector (A^b)

Returns
A^b
Parameters
[in]A[in] vector to raise to the power
[in]b[in] the exponenet

Definition at line 777 of file SPIVec.cpp.

◆ pow() [2/2]

SPIVec SPI::pow ( const SPIVec A,
SPIVec B 
)

take the pow of each element in the vectors

Definition at line 775 of file SPIVec.cpp.

◆ printf()

PetscInt SPI::printf ( std::string  msg,
  ... 
)

print a message to string using PetscPrintf (also adds a newline at end) (note: only prints on rank 0 processor)

Returns
0 if successful
Parameters
[in]msg[in] message to print with formatting such as %g
[in]...[in] scalar or double or int to match formatting string to be output on processor with rank 0

Definition at line 5 of file SPIprint.cpp.

◆ printfc()

PetscInt SPI::printfc ( std::string  msg,
PetscScalar  val 
)

print a message to string using PetscPrintf (also adds a newline at end) (note: only prints on rank 0 processor) with PetscScalars as input and two formats per argument

Returns
0 if successful
Parameters
[in]msg[in] message to print with formatting such as %g
[in]val[in] scalar or double or int to match formatting string to be output on processor with rank 0

Definition at line 29 of file SPIprint.cpp.

◆ proj() [1/2]

SPIVec SPI::proj ( SPIVec u,
SPIVec v,
SPIgrid1D grid 
)
Parameters
[in]u[in] first vector to project
[in]v[in] second vector to project
[in]grid[in] respective grid

Definition at line 1220 of file SPIgrid.cpp.

◆ proj() [2/2]

SPIVec SPI::proj ( SPIVec u,
SPIVec v,
SPIgrid2D grid 
)
Parameters
[in]u[in] first vector to project
[in]v[in] second vector to project
[in]grid[in] respective grid

Definition at line 1267 of file SPIgrid.cpp.

◆ real()

SPIVec SPI::real ( const SPIVec A)

return the real part of the vector

Parameters
[in]A[in] vector to take real part of

Definition at line 631 of file SPIVec.cpp.

◆ save() [1/5]

PetscInt SPI::save ( const SPIMat A,
const std::string  filename 
)

save matrix to filename in binary format (see Petsc documentation for format ) Format is (from Petsc documentation): int MAT_FILE_CLASSID int number of rows int number of columns int total number of nonzeros int *number nonzeros in each row int *column indices of all nonzeros (starting index is zero) PetscScalar *values of all nonzeros

Returns
0 if successful
Parameters
[in]A[in] A to save in
[in]filename[in] filename to save data to

Definition at line 1924 of file SPIMat.cpp.

◆ save() [2/5]

PetscInt SPI::save ( const SPIVec A,
const std::string  filename 
)

save A to hdf5 to filename as variable A.name (note: this will append if filename already exists)

Returns
0 if successful
Parameters
[in]A[in] A to save in hdf5 format under A.name variable
[in]filename[in] filename to save

Definition at line 543 of file SPIVec.cpp.

◆ save() [3/5]

PetscInt SPI::save ( const std::vector< SPIMat > &  As,
const std::string  filename 
)

save matrices to filename in binary format (see Petsc documentation for format

Returns
0 if successful
Parameters
[in]As[in] A to save in
[in]filename[in] filename to save data to

Definition at line 1944 of file SPIMat.cpp.

◆ save() [4/5]

PetscInt SPI::save ( std::vector< PetscScalar >  A,
std::string  variablename,
std::string  filename 
)

save A to hdf5 to filename as variable A.name (note: this will append if filename already exists)

Returns
0 if successful
Parameters
[in]A[in] variable to save to hdf5 file
[in]variablename[in] the name of the variable to save in the file
[in]filename[in] hdf5 filename

Definition at line 562 of file SPIVec.cpp.

◆ save() [5/5]

PetscInt SPI::save ( std::vector< SPIVec A,
std::string  variablename,
std::string  filename 
)

save A to hdf5 to filename as variable A.name (note: this will append if filename already exists)

Returns
0 if successful
Parameters
[in]A[in] variable to save to hdf5 file
[in]variablename[in] the name of the variable to save in the file
[in]filename[in] hdf5 filename

Definition at line 576 of file SPIVec.cpp.

◆ set_Cheby_mapped_y()

SPIVec SPI::set_Cheby_mapped_y ( PetscScalar  a,
PetscScalar  b,
PetscInt  ny 
)

create a stretched Chebyshev grid from [a,b] using default Chebfun like mapping

Returns
grid in [a,b] using Chebyshev collocated points
Parameters
[in]a[in] lower bound of grid domain [a,b]
[in]b[in] upper bound of grid domain [a,b]
[in]ny[in] number of points on grid

Definition at line 318 of file SPIgrid.cpp.

◆ set_Cheby_stretched_y()

SPIVec SPI::set_Cheby_stretched_y ( PetscScalar  y_max,
PetscInt  ny,
PetscScalar  yi 
)

create a stretched Chebyshev grid from [0,y_max] with yi being the midpoint location for the number of points

Returns
stretched grid using Chebyshev collocated points
Parameters
[in]y_max[in] maximum range for [0, y_max] stretching
[in]ny[in] number of points
[in]yi[in] midpoint location for the number of points (i.e. ny/2 are above and below this point)

Definition at line 305 of file SPIgrid.cpp.

◆ set_Cheby_y()

SPIVec SPI::set_Cheby_y ( PetscInt  ny)

create a Chebyshev grid from [-1,1]

Returns
grid using Chebyshev collocated points
Parameters
[in]ny[in] number of points

Definition at line 329 of file SPIgrid.cpp.

◆ set_D()

SPIMat SPI::set_D ( SPIVec y,
PetscInt  d,
PetscInt  order,
PetscBool  uniform 
)

set the derivative operator for the proper y grid if uniform=false. Uses map_D function

Returns
mapped derivative to streched grid
Parameters
[in]y[in] grid points
[in]d[in] order of derivative
[in]order[in] order of accuracy of derivative (default 4)
[in]uniform[in] is this for a uniform grid? (default false)

Definition at line 82 of file SPIgrid.cpp.

◆ set_D_Chebyshev()

SPIMat SPI::set_D_Chebyshev ( SPIVec x,
PetscInt  d,
PetscBool  need_map 
)

set a Chebyshev collocated operator acting with respect to the collocated grid

Returns
Chebyshev collocated derivative operator
Parameters
[in]x[in] grid (created from set_Cheby_stretched_y or set_Cheby_y or similar example)
[in]d[in] order of the derivative (default 1)
[in]need_map[in] need mapping? (default false)

Definition at line 225 of file SPIgrid.cpp.

◆ set_D_Fourier()

SPIMat SPI::set_D_Fourier ( SPIVec  t,
PetscInt  d 
)

create a Fourier derivative operator acting on grid t

Returns
derivative operator using Fourier collocated points
Parameters
[in]t[in] grid point created from set_Fourier_t
[in]d[in] order of derivatives

Definition at line 407 of file SPIgrid.cpp.

◆ set_D_periodic()

SPIMat SPI::set_D_periodic ( SPIVec y,
PetscInt  d,
PetscInt  order 
)

set the derivative operator for the proper periodic grid assuming uniform discretization.

Returns
derivative operator for uniform periodic grid
Parameters
[in]y[in] grid points
[in]d[in] order of derivative
[in]order[in] order of accuracy of derivative (default 4)

Definition at line 172 of file SPIgrid.cpp.

◆ set_D_UltraS()

std::tuple< SPIMat, SPIMat > SPI::set_D_UltraS ( SPIVec x,
PetscInt  d 
)

set a UltraSpherical operator acting with respect to the collocated grid (keeps everything in UltraSpherical space), take in Chebyshev coefficients and outputs coefficients in C(d) coefficient space

Returns
UltraSpherical derivative operator
Parameters
[in]x[in] grid (must be created from set_Cheby_y or set_Cheby_mapped_y)
[in]d[in] order of the derivative (default 1)

Definition at line 338 of file SPIgrid.cpp.

◆ set_FD_stretched_y()

SPIVec SPI::set_FD_stretched_y ( PetscScalar  y_max,
PetscInt  ny,
PetscScalar  delta 
)

set stretched grid from [0,y_max] using tanh stretching for use with finite difference operators

Returns
y stretched grid
Parameters
[in]y_max[in] upper bound of y in [0,y_max]
[in]ny[in] number of points to use in the domain
[in]delta[in] stretching parameter, near zero is no stretching, (default 2.0001)

Definition at line 214 of file SPIgrid.cpp.

◆ set_Fourier_t()

SPIVec SPI::set_Fourier_t ( PetscScalar  T,
PetscInt  nt 
)

create a Fourier grid from [0,T]

Returns
grid using Fourier collocated points without the last point of linspace(0,T,n+1)[:-1]
Parameters
[in]T[in] period or end point for [0,T] domain
[in]nt[in] number of points

Definition at line 392 of file SPIgrid.cpp.

◆ set_T_That()

std::tuple< SPIMat, SPIMat > SPI::set_T_That ( PetscInt  n)

set a T and That operators acting with respect to the Chebyshev collocated grid, That: physical -> Chebyshev coefficient, T: Chebyshev coefficient -> physical

Returns
UltraSpherical derivative operator
Parameters
[in]n[in] number of discretized points in Chebyshev grid

Definition at line 367 of file SPIgrid.cpp.

◆ sin() [1/2]

SPIMat SPI::sin ( const SPIMat A)

take the sin of each element in a matrix

Definition at line 2025 of file SPIMat.cpp.

◆ sin() [2/2]

SPIVec SPI::sin ( const SPIVec A)

take the sin of each element in a vector

Definition at line 706 of file SPIVec.cpp.

◆ sinh()

SPIVec SPI::sinh ( const SPIVec A)

take the sinh of each element in a vector

Definition at line 723 of file SPIVec.cpp.

◆ solve()

SPIVec SPI::solve ( const SPIMat A,
const SPIVec b 
)

Solve linear system, Ax=b using solve(A,b) notation.

Returns
x vector in Ax=b solved linear system
Parameters
[in]A[in] A in Ax=b
[in]b[in] b in Ax=b

Definition at line 596 of file SPIMat.cpp.

◆ SPIVec1Dto2D()

SPIVec SPI::SPIVec1Dto2D ( SPIgrid2D grid2D,
SPIVec u 
)

expand a 1D vector to a 2D vector copying data along time dimension

Parameters
[in]grid2D[in] 3D grid containing wall-normal and time dimensions
[in]u[in] vector to inflate to match the 3D operator

Definition at line 850 of file SPIgrid.cpp.

◆ sqrt()

SPIVec SPI::sqrt ( const SPIVec A)

take the atanh of each element in a vector

Definition at line 741 of file SPIVec.cpp.

◆ sum()

PetscScalar SPI::sum ( SPIVec  x1)

take the sum of a vector

Parameters
[in]x1[in] vector to sum

Definition at line 805 of file SPIVec.cpp.

◆ svd()

std::tuple< std::vector< PetscReal >, std::vector< SPIVec >, std::vector< SPIVec > > SPI::svd ( const SPIMat A)

solve SVD problem of A=U*E*V^H for skinny A matrix

Returns
tuple of sigma,u,v std::tie(sigma, u, v) = svd(A)
Parameters
[in]A[in] A to perform SVD decomposition A=U*E*V^H problem

Definition at line 851 of file SPIMat.cpp.

◆ tan() [1/2]

SPIMat SPI::tan ( const SPIMat A)

take the tan of each element in a matrix

Definition at line 2031 of file SPIMat.cpp.

◆ tan() [2/2]

SPIVec SPI::tan ( const SPIVec A)

take the tan of each element in a vector

Definition at line 710 of file SPIVec.cpp.

◆ tanh()

SPIVec SPI::tanh ( const SPIVec A)

take the tanh of each element in a vector

Definition at line 727 of file SPIVec.cpp.

◆ trapz() [1/2]

PetscScalar SPI::trapz ( SPIVec  y)

trapezoidal integration of y with unity coordinate spacing, \(\int y dx \)

Returns
integrated value
Parameters
[in]y[in] vector to integrate, assuming default spacing of one

Definition at line 898 of file SPIVec.cpp.

◆ trapz() [2/2]

PetscScalar SPI::trapz ( SPIVec  y,
SPIVec  x 
)

trapezoidal integration of y with x coordinates, \(\int y dx \)

Returns
integrated value
Parameters
[in]y[in] vector to integrate
[in]x[in] optional, coordinates to integrate over, must be same size as y, and defaults to spacing of one if not given

Definition at line 915 of file SPIVec.cpp.

◆ zeros() [1/2]

SPIMat SPI::zeros ( const PetscInt  m,
const PetscInt  n 
)

create, form, and return zeros matrix of size mxn

Returns
zero matrix of size mxn
Parameters
[in]m[in] m size of zero matrix
[in]n[in] n size of zero matrix

Definition at line 708 of file SPIMat.cpp.

◆ zeros() [2/2]

SPIVec SPI::zeros ( const PetscInt  rows)

create and return a vector of size rows full of zeros

Returns
vector of size zeros
Parameters
[in]rows[in] size of vector to create

Definition at line 614 of file SPIVec.cpp.