PB3D [2.47]
Ideal linear high-n MHD stability in 3-D
Loading...
Searching...
No Matches
sol_utilities::calc_xuq Interface Reference

Calculates the normal ( \(\cdot_n\)) or geodesic ( \(\cdot_g\)) components of the plasma perturbation \(\vec{\xi}\) or the magnetic perturbation \(\vec{Q} = \nabla \times \left(\vec{x} \times \vec{B}\right)\). More...

Public Member Functions

integer function calc_xuq_ind (mds_x, mds_sol, grid_eq, grid_x, grid_sol, eq_1, eq_2, x, sol, x_id, xuq_style, time, xuq, deriv)
 (time) individual version
integer function calc_xuq_arr (mds_x, mds_sol, grid_eq, grid_x, grid_sol, eq_1, eq_2, x, sol, x_id, xuq_style, time, xuq, deriv)
 (time) array version

Detailed Description

Calculates the normal ( \(\cdot_n\)) or geodesic ( \(\cdot_g\)) components of the plasma perturbation \(\vec{\xi}\) or the magnetic perturbation \(\vec{Q} = \nabla \times \left(\vec{x} \times \vec{B}\right)\).

Which variables are plot depends on XUQ_style [17] :

  • XUQ_style = 1: \(X\) (supports parallel derivative)
  • XUQ_style = 2: \(U\) (supports parallel derivative)
  • XUQ_style = 3: \(Q_n\)
  • XUQ_style = 4: \(Q_g\)

\(X\) and \(U\) support the calculation of the parallel derivative.

For \(Q_n\) and \(Q_g\), the metric variables have to be provided as well.

The perturbation grid is assumed to have the same angular coordinates as the equilibrium grid, and the normal coordinates correspond to either the equilibrium grid (X_grid_style 1), the solution grid (X_grid_style 2) or the enriched equilibrium grid (X_grid_style 3).

The output grid, furthermore, has the angular part corresponding to the equilibrium grid, and the normal part to the solution grid.

Returns
ierr

Definition at line 56 of file sol_utilities.f90.

Member Function/Subroutine Documentation

◆ calc_xuq_arr()

integer function sol_utilities::calc_xuq::calc_xuq_arr ( type(modes_type), intent(in), target mds_x,
type(modes_type), intent(in), target mds_sol,
type(grid_type), intent(in) grid_eq,
type(grid_type), intent(in) grid_x,
type(grid_type), intent(in) grid_sol,
type(eq_1_type), intent(in) eq_1,
type(eq_2_type), intent(in) eq_2,
type(x_1_type), intent(in), target x,
type(sol_type), intent(in) sol,
integer, intent(in) x_id,
integer, intent(in) xuq_style,
real(dp), dimension(:), intent(in) time,
complex(dp), dimension(:,:,:,:), intent(inout) xuq,
logical, intent(in), optional deriv )

(time) array version

Parameters
[in]mds_xgeneral modes variables in perturbation grid
[in]mds_solgeneral modes variables in solution grid
[in]grid_eqequilibrium grid
[in]grid_xperturbation grid
[in]grid_solsolution grid
[in]eq_1flux equilibrium
[in]eq_2metric equilibrium
[in]xperturbation variables
[in]solsolution variables
[in]x_idnr. of Eigenvalue
[in]xuq_stylewhether to calculate \(X\), \(U\), \(Q_n\) or \(Q_g\)
[in]timetime range
[in,out]xuq\(X\), \(U\), \(Q_n\) or \(Q_g\)
[in]derivreturn parallel derivative

Definition at line 65 of file sol_utilities.f90.

Here is the call graph for this function:

◆ calc_xuq_ind()

integer function sol_utilities::calc_xuq::calc_xuq_ind ( type(modes_type), intent(in) mds_x,
type(modes_type), intent(in) mds_sol,
type(grid_type), intent(in) grid_eq,
type(grid_type), intent(in) grid_x,
type(grid_type), intent(in) grid_sol,
type(eq_1_type), intent(in) eq_1,
type(eq_2_type), intent(in) eq_2,
type(x_1_type), intent(in) x,
type(sol_type), intent(in) sol,
integer, intent(in) x_id,
integer, intent(in) xuq_style,
real(dp), intent(in) time,
complex(dp), dimension(:,:,:), intent(inout) xuq,
logical, intent(in), optional deriv )

(time) individual version

Parameters
[in]mds_xgeneral modes variables in perturbation grid
[in]mds_solgeneral modes variables in solution grid
[in]grid_eqequilibrium grid
[in]grid_xperturbation grid
[in]grid_solsolution grid
[in]eq_1flux equilibrium
[in]eq_2metric equilibrium
[in]xperturbation variables
[in]solsolution variables
[in]x_idnr. of Eigenvalue
[in]xuq_stylewhether to calculate X, U, Qn or Qg
[in,out]xuqnormal component of perturbation
[in]derivreturn parallel derivative

Definition at line 518 of file sol_utilities.f90.


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