:orphan: # PCApplyTranspose Applies the transpose of preconditioner to a vector. ## Synopsis ``` #include "petscksp.h" PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y) ``` Collective ## Input Parameters - ***pc -*** the preconditioner context - ***x -*** input vector ## Output Parameter - ***y -*** output vector ## Note For complex numbers this applies the non-Hermitian transpose. ## Developer Note We need to implement a `PCApplyHermitianTranspose()` ## See Also `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()` ## Level developer ## Location src/ksp/pc/interface/precon.c ## Implementations PCApplyTranspose_ASM in src/ksp/pc/impls/asm/asm.c
PCApplyTranspose_BDDC in src/ksp/pc/impls/bddc/bddc.c
PCApplyTranspose_BDDCIPC in src/ksp/pc/impls/bddc/bddc.c
PCApplyTranspose_Eisenstat in src/ksp/pc/impls/eisens/eisen.c
PCApplyTranspose_Cholesky in src/ksp/pc/impls/factor/cholesky/cholesky.c
PCApplyTranspose_ILU in src/ksp/pc/impls/factor/ilu/ilu.c
PCApplyTranspose_LU in src/ksp/pc/impls/factor/lu/lu.c
PCApplyTranspose_QR in src/ksp/pc/impls/factor/qr/qr.c
PCApplyTranspose_FieldSplit in src/ksp/pc/impls/fieldsplit/fieldsplit.c
PCApplyTranspose_GASM in src/ksp/pc/impls/gasm/gasm.c
PCApplyTranspose_H2OPUS in src/ksp/pc/impls/h2opus/pch2opus.c
PCApplyTranspose_KSP in src/ksp/pc/impls/ksp/pcksp.c
PCApplyTranspose_Mat in src/ksp/pc/impls/mat/pcmat.c
PCApplyTranspose_MG in src/ksp/pc/impls/mg/mg.c
PCApplyTranspose_PBJacobi in src/ksp/pc/impls/pbjacobi/pbjacobi.c
PCApplyTranspose_Redundant in src/ksp/pc/impls/redundant/redundant.c
PCApplyTranspose_Shell in src/ksp/pc/impls/shell/shellpc.c
PCApplyTranspose_SOR in src/ksp/pc/impls/sor/sor.c
PCApplyTranspose_SVD in src/ksp/pc/impls/svd/svd.c
PCApplyTranspose_VPBJacobi in src/ksp/pc/impls/vpbjacobi/vpbjacobi.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/pc/interface/precon.c) [Index of all PC routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)