PetscDTAltVStar#

Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form

Synopsis#

#include "petscdt.h" 
PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw)

Input Parameters#

  • N - the dimension of the vector space, N >= 0

  • k - the degree k of the k-form w, 0 <= k <= N

  • pow - the number of times to apply the Hodge star operator: pow < 0 indicates that the inverse of the Hodge star operator should be applied |pow| times.

  • w - a k-form, size [N choose k]

Output Parameter#

  • starw = (star)^pow w. Each degree of freedom of a k- form is associated with a subset S of k coordinates of the N dimensional vector space: the Hodge start operator (star) maps that degree of freedom to the degree of freedom associated with S’, the complement of S, with a sign change if the permutation of coordinates {S[0], … S[k-1], S’[0], … S’[N-k- 1]} is an odd permutation. This implies (star)^2 w = (-1)^{k(N-k)} w, and (star)^4 w = w.

See Also#

PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()

Level#

intermediate

Location#

src/dm/dt/interface/dtaltv.c


Edit on GitLab

Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages