XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Functions
velocity_interpolate.hpp File Reference
#include <petscversion.h>
#include <petscsystypes.h>
#include <petscsys.h>
#include <petscdmplex.h>
#include <petscds.h>
#include "petscdmswarm.h"
#include "petscksp.h"
Include dependency graph for velocity_interpolate.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

void initialize_pseudo_inverse_interpolation (DM &dm, PetscInt order, const PetscInt Nspecies, const PetscReal vPara_min, const PetscReal vPara_max, const PetscInt num_vPara, const PetscReal vPerp_min, const PetscReal vPerp_max, const PetscInt num_vPerp)
 
void pseudo_inverse_interpolate_particles_to_velocity_grid (DM &dm, PetscReal *particle_coords_fptr, PetscReal *particle_weights_fptr, const PetscInt Nparticles, PetscInt order, const PetscInt Nspecies, PetscReal *grid_coords_fptr, PetscReal *grid_weights_fptr, PetscInt *Ngrid)
 
void pseudo_inverse_interpolate_velocity_grid_to_particles (DM &dm, PetscReal *grid_coords_fptr, PetscReal *grid_weights_fptr, const PetscInt Ngrid, PetscInt order, const PetscInt Nspecies, PetscReal *particle_coords_fptr, PetscReal *particle_weights_fptr, const PetscInt Nparticles, bool *ksp_converged)
 

Detailed Description

1.0 A. Mollen, M. Adams, M. Knepley Feb 2021

DESCRIPTION

Header file to velocity_interpolate.cpp

Function Documentation

void initialize_pseudo_inverse_interpolation ( DM &  dm,
PetscInt  order,
const PetscInt  Nspecies,
const PetscReal  vPara_min,
const PetscReal  vPara_max,
const PetscInt  num_vPara,
const PetscReal  vPerp_min,
const PetscReal  vPerp_max,
const PetscInt  num_vPerp 
)

Creates the PETSc DM object used for the pseudo-inverse velocity interpolation, based on the 2D velocity grid parameters.

Parameters
[in,out]dmis the PETSc DM object created.
[in]orderis the interpolation order (1 = linear), (2 = quadratic).
[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
void

Here is the call graph for this function:

void pseudo_inverse_interpolate_particles_to_velocity_grid ( DM &  dm,
PetscReal *  particle_coords_fptr,
PetscReal *  particle_weights_fptr,
const PetscInt  Nparticles,
PetscInt  order,
const PetscInt  Nspecies,
PetscReal *  grid_coords_fptr,
PetscReal *  grid_weights_fptr,
PetscInt *  Ngrid 
)

Performs the particles -> velocity grid interpolation of marker particle weights to 2D velocity grid weights.

Parameters
[in]dmis the PETSc DM object created by cpp_initialize_pseudo_inverse_interpolation.
[in]particle_coords_fptris a 2*Nparticles-sized array containing the particle coordinates.
[in]particle_weights_fptris a Nparticles-sized array containing the particle weights.
[in]Nparticlesis the number of particles being interpolated.
[in]orderis the interpolation order (1 = linear), (2 = quadratic).
[in]Nspeciesis the number of species being interpolated, typically 1.
[in,out]grid_coords_fptris an array containing the velocity grid coordinates.
[in,out]grid_weights_fptris an array containing the velocity grid weights.
[in,out]Ngridis the number of grid points.
Returns
void

Here is the call graph for this function:

void pseudo_inverse_interpolate_velocity_grid_to_particles ( DM &  dm,
PetscReal *  grid_coords_fptr,
PetscReal *  grid_weights_fptr,
const PetscInt  Ngrid,
PetscInt  order,
const PetscInt  Nspecies,
PetscReal *  particle_coords_fptr,
PetscReal *  particle_weights_fptr,
const PetscInt  Nparticles,
bool *  ksp_converged 
)

Performs the velocity grid -> particles interpolation of 2D velocity grid weights to marker particle weights.

Parameters
[in]dmis the PETSc DM object created by cpp_initialize_pseudo_inverse_interpolation.
[in]grid_coords_fptris an array containing the velocity grid coordinates.
[in]grid_weights_fptris an array containing the velocity grid weights.
[in]Ngridis the number of grid points.
[in]orderis the interpolation order (1 = linear), (2 = quadratic).
[in]Nspeciesis the number of species being interpolated, typically 1.
[in,out]particle_coords_fptris a 2*Nparticles-sized array containing the particle coordinates.
[in,out]particle_weights_fptris a Nparticles-sized array containing the particle weights.
[in]Nparticlesis the number of particles being interpolated.
[out]ksp_convergedis true if the pseudo-inverse calculation converged, otherwise false.
Returns
void

Here is the call graph for this function: