V3FIT
Functions/Subroutines | Variables
metrics Module Reference

Module for reading in VMEC input and computing metric elements on the SIESTA (sqrt-flux) radial grid. Computes and stores the real-space metric elements, Jacobian based on square root of fluc splined vmec coordinate system. More...

Functions/Subroutines

subroutine init_metric_elements ()
 Initialize the metric elements. More...
 
subroutine set_grid_sizes (mpol_in, ntor_in)
 Set grid sizes. More...
 
subroutine surfavg (average, q3d, nsmin, nsmax)
 Surface average a quantity. More...
 
subroutine loadgrid (istat)
 Load the R, Z and lambda arrays from VMEC. More...
 
subroutine rz_to_ijsp (rmn, zmn, lasym)
 Transform from fourier R Z to real space qauntities. More...
 
subroutine half_mesh_metrics (r1_i, ru_i, rv_i, z1_i, zu_i, zv_i)
 Compute the metric elements on the half mesh. More...
 
subroutine full_mesh_metrics
 Gets the full grid metric elements. More...
 
subroutine toupper (xsubsij, xsubuij, xsubvij,
 Converts to contravariant component full grid. More...
 
subroutine tolowerh (xsupsij, xsupuij, xsupvij,
 Converts to covariant component half grid. More...
 
subroutine tolowerf (xsupsij, xsupuij, xsupvij,
 Converts to covariant component full grid. More...
 
subroutine cleanup_metric_elements
 Deallocate memory containing metric elements on the half mesh. More...
 
subroutine dealloc_full_lower_metrics
 Deallocate memory containing metric elements on the full mesh.
 

Variables

real(dp), dimension(:,:), allocatable sqrtg
 Jacobian on half grid.
 
real(dp), dimension(:,:), allocatable gss
 Lower metric tensor half grid. e_s . e_s.
 
real(dp), dimension(:,:), allocatable gsu
 Lower metric tensor half grid. e_s . e_u.
 
real(dp), dimension(:,:), allocatable gsv
 Lower metric tensor half grid. e_s . e_v.
 
real(dp), dimension(:,:), allocatable guu
 Lower metric tensor half grid. e_u . e_u.
 
real(dp), dimension(:,:), allocatable guv
 Lower metric tensor half grid. e_u . e_v.
 
real(dp), dimension(:,:), allocatable gvv
 Lower metric tensor half grid. e_v . e_v.
 
real(dp), dimension(:,:), allocatable hss
 Upper metric tensor half grid. e^s . e^s.
 
real(dp), dimension(:,:), allocatable hsu
 Upper metric tensor. e^s . e^u.
 
real(dp), dimension(:,:), allocatable hsv
 Upper metric tensor. e^s . e^v.
 
real(dp), dimension(:,:), allocatable huu
 Upper metric tensor. e^u . e^u.
 
real(dp), dimension(:,:), allocatable huv
 Upper metric tensor. e^u . e^v.
 
real(dp), dimension(:,:), allocatable hvv
 Upper metric tensor. e^v . e^v.
 
real(dp), dimension(:,:), allocatable gssf
 Lower metric tensor full grid. e_s . e_s.
 
real(dp), dimension(:,:), allocatable gsuf
 Lower metric tensor full grid. e_s . e_u.
 
real(dp), dimension(:,:), allocatable gsvf
 Lower metric tensor full grid. e_s . e_v.
 
real(dp), dimension(:,:), allocatable guuf
 Lower metric tensor full grid. e_u . e_u.
 
real(dp), dimension(:,:), allocatable guvf
 Lower metric tensor full grid. e_u . e_v.
 
real(dp), dimension(:,:), allocatable gvvf
 Lower metric tensor full grid. e_v . e_v.
 
real(dp) rmax
 Maximum of the grid R inside the last closed flux surface.
 
real(dp) rmin
 Minimum of the grid R inside the last closed flux surface.
 
real(dp) zmax
 Maximum of the grid Z inside the last closed flux surface.
 
real(dp) zmin
 Minimum of the grid Z inside the last closed flux surface.
 
real(dp) skston
 Start timer.
 
real(dp) skstoff
 Stop timer.
 

Detailed Description

Module for reading in VMEC input and computing metric elements on the SIESTA (sqrt-flux) radial grid. Computes and stores the real-space metric elements, Jacobian based on square root of fluc splined vmec coordinate system.

Function/Subroutine Documentation

◆ cleanup_metric_elements()

subroutine metrics::cleanup_metric_elements

Deallocate memory containing metric elements on the half mesh.

Also deallocates the grid variables.

Note
sqrtg is deallocated in init_bcovar and stored in jacob variable

Definition at line 642 of file metrics.f90.

◆ full_mesh_metrics()

subroutine metrics::full_mesh_metrics

Gets the full grid metric elements.

This subroutine gets the lower metric elements on the full mesh. Preserves positive definiteness of metric tensor.

Definition at line 401 of file metrics.f90.

◆ half_mesh_metrics()

subroutine metrics::half_mesh_metrics ( real(dp), dimension(nuv,ns_i), intent(in)  r1_i,
real(dp), dimension(nuv,ns_i), intent(in)  ru_i,
real(dp), dimension(nuv,ns_i), intent(in)  rv_i,
real(dp), dimension(nuv,ns_i), intent(in)  z1_i,
real(dp), dimension(nuv,ns_i), intent(in)  zu_i,
real(dp), dimension(nuv,ns_i), intent(in)  zv_i 
)

Compute the metric elements on the half mesh.

Parameters
[in]r1_iR on the full mesh.
[in]ru_idRdu on the full mesh.
[in]rv_idRdv on the full mesh.
[in]z1_iZ on the full mesh.
[in]zu_idZdu on the full mesh.
[in]zv_idZdv on the full mesh.

Definition at line 308 of file metrics.f90.

◆ init_metric_elements()

subroutine metrics::init_metric_elements

Initialize the metric elements.

Loads values on to the SIESTA meshes. This involves splining and interpolating the VMEC quantities from the VMEC radial grid to the SIESTA radial grid.

Definition at line 103 of file metrics.f90.

◆ loadgrid()

subroutine metrics::loadgrid ( integer, intent(out)  istat)

Load the R, Z and lambda arrays from VMEC.

R, Z and Lambda are recast onto the SIESTA grid.

Parameters
[out]istatError status.

Definition at line 218 of file metrics.f90.

◆ rz_to_ijsp()

subroutine metrics::rz_to_ijsp ( real(dp), dimension(:,:,:), intent(in)  rmn,
real(dp), dimension(:,:,:), intent(in)  zmn,
logical, intent(in)  lasym 
)

Transform from fourier R Z to real space qauntities.

Computes R, Z, dr/du, dr/dv, dz/du and dz/dv from the fourier representation.

Parameters
[in]rmnRadial fourier amplitudes.
[in]zmnVertical fourier amplitudes.
[in]lasymA symmetric flag.

Definition at line 257 of file metrics.f90.

◆ set_grid_sizes()

subroutine metrics::set_grid_sizes ( integer, intent(in)  mpol_in,
integer, intent(in)  ntor_in 
)

Set grid sizes.

The real space grid is determined from the number of toroidal and poloidal modes.

Parameters
[in]mpol_inNumber of SIESTA poloidal modes.
[in]ntor_inNumber of SIESTA toroidal modes.

Definition at line 127 of file metrics.f90.

◆ surfavg()

subroutine metrics::surfavg ( real(dp), dimension(nsmin:nsmax), intent(out)  average,
real(dp), dimension(nu_i,nv_i,nsmin:nsmax), intent(in)  q3d,
integer, intent(in)  nsmin,
integer, intent(in)  nsmax 
)

Surface average a quantity.

Parameters
[out]averageSurface average
[in]q3d3D quantity is real space.
[in]nsminMinimum radial index.
[in]nsmaxMaximum radial index.

Definition at line 183 of file metrics.f90.

◆ tolowerf()

subroutine metrics::tolowerf ( real(dp), dimension(nuv,ns), intent(in)  xsupsij,
real(dp), dimension(nuv,ns), intent(in)  xsupuij,
real(dp), dimension(nuv,ns), intent(in)  xsupvij 
)

Converts to covariant component full grid.

Computes covariant (lower) from contravariant (upper) magnetic field components on the full mesh.

xsubsij = gssf*xsupsij + gsuf*xsupuij + gsvf*xsupvij xsubuij = gsuf*xsupsij + guuf*xsupuij + guvf*xsupvij xsubvij = gsvf*xsupsij + guvf*xsupuij + gvvf*xsupvij

Parameters
[in]xsupsijContravariant s component.
[in]xsupuijContravariant u component.
[in]xsupvijContravariant v component.
[out]xsubsijCovariant s component.
[out]xsubuijCovariant u component.
[out]xsubvijCovariant v component.
[in]nsminMin radial index.
[in]nsmaxMax radial index.

Definition at line 590 of file metrics.f90.

◆ tolowerh()

subroutine metrics::tolowerh ( real(dp), dimension(nuv,nsmin:nsmax), intent(in)  xsupsij,
real(dp), dimension(nuv,nsmin:nsmax), intent(in)  xsupuij,
real(dp), dimension(nuv,nsmin:nsmax), intent(in)  xsupvij 
)

Converts to covariant component half grid.

Computes covariant (lower) from contravariant (upper) magnetic field components on the half mesh.

xsubsij = gss*xsupsij + gsu*xsupuij + gsv*xsupvij xsubuij = gsu*xsupsij + guu*xsupuij + guv*xsupvij xsubvij = gsv*xsupsij + guv*xsupuij + gvv*xsupvij

Parameters
[in]xsupsijContravariant s component.
[in]xsupuijContravariant u component.
[in]xsupvijContravariant v component.
[out]xsubsijCovariant s component.
[out]xsubuijCovariant u component.
[out]xsubvijCovariant v component.
[in]nsminMin radial index.
[in]nsmaxMax radial index.

Definition at line 539 of file metrics.f90.

◆ toupper()

subroutine metrics::toupper ( real(dp), dimension(nuv,nsmin:nsmax), intent(in)  xsubsij,
real(dp), dimension(nuv,nsmin:nsmax), intent(in)  xsubuij,
real(dp), dimension(nuv,nsmin:nsmax), intent(in)  xsubvij 
)

Converts to contravariant component full grid.

Computes contravariant (upper) from covariant (lower) components on the full mesh.

xsupsij = hss*xsubsij + hsu*xsupuij + hsv*xsubvij xsupuij = hsu*xsubsij + huu*xsupuij + huv*xsubvij xsupvij = hsv*xsubsij + huv*xsupuij + hvv*xsubvij

Parameters
[in]xsubsijCovariant s component.
[in]xsubuijCovariant u component.
[in]xsubvijCovariant v component.
[out]xsupsijContravariant s component.
[out]xsupuijContravariant u component.
[out]xsupvijContravariant v component.
[in]nsminMin radial index.
[in]nsmaxMax radial index.

Definition at line 485 of file metrics.f90.