:orphan: # PetscDTAltVPullbackMatrix Compute the pullback matrix for k-forms under a linear transformation ## Synopsis ``` #include "petscdt.h" PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar) ``` ## Input Parameters - ***N -*** the dimension of the origin vector space of the linear transformation, N >= 0 - ***M -*** the dimension of the image vector space of the linear transformation, M >= 0 - ***L -*** a linear transformation, an [M x N] matrix in row-major format - ***k -*** the *signed* degree k of the |k|-forms on which Lstar acts, -(min(M,N)) <= k <= min(M,N). A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note in `PetscDTAltvPullback()`) ## Output Parameter - ***Lstar -*** the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w ## See Also `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVStar()` ## 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)