:orphan: # PetscDTAltVWedge Compute the wedge product of a j-form and a k-form, giving a (j+k) form ## Synopsis ``` #include "petscdt.h" PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb) ``` ## Input Parameters - ***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-form b, 0 <= k <= N and 0 <= j+k <= N - ***a -*** a j-form, size [N choose j] - ***b -*** a k-form, size [N choose k] ## Output Parameter - ***awedgeb -*** the (j+k)-form a wedge b, size [N choose (j+k)]: (a wedge b)(v_1,...,v_{j+k}) = \sum_{s} sign(s) a(v_{s_1},...,v_{s_j}) b(v_{s_{j+1}},...,v_{s_{j+k}}), where the sum is over permutations s such that s_1 < s_2 < ... < s_j and s_{j+1} < s_{j+2} < ... < s_{j+k}. ## See Also `PetscDTAltV`, `PetscDTAltVWedgeMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` ## Level intermediate ## Location src/dm/dt/interface/dtaltv.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/dt/interface/dtaltv.c) [Index of all DT routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)