#include <DM_wrapper.hpp>
|
| DMWrapper () |
|
| ~DMWrapper () |
|
| DMWrapper (const DMWrapper &original) |
|
DMWrapper & | operator= (const DMWrapper &original) |
|
PetscErrorCode | initialize_pseudo_inverse_interpolation_petsc (int nthreads, PetscInt order, const PetscInt Nspecies, const PetscReal vPara_min, const PetscReal vPara_max, const PetscInt nvp, const PetscReal vPerp_min, const PetscReal vPerp_max, const PetscInt nmu) |
|
PetscErrorCode | destroy_DM_petsc () |
|
void | initialize_pseudo_inverse_interpolation (int nthreads, int order, const int Nspecies, const double vPara_min, const double vPara_max, const int nvp, const double vPerp_min, const double vPerp_max, const int nmu) |
|
void | destroy_DM () |
|
◆ DMWrapper() [1/2]
◆ ~DMWrapper()
DMWrapper::~DMWrapper |
( |
| ) |
|
|
inline |
◆ DMWrapper() [2/2]
DMWrapper::DMWrapper |
( |
const DMWrapper & |
original | ) |
|
|
inline |
◆ destroy_DM()
void DMWrapper::destroy_DM |
( |
| ) |
|
|
inline |
◆ destroy_DM_petsc()
PetscErrorCode DMWrapper::destroy_DM_petsc |
( |
| ) |
|
|
inline |
Destroys the PETSc DM and mass matrix objects used for the pseudo-inverse velocity interpolation.
- Returns
- PetscErrorCode
◆ initialize_pseudo_inverse_interpolation()
void DMWrapper::initialize_pseudo_inverse_interpolation |
( |
int |
nthreads, |
|
|
int |
order, |
|
|
const int |
Nspecies, |
|
|
const double |
vPara_min, |
|
|
const double |
vPara_max, |
|
|
const int |
nvp, |
|
|
const double |
vPerp_min, |
|
|
const double |
vPerp_max, |
|
|
const int |
nmu |
|
) |
| |
|
inline |
◆ initialize_pseudo_inverse_interpolation_petsc()
PetscErrorCode DMWrapper::initialize_pseudo_inverse_interpolation_petsc |
( |
int |
nthreads, |
|
|
PetscInt |
order, |
|
|
const PetscInt |
Nspecies, |
|
|
const PetscReal |
vPara_min, |
|
|
const PetscReal |
vPara_max, |
|
|
const PetscInt |
nvp, |
|
|
const PetscReal |
vPerp_min, |
|
|
const PetscReal |
vPerp_max, |
|
|
const PetscInt |
nmu |
|
) |
| |
|
inline |
Creates the PETSc DM and mass matrix objects used for the pseudo-inverse velocity interpolation, based on the 2D velocity grid parameters.
- Parameters
-
[in] | nthreads | is the number of threads which sets the number of PETSc objects being created. |
[in] | order | is the interpolation order (0 = nearest neighbor), (1 = linear), (2 = quadratic), (3 = cubic), .... |
[in] | Nspecies | is the number of species being interpolated, typically 1. |
[in] | vPara_min | is the minimum v|| of the velocity grid. |
[in] | vPara_max | is the maximum v|| of the velocity grid. |
[in] | num_vPara | is the number of grid points in v||. |
[in] | vPerp_min | is the minimum v⊥ of the velocity grid. |
[in] | vPerp_max | is the maximum v⊥ of the velocity grid. |
[in] | num_vPerp | is the number of grid points in v⊥. |
- Returns
- PetscErrorCode
◆ operator=()
◆ dm_len
Number of PETSc mesh objects to store.
◆ is_owner
◆ obj
The documentation for this struct was generated from the following file:
- /p/test_ssd/builds/t3_84szKM/0/xgc/XGC-Devel/XGC_core/cpp/DM_wrapper.hpp