petsc-3.14.6 2021-03-30
Report Typos and Errors

PetscDTAltVWedgeMatrix

Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms

Synopsis

#include "petscdt.h" 
PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat)

Input Arguments

N - the dimension of the vector space, N >= 0
j - the degree j of the j-form a, 0 <= j <= N
k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N
a - a j-form, size [N choose j]

Output Arguments

awedge - (a wedge), an [(N choose j+k) x (N choose k)] matrix in row-major order, such that (a wedge) * b = a wedge b

See Also

PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()

Level

intermediate

Location

src/dm/dt/interface/dtaltv.c
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages