The Open FUSION Toolkit 1.0.0-beta6
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
taylor Module Reference

Detailed Description

Subroutines and fields for Taylor state calculations using mimetic operators.

  • Force-Free eigenmodes
  • Vacuum fields for geometries with cut planes
  • Inhomogeneous force free states from vacuum fields
Author
Chris Hansen
Date
June 2010

Data Types

type  oft_taylor_hmodes
 Force-free, uniform \( \lambda \) eigenmode object. More...
type  oft_taylor_ifield
 Force-free, uniform \( \lambda \) field object (inhomogeneous BCs). More...
type  oft_taylor_rinterp
 Interpolate a Taylor state field. More...

Functions/Subroutines

subroutine ff_delete (self, storage_only)
 Setup eigenmodes object.
subroutine ff_setup (self, nh, hcpc, hcpv, htags, ml_hcurl, ml_h1, ml_hcurl_grad, ml_h1grad, ml_lagrange, minlev)
 Setup eigenmodes object.
subroutine hmodes_delete (self, storage_only)
 Setup eigenmodes object.
subroutine hmodes_setup (self, ml_hcurl, ml_lagrange, minlev, htor_axis)
 Setup eigenmodes object.
subroutine taylor_hmodes (self, nm, rst_filename)
 Compute 'taylor_nm' Force-Free eignemodes.
subroutine taylor_injector_single (self, hmodes, lambda, fluxes, gffa)
 Compute force-free plasma response to external fields generated by taylor_vacuum.
subroutine taylor_injectors (self, hmodes, lambda, rst_filename)
 Compute force-free plasma response to external fields generated by taylor_injectors.
subroutine taylor_rinterp (self, cell, f, gop, val)
 Reconstruct a composite Taylor state field.
subroutine taylor_rinterp_setup1 (self, hcurl_grad_rep)
 Setup interpolator for composite Taylor state fields.
subroutine taylor_rinterp_setup2 (self, hcurl_rep, hgrad_rep)
 Setup interpolator for composite Taylor state fields.
subroutine taylor_vac_curr (self, hmodes, rst_filename)
 Generate vector potential whose corresponding current matches the vacuum fields stored in taylor::taylor_hvac.
subroutine taylor_vacuum (self, energy, hmodes, rst_filename)
 Generate vacuum fields for specified handle (cut planes).

Variables

integer(i4), parameter taylor_tag_size = 4
 Size of jump planes character tags.

Function/Subroutine Documentation

◆ ff_delete()

subroutine ff_delete ( class(oft_taylor_ifield), intent(inout) self,
logical, intent(in), optional storage_only )

Setup eigenmodes object.

Parameters
[in,out]selfForce-free eigenmode object
[in]storage_onlyOnly reset storage, but do not clear references

◆ ff_setup()

subroutine ff_setup ( class(oft_taylor_ifield), intent(inout) self,
integer(i4), intent(in) nh,
real(r8), dimension(3,nh), intent(in) hcpc,
real(r8), dimension(3,nh), intent(in) hcpv,
character(len=taylor_tag_size), dimension(nh), intent(in), optional htags,
type(oft_ml_fem_type), intent(in), optional, target ml_hcurl,
type(oft_ml_fem_type), intent(in), optional, target ml_h1,
type(oft_ml_fem_comp_type), intent(in), optional, target ml_hcurl_grad,
type(oft_ml_fem_type), intent(in), optional, target ml_h1grad,
type(oft_ml_fem_type), intent(in), optional, target ml_lagrange,
integer(i4), intent(in), optional minlev )

Setup eigenmodes object.

Parameters
[in,out]selfForce-free field object
[in]nhNumber of jump planes
[in]hcpcJump plane center possitions [3,nh]
[in]hcpvJump plane normal vectors [3,nh]
[in]htagsNames for each jump plane [LEN=taylor_tag_size,nh] (optional)
[in]ml_hcurlMulti-level H(Curl) FE representation
[in]ml_h1Multi-level H^1 FE representation
[in]ml_hcurl_gradMulti-level H(Curl) + grad(H^1) FE representation
[in]ml_h1gradMulti-level grad(H^1) FE representation
[in]ml_lagrangeMulti-level Lagrange FE representation

◆ hmodes_delete()

subroutine hmodes_delete ( class(oft_taylor_hmodes), intent(inout) self,
logical, intent(in), optional storage_only )

Setup eigenmodes object.

Parameters
[in,out]selfForce-free eigenmode object
[in]storage_onlyOnly reset storage, but do not clear references

◆ hmodes_setup()

subroutine hmodes_setup ( class(oft_taylor_hmodes), intent(inout) self,
type(oft_ml_fem_type), intent(in), optional, target ml_hcurl,
type(oft_ml_fem_type), intent(in), optional, target ml_lagrange,
integer(i4), intent(in), optional minlev,
integer(i4), intent(in), optional htor_axis )

Setup eigenmodes object.

Parameters
[in,out]selfForce-free eigenmode object
[in]ml_hcurlMulti-level H(Curl) FE representation
[in]ml_lagrangeMulti-level Lagrange FE representation

◆ taylor_hmodes()

subroutine taylor_hmodes ( type(oft_taylor_hmodes), intent(inout) self,
integer(i4), intent(in), optional nm,
character(len=*), intent(in), optional rst_filename )

Compute 'taylor_nm' Force-Free eignemodes.

Parameters
[in,out]selfForce-free eigenmode object
[in]nmNumber of modes to compute (optional: 1)
[in]rst_filenameFile name to store/load restart information

◆ taylor_injector_single()

subroutine taylor_injector_single ( type(oft_taylor_ifield), intent(inout) self,
type(oft_taylor_hmodes), intent(inout) hmodes,
real(r8), intent(in) lambda,
real(r8), dimension(:), intent(in) fluxes,
class(oft_vector), intent(inout), pointer gffa )

Compute force-free plasma response to external fields generated by taylor_vacuum.

Parameters
[in,out]selfForce-free field object
[in,out]hmodesForce-free eigenmode object
[in]lambdaDesired lambda for force-free state
[in]fluxesFlux for each handle
[in,out]gffaPlasma component (non-vacuum) of handle field

◆ taylor_injectors()

subroutine taylor_injectors ( type(oft_taylor_ifield), intent(inout) self,
type(oft_taylor_hmodes), intent(inout) hmodes,
real(r8), intent(in) lambda,
character(len=*), intent(in), optional rst_filename )

Compute force-free plasma response to external fields generated by taylor_injectors.

Parameters
[in,out]selfForce-free field object
[in,out]hmodesForce-free eigenmode object
[in]lambdaDesired lambda for force-free state
[in]rst_filenameFile name to store/load restart information (optional)

◆ taylor_rinterp()

subroutine taylor_rinterp ( class(oft_taylor_rinterp), intent(inout) self,
integer(i4), intent(in) cell,
real(r8), dimension(:), intent(in) f,
real(r8), dimension(3,4), intent(in) gop,
real(r8), dimension(:), intent(out) val )

Reconstruct a composite Taylor state field.

Parameters
[in]cellCell for interpolation
[in]fPosition in cell in logical coord [4]
[in]gopLogical gradient vectors at f [3,4]
[out]valReconstructed field at f [3]

◆ taylor_rinterp_setup1()

subroutine taylor_rinterp_setup1 ( class(oft_taylor_rinterp), intent(inout) self,
class(oft_fem_comp_type), intent(inout), target hcurl_grad_rep )

Setup interpolator for composite Taylor state fields.

Fetches local representation used for interpolation from vector object

◆ taylor_rinterp_setup2()

subroutine taylor_rinterp_setup2 ( class(oft_taylor_rinterp), intent(inout) self,
class(oft_afem_type), intent(inout), target hcurl_rep,
class(oft_afem_type), intent(inout), target hgrad_rep )

Setup interpolator for composite Taylor state fields.

Fetches local representation used for interpolation from vector object

◆ taylor_vac_curr()

subroutine taylor_vac_curr ( type(oft_taylor_ifield), intent(inout) self,
type(oft_taylor_hmodes), intent(inout) hmodes,
character(len=*), intent(in), optional rst_filename )

Generate vector potential whose corresponding current matches the vacuum fields stored in taylor::taylor_hvac.

Parameters
[in,out]selfForce-free field object
[in,out]hmodesForce-free eigenmode object
[in]rst_filenameFile name to store/load restart information (optional)

◆ taylor_vacuum()

subroutine taylor_vacuum ( type(oft_taylor_ifield), intent(inout) self,
real(r8), dimension(:), intent(out), optional energy,
type(oft_taylor_hmodes), intent(inout), optional hmodes,
character(len=*), intent(in), optional rst_filename )

Generate vacuum fields for specified handle (cut planes).

Parameters
[in,out]selfForce-free field object
[out]energyVacuum energy for each handle (optional)
[in,out]hmodesForce-free eigenmode object (optional)
[in]rst_filenameFile name to store/load restart information (optional)

Variable Documentation

◆ taylor_tag_size

integer(i4), parameter taylor_tag_size = 4

Size of jump planes character tags.