CeedBasis

static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col)

Compute Householder reflection.

Computes A = (I - b v v^T) A where A is an mxn matrix indexed as A[i*row + j*col]

Library Developer Functions

Parameters
  • A[inout] Matrix to apply Householder reflection to, in place

  • v – Householder vector

  • b – Scaling factor

  • m – Number of rows in A

  • n – Number of columns in A

  • row – Row stride

  • col – Col stride

Returns

An error code: 0 - success, otherwise - failure

int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, CeedInt k, CeedInt row, CeedInt col)

Apply Householder Q matrix.

Compute A = Q A where Q is mxm and A is mxn.

Library Developer Functions

Parameters
  • A[inout] Matrix to apply Householder Q to, in place

  • Q – Householder Q matrix

  • tau – Householder scaling factors

  • t_mode – Transpose mode for application

  • m – Number of rows in A

  • n – Number of columns in A

  • k – Number of elementary reflectors in Q, k<m

  • row – Row stride in A

  • col – Col stride in A

Returns

An error code: 0 - success, otherwise - failure

static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n)

Compute Givens rotation.

Computes A = G A (or G^T A in transpose mode) where A is an mxn matrix indexed as A[i*n + j*m]

Library Developer Functions

Parameters
  • A[inout] Row major matrix to apply Givens rotation to, in place

  • c – Cosine factor

  • s – Sine factor

  • t_modeCEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of A clockwise; CEED_TRANSPOSE for the opposite rotation

  • i – First row/column to apply rotation

  • k – Second row/column to apply rotation

  • m – Number of rows in A

  • n – Number of columns in A

Returns

An error code: 0 - success, otherwise - failure

static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream)

View an array stored in a CeedBasis.

Library Developer Functions

Parameters
  • name[in] Name of array

  • fp_fmt[in] Printing format

  • m[in] Number of rows in array

  • n[in] Number of columns in array

  • a[in] Array to be viewed

  • stream[in] Stream to view to, e.g., stdout

Returns

An error code: 0 - success, otherwise - failure

static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project)

Create the interpolation and gradient matrices for projection from the nodes of basis_from to the nodes of basis_to.

The interpolation is given by interp_project = interp_to^+ * interp_from, where the pesudoinverse interp_to^+ is given by QR factorization. The gradient is given by grad_project = interp_to^+ * grad_from. Note: basis_from and basis_to must have compatible quadrature spaces.

Library Developer Functions

Parameters
  • basis_from[in] CeedBasis to project from

  • basis_to[in] CeedBasis to project to

  • interp_project[out] Address of the variable where the newly created interpolation matrix will be stored.

  • grad_project[out] Address of the variable where the newly created gradient matrix will be stored.

Returns

An error code: 0 - success, otherwise - failure