XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
current_drive.cpp File Reference
#include "current_drive.hpp"
#include "electric_field.hpp"
#include <mpi.h>
Include dependency graph for current_drive.cpp:

Functions

double * get_jpar0_loc ()
 
KOKKOS_INLINE_FUNCTION double parallel_spitzer_resistivity (double den, double te)
 
void f_current_drive (const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, Plasma &plasma, const Moments &f0_moments, const DomainDecomposition< DeviceType > &pol_decomp, ElectricField< DeviceType > &electric_field, const MagneticField< DeviceType > &magnetic_field, LoopVolDiagnostics &loop_vol_diag)
 

Function Documentation

void f_current_drive ( const Simulation< DeviceType > &  sml,
const Grid< DeviceType > &  grid,
Plasma plasma,
const Moments f0_moments,
const DomainDecomposition< DeviceType > &  pol_decomp,
ElectricField< DeviceType > &  electric_field,
const MagneticField< DeviceType > &  magnetic_field,
LoopVolDiagnostics loop_vol_diag 
)

Drives electron current such that the plasma current required by the Grad-Shafranov equilibrium magnetic field is maintained. The plasma current is \(j_{\parallel,0} = \hat{\boldsymbol{b}}\cdot\left(\nabla\times\boldsymbol{B}_0\right)/\mu_0 \), or \(j_{\parallel,0} = \sum_s (q_s n_s u_{\parallel,s}) \). Therefore, \(u_{\parallel,e}=j_{\parallel,0}-\sum_{ions}(q_s n_s u_{\parallel,s}) \). The actual electron flow is stored in f0_moments.u_local). Any difference between the targeted and actual flow is damped away with an exponential damping \(\partial u/\partial t = \gamma (u_{target}-u_{actual}) \).

Parameters
[in]smlSimulation object (runtime parameters), class Simulation
[in]gridGrid object, class Grid
[in,out]plasmaXGC plasma object; class plasma
[in]f0_momentsCurrent low-order plasma moments, struct Moments
[in]pol_decompDomain decomposition object (patch boundaries, etc.), class DomainDecomposition
[in]loop_vol_diagDiagnostic ouput object for the loop voltage, class LoopVolDiagnostics

Here is the call graph for this function:

Here is the caller graph for this function:

double* get_jpar0_loc ( )

Here is the caller graph for this function:

KOKKOS_INLINE_FUNCTION double parallel_spitzer_resistivity ( double  den,
double  te 
)

Evaluates the parallel Spitzer resistivity in \(\Omega m\) for given density and temperatures. Uses a simplified formula (NRL Formulary 2019, assuming \(T_e>10\) eV) because the logic to adjust V_loop is heuristic anyway.

Parameters
[in]denPlasma (electron) density in \(m^{-3}\), double
[in]teElectron temperature in eV, double

Here is the caller graph for this function: