XGC1
Public Member Functions | Public Attributes | List of all members
DMWrapper Struct Reference

#include <DM_wrapper.hpp>

Public Member Functions

 DMWrapper ()
 
 ~DMWrapper ()
 
 DMWrapper (const DMWrapper &original)
 
DMWrapperoperator= (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 ()
 

Public Attributes

bool is_owner
 
int dm_len
 Number of PETSc mesh objects to store. More...
 
Kokkos::View< DMObjects *, Kokkos::LayoutRight, HostTypeobj
 Stores number of threads DMObjects. More...
 

Constructor & Destructor Documentation

◆ DMWrapper() [1/2]

DMWrapper::DMWrapper ( )
inline

◆ ~DMWrapper()

DMWrapper::~DMWrapper ( )
inline

◆ DMWrapper() [2/2]

DMWrapper::DMWrapper ( const DMWrapper original)
inline

Member Function Documentation

◆ destroy_DM()

void DMWrapper::destroy_DM ( )
inline
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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
Here is the caller graph for this function:

◆ 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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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]nthreadsis the number of threads which sets the number of PETSc objects being created.
[in]orderis the interpolation order (0 = nearest neighbor), (1 = linear), (2 = quadratic), (3 = cubic), ....
[in]Nspeciesis the number of species being interpolated, typically 1.
[in]vPara_minis the minimum v|| of the velocity grid.
[in]vPara_maxis the maximum v|| of the velocity grid.
[in]num_vParais the number of grid points in v||.
[in]vPerp_minis the minimum v⊥ of the velocity grid.
[in]vPerp_maxis the maximum v⊥ of the velocity grid.
[in]num_vPerpis the number of grid points in v⊥.
Returns
PetscErrorCode
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

DMWrapper& DMWrapper::operator= ( const DMWrapper original)
inline

Member Data Documentation

◆ dm_len

int DMWrapper::dm_len

Number of PETSc mesh objects to store.

◆ is_owner

bool DMWrapper::is_owner

◆ obj

Kokkos::View<DMObjects*, Kokkos::LayoutRight, HostType> DMWrapper::obj

Stores number of threads DMObjects.


The documentation for this struct was generated from the following file: