V3FIT
Functions/Subroutines | Variables
pchelms Module Reference

This file solves the helmholtz equation to set inital fields that match vmec and vacuum currents from the vector potential. Initial vector potential on the edge is supplied by the BMW code. More...

Functions/Subroutines

subroutine run_pchelms
 Run the pchelms code to solve for inital jbsup values.
 
subroutine curla_pchelms (Asubsmnf, Asubumnf, Asubvmnf,
 Compute real-space components of contravariant jac*B^i on half radial grid from curl(A). More...
 
subroutine curlb_pchelms (ksupsmnf, ksupumnf, ksupvmnf, asubsmnf,
 Compute real-space components of contravariant jac*J^i on full radial grid from curl(B). More...
 
subroutine cyl2vmec_a (A_r, A_p, A_z, cA_s, cA_u, cA_v)
 Convert cylindical vector potential to contravariant components. More...
 
subroutine init_a (cA_s, cA_u, cA_v, A_s, A_u, A_v, parity)
 Initialize vector potential. More...
 
subroutine init_f (Fsupsmn, Fsupumn, Fsupvmn, jcurrumn, jcurrvmn)
 Initialize expected currents. More...
 
subroutine compare_current (Fsupumn, Fsupvmn, jcurrumn, jcurrvmn)
 Compare Curl(Curl(A)) with the expected VMEC and vacuum currents. More...
 
subroutine dump_a (js, iunit)
 Write out vector potential to file. More...
 
subroutine boundaryconditions (Fsupsmn, Fsupumn, Fsupvmn, iparity)
 Apply boundary conditions. More...
 
subroutine gmres_pchelms
 Setup and run GMRES solver for the helmholtz problem. More...
 
subroutine matvec_pchelms (p, Ap, ndim)
 Serial callback function for MatVec GMRES routine. More...
 
subroutine matvec_par_pchelms (ploc, Ap, nloc)
 Parallel callback function for MatVec GMRES routine. More...
 
subroutine getnlforce_pchelms (xcstate, fsq_nl, bnorm)
 Non linear force callback. More...
 

Variables

logical, parameter, private l_asedge = .TRUE.
 Controls if the edge Aubs values are evolved or fixed. More...
 
integer, private nsmin
 
integer, private nsmax
 
integer, private ns_match
 
logical, private linit
 
logical, private lhessian
 
real(dp), private bnorm = -1
 
real(dp), private line_bu
 
real(dp), dimension(:,:,:), allocatable, private bsupsijh
 
real(dp), dimension(:,:,:), allocatable, private bsupuijh
 
real(dp), dimension(:,:,:), allocatable, private bsupvijh
 
real(dp), dimension(:,:,:), allocatable, private bsubsijh
 
real(dp), dimension(:,:,:), allocatable, private bsubuijh
 
real(dp), dimension(:,:,:), allocatable, private bsubvijh
 
real(dp), dimension(:,:,:), allocatable, private jacobmnch
 
real(dp), dimension(:,:,:), allocatable, private jacobmnsh
 

Detailed Description

This file solves the helmholtz equation to set inital fields that match vmec and vacuum currents from the vector potential. Initial vector potential on the edge is supplied by the BMW code.

Function/Subroutine Documentation

◆ boundaryconditions()

subroutine pchelms::boundaryconditions (   Fsupsmn,
  Fsupumn,
  Fsupvmn,
  iparity 
)

Apply boundary conditions.

Parameters
(inout)Fsupsmn Contravariant current in the s direction.
(inout)Fsupumn Contravariant current in the u direction.
(inout)Fsupvmn Contravariant current in the v direction.
(in)iparity Fourier parity of the residuals.

Definition at line 1353 of file pchelms.f90.

◆ compare_current()

subroutine pchelms::compare_current ( real (dp), dimension(:,:,:), intent(in)  Fsupumn,
real (dp), dimension(:,:,:), intent(in)  Fsupvmn,
real (dp), dimension(:,:,:), intent(in)  jcurrumn,
real (dp), dimension(:,:,:), intent(in)  jcurrvmn 
)

Compare Curl(Curl(A)) with the expected VMEC and vacuum currents.

Parameters
[out]FsupsmnContravariant current in the s direction on the full grid.
[out]FsupumnContravariant current in the u direction on the full grid.
[out]FsupvmnContravariant current in the v direction on the full grid.
[in]jcurrumnContravariant VMEC current in the u direction on the full grid.
[in]jcurrvmnContravariant VMEC current in the v direction on the full grid.
[in]iparityParity flag to indicate stellaratory symmetry.

Definition at line 650 of file pchelms.f90.

◆ curla_pchelms()

subroutine pchelms::curla_pchelms (   Asubsmnf,
  Asubumnf,
  Asubvmnf 
)

Compute real-space components of contravariant jac*B^i on half radial grid from curl(A).

Parameters
[in]AsubsmnfCovariant vector potential in the s direction.
[in]AsubumnfCovariant vector potential in the u direction.
[in]AsubvmnfCovariant vector potential in the v direction.
[in,out]jbsupsmnhContravariant B field and jacobian in the s direction.
[in,out]jbsupumnhContravariant B field and jacobian in the u direction.
[in,out]jbsupvmnhContravariant B field and jacobian in the v direction.
[in]iparityParity flag to indicate stellaratory symmetry.

Definition at line 187 of file pchelms.f90.

◆ curlb_pchelms()

subroutine pchelms::curlb_pchelms ( real(dp), dimension(:,:,:), pointer  ksupsmnf,
real(dp), dimension(:,:,:), pointer  ksupumnf,
real(dp), dimension(:,:,:), pointer  ksupvmnf,
real(dp), dimension(:,:,:), pointer  asubsmnf 
)

Compute real-space components of contravariant jac*J^i on full radial grid from curl(B).

Parameters
[in,out]ksupsmnfContravariant current and jacobian in the s direction on the full grid.
[in,out]ksupumnfContravariant current and jacobian in the u direction on the full grid.
[in,out]ksupvmnfContravariant current and jacobian in the v direction on the full grid.
[in,out]asubsmnfCovariant vector potential in the s direction on the full grid.
[in]bsupsmnhContravariant current and jacobian in the v direction on the half grid.
[in]iparityParity flag to indicate stellaratory symmetry.
[in,out]curtorCurrent enclosed at the boundary.

Definition at line 285 of file pchelms.f90.

◆ cyl2vmec_a()

subroutine pchelms::cyl2vmec_a ( real (dp), dimension(:,:), intent(in)  A_r,
real (dp), dimension(:,:), intent(in)  A_p,
real (dp), dimension(:,:), intent(in)  A_z,
real (dp), dimension(:,:), intent(out)  cA_s,
real (dp), dimension(:,:), intent(out)  cA_u,
real (dp), dimension(:,:), intent(out)  cA_v 
)

Convert cylindical vector potential to contravariant components.

Parameters
[in]A_rCylindical r component of the vector potential at the last surface.
[in]A_pCylindical p component of the vector potential at the last surface.
[in]A_zCylindical z component of the vector potential at the last surface.
[out]cA_sContravariant vector potential in the s direction on the last surface grid.
[out]cA_uContravariant vector potential in the u direction on the last surface grid.
[out]cA_vContravariant vector potential in the v direction on the last surface grid.

Definition at line 441 of file pchelms.f90.

◆ dump_a()

subroutine pchelms::dump_a ( integer, intent(in)  js,
integer, intent(in)  iunit 
)

Write out vector potential to file.

Parameters
[in]jsRadial index.
[in]iunitFile unit to write to.

Definition at line 707 of file pchelms.f90.

◆ getnlforce_pchelms()

subroutine pchelms::getnlforce_pchelms ( real (dp), dimension(neqs), intent(in)  xcstate,
real (dp), intent(out)  fsq_nl,
real (dp), intent(in)  bnorm 
)

Non linear force callback.

Parameters
[in]xcstateDisplacement parameters.
[out]fsq_nl|F^2| non-linear residual.
[in]bnormInternal GMRES normalization

Definition at line 1578 of file pchelms.f90.

◆ gmres_pchelms()

subroutine pchelms::gmres_pchelms

Setup and run GMRES solver for the helmholtz problem.

GMRES solves r == b + Ax = 0

Definition at line 1396 of file pchelms.f90.

◆ init_a()

subroutine pchelms::init_a ( real (dp), dimension(:,:,:), intent(in)  cA_s,
real (dp), dimension(:,:,:), intent(in)  cA_u,
real (dp), dimension(:,:,:), intent(in)  cA_v,
real (dp), dimension(0:mpol,-ntor:ntor,ns), intent(out)  A_s,
real (dp), dimension(0:mpol,-ntor:ntor,ns), intent(out)  A_u,
real (dp), dimension(0:mpol,-ntor:ntor,ns), intent(out)  A_v,
integer, intent(in)  parity 
)

Initialize vector potential.

This subroutine is only called once so no need to parallelize. Linearly extrapolat to the center or match surface.

Parameters
[in,out]cA_sContravariant vector potential in the s direction on the full grid.
[in,out]cA_uContravariant vector potential in the u direction on the full grid.
[in,out]cA_vContravariant vector potential in the v direction on the full grid.
[out]A_sContravariant fourier vector potential in the s direction on the full grid.
[out]A_uContravariant fourier vector potential in the u direction on the full grid.
[out]A_vContravariant fourier vector potential in the v direction on the full grid.
[in]parityParity flag to indicate stellaratory symmetry.

Definition at line 495 of file pchelms.f90.

◆ init_f()

subroutine pchelms::init_f ( real (dp), dimension(0:mpol,-ntor:ntor,ns), intent(out)  Fsupsmn,
real (dp), dimension(0:mpol,-ntor:ntor,ns), intent(out)  Fsupumn,
real (dp), dimension(0:mpol,-ntor:ntor,ns), intent(out)  Fsupvmn,
  jcurrumn,
  jcurrvmn 
)

Initialize expected currents.

Up to the vmec last closed flux surface, initialize with the VMEC currents. outside, initalize to zero.

Parameters
[out]FsupsmnContravariant current in the s direction on the full grid.
[out]FsupumnContravariant current in the u direction on the full grid.
[out]FsupvmnContravariant current in the v direction on the full grid.
[in]jcurrumnContravariant VMEC current in the u direction on the full grid.
[in]jcurrvmnContravariant VMEC current in the v direction on the full grid.

Definition at line 586 of file pchelms.f90.

◆ matvec_par_pchelms()

subroutine pchelms::matvec_par_pchelms ( real(dp), dimension(nloc), intent(in)  ploc,
real(dp), dimension(nloc), intent(out)  Ap,
integer, intent(in)  nloc 
)

Parallel callback function for MatVec GMRES routine.

Parameters
[in]plocLocal displacement parameters.
[out]ApLocal Matrix element.
[in]nlocLocal number of dimensiona.

Definition at line 1532 of file pchelms.f90.

◆ matvec_pchelms()

subroutine pchelms::matvec_pchelms ( real(dp), dimension(ndim), intent(in)  p,
real(dp), dimension(ndim), intent(out)  Ap,
integer, intent(in)  ndim 
)

Serial callback function for MatVec GMRES routine.

Parameters
[in]pDisplacement parameters.
[out]ApMatrix element.
[in]ndimNumber of dimensiona.

Definition at line 1506 of file pchelms.f90.

Variable Documentation

◆ l_asedge

logical, parameter, private pchelms::l_asedge = .TRUE.
private

Controls if the edge Aubs values are evolved or fixed.

  • True Evolve the edge Asubs values
  • False Fix edge Asubs values.

Definition at line 46 of file pchelms.f90.