Actual source code: petsccharacteristic.h

  1: /*
  2:    Defines the interface functions for the method of characteristics solvers
  3: */
  4: #pragma once

  6: #include <petscvec.h>
  7: #include <petscdmdatypes.h>

  9: /* SUBMANSEC = Characteristic */

 11: PETSC_EXTERN PetscErrorCode CharacteristicInitializePackage(void);
 12: PETSC_EXTERN PetscErrorCode CharacteristicFinalizePackage(void);

 14: /*S
 15:      Characteristic - Abstract PETSc object that manages method of characteristics solves

 17:    Level: beginner

 19: .seealso: `CharacteristicCreate()`, `CharacteristicSetType()`, `CharacteristicType`, `SNES`, `TS`, `PC`, `KSP`
 20: S*/
 21: typedef struct _p_Characteristic *Characteristic;

 23: /*J
 24:     CharacteristicType - String with the name of a characteristics method.

 26:    Level: beginner

 28: .seealso: `CharacteristicSetType()`, `Characteristic`
 29: J*/
 30: #define CHARACTERISTICDA "da"
 31: typedef const char *CharacteristicType;

 33: PETSC_EXTERN PetscErrorCode CharacteristicCreate(MPI_Comm, Characteristic *);
 34: PETSC_EXTERN PetscErrorCode CharacteristicSetType(Characteristic, CharacteristicType);
 35: PETSC_EXTERN PetscErrorCode CharacteristicSetUp(Characteristic);
 36: PETSC_EXTERN PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic, DM, Vec, Vec, PetscInt, PetscInt[], PetscErrorCode (*)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *);
 37: PETSC_EXTERN PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic, DM, Vec, Vec, PetscInt, PetscInt[], PetscErrorCode (*)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *);
 38: PETSC_EXTERN PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic, DM, Vec, PetscInt, PetscInt[], PetscErrorCode (*)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *);
 39: PETSC_EXTERN PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic, DM, Vec, PetscInt, PetscInt[], PetscErrorCode (*)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *);
 40: PETSC_EXTERN PetscErrorCode CharacteristicSolve(Characteristic, PetscReal, Vec);
 41: PETSC_EXTERN PetscErrorCode CharacteristicDestroy(Characteristic *);

 43: PETSC_EXTERN PetscFunctionList CharacteristicList;

 45: PETSC_EXTERN PetscErrorCode CharacteristicRegister(const char[], PetscErrorCode (*)(Characteristic));