|
XGC1
|
#include "petsc_version_defs.h"#include <petsc/finclude/petsc.h>
Data Types | |
| interface | f90moduleinterfaces::pcsetapplicationcontext |
| interface | f90moduleinterfaces::pcgetapplicationcontext |
| interface | f90moduleinterfaces::matcreateshell |
| interface | f90moduleinterfaces::matshellsetcontext |
| interface | f90moduleinterfaces::matshellgetcontext |
Modules | |
| module | f90moduleinterfaces |
| Explicit interfaces to somve PETSc function used by the FSA solver. | |
Functions/Subroutines | |
| subroutine | clear_monitors () |
| integer function | get_ntor_num_from_ntor_real (ntor_real, nphi, wedge_n) |
| Calculates the numerical toroidal mode number for a given real toroidal mode number. More... | |
| subroutine | set_petsc_rhs_vec_axisymmetric_poisson (xgc_rhs, xgc_rhs_bd, xgc_sol_bd) |
| Places rhs and (optionally) boundary condition data into a PETSc vector. More... | |
| subroutine | solve_axisymmetric_poisson_iter (pot0m, is_full_solve_int) |
| Poisson solver that uses an iterative method to solve for the axisymmetric electrostatic potential: The Poisson equation for the axisymmetric mode is. More... | |
| subroutine | axisym_mat_mult (A, X, Y, ierr) |
| Subroutine used to apply the shell matrix representing the left hand side of the axisymmetric Poisson equation. More... | |
| subroutine | axisym_ad_mat_mult (A, X, Y, ierr) |
| Subroutine used to apply the shell matrix representing the adiabatic term in the axisymmetric Poisson equation. More... | |
| subroutine | axisym_pc_apply (pc, X, Y, ierr) |
| Subroutine used to apply the shell preconditioner used in the iterative solution of the axisymmetric Poisson equation. More... | |
| subroutine | solve_axisym_poisson_step_one (xgc_rhs00, xgc_rhs_bd00, xgc_sol_bd00, pot0m) |
| Simple two-step 2D solver driver for the axisymmetric potential and one-step solver ("FSA-solver") with constraint equation for the flux-surface averaged potential. For the two-step solver, the RHS is the flux-surface averaged charge density, the LHS consists only of the polarization operator. The result of the 2D solve is flux-surface averaged; i.e., the result is \(\langle\phi\rangle\) The Poisson equation for the flux-surface averaged mode is. More... | |
| subroutine | positive_phi00_sol (tmp00_surf, pot0m) |
| subroutine | get_xpt_i (i_x1, i_x2) |
| integer(c_int) function | get_solver00_nbndry () |
| integer(c_int) function | solver00_use_this_rank_int () |
| integer(c_int) function | solvera_use_this_rank_int (CV_int) |
| subroutine | solve_ampere_cv (xgc_rhs, xgc_rhs2, Ah) |
| subroutine | solve_ampere_cv_spec (xgc_rhs_spectral, xgc_rhs2_spectral, xgc_field_spectral) |
| subroutine | solve_ampere (xgc_rhs, xgc_rhs2, Ah) |
| subroutine | solve_ampere_spec (xgc_rhs_spectral, xgc_rhs2_spectral, xgc_field_spectral) |
| subroutine axisym_ad_mat_mult | ( | A, | |
| X, | |||
| Y, | |||
| ierr | |||
| ) |
Subroutine used to apply the shell matrix representing the adiabatic term in the axisymmetric Poisson equation.
\[ \frac{e n_0/T_{e,0}} \left( \overline{\phi} - \alpha(\psi) \langle \overline{\phi} \rangle_{fs} \right) \]
| [in] | A | Shell matrix which multiplies X, Petsc Mat |
| [in] | X | Vector to be multiplied, Petsc Vec |
| [out] | Y | Vector to store product A*X, Petsc Vec |
| [out] | ierr | PETSc Error code |

| subroutine axisym_mat_mult | ( | A, | |
| X, | |||
| Y, | |||
| ierr | |||
| ) |
Subroutine used to apply the shell matrix representing the left hand side of the axisymmetric Poisson equation.
\[ -\nabla_\perp \cdot \xi \nabla_\perp \overline{\phi} + \frac{e n_0/T_{e,0}} \left( \overline{\phi} - \langle \overline{\phi} \rangle_{fs} \right) \]
| [in] | A | Shell matrix which multiplies X, Petsc Mat |
| [in] | X | Vector to be multiplied, Petsc Vec |
| [out] | Y | Vector to store product A*X, Petsc Vec |
| [out] | ierr | PETSc Error code |

| subroutine axisym_pc_apply | ( | pc, | |
| X, | |||
| Y, | |||
| ierr | |||
| ) |
Subroutine used to apply the shell preconditioner used in the iterative solution of the axisymmetric Poisson equation.
| [in] | pc | Shell preconditioner which is applied to residual X, Petsc PC |
| [in] | X | Residual vector, Petsc Vec |
| [out] | Y | Preconditioned residual vector, Petsc Vec |
| [out] | ierr | PETSc Error code |

| subroutine clear_monitors |
| integer function get_ntor_num_from_ntor_real | ( | integer, intent(in) | ntor_real, |
| integer, intent(in) | nphi, | ||
| integer, intent(in) | wedge_n | ||
| ) |
Calculates the numerical toroidal mode number for a given real toroidal mode number.
| [in] | ntor_real | Real toroidal mode number, integer |
| [in] | nphi | Number of toroidal planes per wedge, integer |
| [in] | wedge_n | Number of wedges per toroidal circuit, integer |

| integer(c_int) function get_solver00_nbndry |
| subroutine get_xpt_i | ( | integer(c_int) | i_x1, |
| integer(c_int) | i_x2 | ||
| ) |
| subroutine positive_phi00_sol | ( | real (kind=8), dimension(grid%npsi_surf) | tmp00_surf, |
| real(8), dimension(grid%nnode) | pot0m | ||
| ) |

| subroutine set_petsc_rhs_vec_axisymmetric_poisson | ( | real(8), dimension(grid%nnode) | xgc_rhs, |
| real(8), dimension(psn%solver00%n_boundary) | xgc_rhs_bd, | ||
| real(8), dimension(psn%solver00%n_boundary) | xgc_sol_bd | ||
| ) |
Places rhs and (optionally) boundary condition data into a PETSc vector.
| [in] | grid | Solver grid data, type(grid_type) |
| [in,out] | psn | Field data; the rhs XGC data is stored there, type(psn_type) |
| [in,out] | solver | XGC solver object; the rhs PETSc vector is stored there, type(xgc_solver) |

| subroutine solve_ampere | ( | real (8), dimension(grid%nnode) | xgc_rhs, |
| real (8), dimension(grid%nnode) | xgc_rhs2, | ||
| real (8), dimension(grid%nnode), intent(out) | Ah | ||
| ) |

| subroutine solve_ampere_cv | ( | real (8), dimension(grid%nnode) | xgc_rhs, |
| real (8), dimension(grid%nnode) | xgc_rhs2, | ||
| real (8), dimension(grid%nnode), intent(out) | Ah | ||
| ) |

| subroutine solve_ampere_cv_spec | ( | real (8), dimension(grid%nnode,2) | xgc_rhs_spectral, |
| real (8), dimension(grid%nnode,2) | xgc_rhs2_spectral, | ||
| real (8), dimension(grid%nnode,2) | xgc_field_spectral | ||
| ) |

| subroutine solve_ampere_spec | ( | real (8), dimension(grid%nnode,2) | xgc_rhs_spectral, |
| real (8), dimension(grid%nnode,2) | xgc_rhs2_spectral, | ||
| real (8), dimension(grid%nnode,2) | xgc_field_spectral | ||
| ) |

| subroutine solve_axisym_poisson_step_one | ( | real (8), dimension(grid%nnode) | xgc_rhs00, |
| real (8), dimension(psn%solver00%n_boundary) | xgc_rhs_bd00, | ||
| real (8), dimension(psn%solver00%n_boundary) | xgc_sol_bd00, | ||
| real(8), dimension(grid%nnode) | pot0m | ||
| ) |
Simple two-step 2D solver driver for the axisymmetric potential and one-step solver ("FSA-solver") with constraint equation for the flux-surface averaged potential. For the two-step solver, the RHS is the flux-surface averaged charge density, the LHS consists only of the polarization operator. The result of the 2D solve is flux-surface averaged; i.e., the result is \(\langle\phi\rangle\) The Poisson equation for the flux-surface averaged mode is.
\[ -\nabla_\perp \cdot \xi \nabla_\perp \overline{\phi} = e \left\langle \left( \langle \overline{\delta n_i} \rangle_g - \overline{\delta n_{e,NA}} \right) \right\rangle_{fs} \]
where \(xi\) is the electric susceptibility, \(e\) is the elementary charge, \(n_0\) is the background density, \(\langle \overline{\delta n_i} \rangle_g\) is the gyroaveraged ion (gyrocenter) charge density, \(\overline{\delta n_{e,NA}}\) is the non-adiabatic electron charge density, and \(\langle\dots\rangle_{fs}\) is the the flux-surface average. \(\overline{\dots}\) is the toroidal average.
NOTE: The two-step solver does not support non-zero Dirichlet boundary conditions!
The FSA solver solves the same equation as the outer-iterative solver (solve_axisymmetric_poisson_iter), but with a twist. The flux-surface averaged potential is treated as an independent field \(\lambda\) in the Poisson equation. In order to make the equation consistent, a second (constraint) equation is added that enforces \(\lambda=\langle\phi\rangle_{fs}\). So the FSA solver solves a 2x2 block system.
| [in] | grid | Solver grid data, type(grid_type) |
| [in,out] | psn | Field data; the potential is stored there, type(psn_type) |

| subroutine solve_axisymmetric_poisson_iter | ( | real(8), dimension(grid%nnode) | pot0m, |
| integer, intent(in), value | is_full_solve_int | ||
| ) |
Poisson solver that uses an iterative method to solve for the axisymmetric electrostatic potential: The Poisson equation for the axisymmetric mode is.
\[ -\nabla_\perp \cdot \xi \nabla_\perp \overline{\phi}_{k+1} + \frac{e n_0/T_{e,0}} \overline{\phi}_{k+1} = e \left( \langle \overline{\delta n_i} \rangle_g - \overline{\delta n_{e,NA}} \right) + \frac{e n_0/T_{e,0}} \langle \phi_{k} \rangle_{fs} \]
where \(k\) is the iteration index, \(xi\) is the electric susceptibility, \(e\) is the elementary charge, \(n_0\) and \(T_0\) are the background density and temperature, \(\langle \overline{\delta n_i} \rangle_g\) is the gyroaveraged ion (gyrocenter) charge density, \(\overline{\delta n_{e,NA}}\) is the non-adiabatic electron charge density, and \(\langle\dots\rangle_{fs}\) is the the flux-surface average. \(\overline{\dots}\) is the toroidal average.
The RHS and LHS of the axisymmetric equation can be Fourier-filtered (poloidal modes). The LHS of the non-axisymmetric equation can be Fourier-filtered in the toroidal and poloidal direction. The solver in both cases is a linear finite-element solver on an unstructured triangular grid.
| [in] | grid | Solver grid data, type(grid_type) |
| [in,out] | psn | Field data; the potential is stored there, type(psn_type) |

| integer(c_int) function solver00_use_this_rank_int |
| integer(c_int) function solvera_use_this_rank_int | ( | integer(c_int), intent(in), value | CV_int | ) |