Deals with most types of equilibria to generate GS2's expected geometrical parameters Input: rhoc :: minor radius of surface of interest (see irho below) use_large_aspect :: default false. Choose true to get "s-alpha" equilibria (see local_eq) Equilibrium options : Choose 1 of: local_eq : Analytic parameterization of local equilibrium ppl_eq = .true. :: use PPL style NetCDF equilibrium (Menard) transp_eq = .true. :: use TRXPL style NetCDF equilibrium (Menard) gen_eq = .true. :: use GA style NetCDF equilibrium (TOQ) chs_eq = .true. :: Use CHEASE generated equilibrium efit_eq = .true. :: use EFIT equilibrium dfit_eq = .true. :: use dipole equilibrium -------------------> [note: if any of the above except local_eq are true set eqfile = input file name] irho :: choose flux surface label case (1) :: sqrt(toroidal flux)/sqrt(toroidal flux of LCFS) case (2) :: diameter/(diameter of LCFS). recommended case (3) :: poloidal flux/(poloidal flux of LCFS) case (4) :: rho_mid? Only used for the vacuum ring dipole, dfit_eq = T
equal_arc :: .true. makes grad_par coefficient a constant.
New feature: Given an actual equilibrium, allow s_hat and d pressure/d rho to vary self-consistently (a la Greene and Chance, Bishop). We use the Bishop formalism [C. M. Bishop, et al., NF, Vol. 24, 1579 (1984).]
bishop :: uses Bishop relations to generate coefficients. case(0) :: do not use Bishop relations -- primarily for testing case(1) :: use Bishop relations with actual equilibrium shat, p' case(2) :: use Bishop relations with new shat, and p' from p_prime_input case(3) :: use Bishop relations with new shat, L_p from below s_hat_input :: new magnetic shear, rho/qdq/drho invLp_input :: new -1/pressure * d pressure/d rho case(4) :: use Bishop relations with new shat, beta' from below s_hat_input :: new magnetic shear, rho/qdq/drho beta_prime_input :: new d beta/d rho case(5) :: use Bishop relations with new shat, alpha from below s_hat_input :: new magnetic shear, rho/qdq/drho alpha_input :: new alpha = -q*2 * rmaj * d pressure/d rho case(6) :: use Bishop relations with equilibrium shat, beta' from below beta_prime_input :: new d beta/d rho case(7) :: use Bishop relations with equilibrium shat, beta' from below dp_mult :: d pressure/d rho = equilibrium gradient * dp_mult (In each case, the definition of rho is determined by irho above)
nperiod :: number of cells of width 2 pi in theta.
geoType :: determines the analytic geometry specification case (0) :: the tilted Miller geometry: generalization of Miller et al. PoP 1998 (see Chapter 3.2.1 of J. Ball MIT Masters thesis) case (1) :: the global MHD solution geometry: calculated by inverting the large aspect ratio, constant-current MHD solution (see J. Ball "Optimized up-down asymmetry to drive fast intrinsic rotation in tokamak reactors") case (2) :: the generalized ellipticity geometry: the polar formula for an ellipse generalized to arbitrary shaping effect mode number (see Chapter 5.1.2 of J. Ball Oxford PhD thesis) case (3) :: the Fourier series geometry: Fourier decomposition of the minor radius as a function of poloidal angle (see Chapter 5.1.1 of J. Ball Oxford PhD thesis) rmaj :: major radius of magnetic axis (in units of aminor) if not local_eq or major radius of center of flux surface of interest (units of aminor) otherwise r_geo :: major radius of last closed flux surface or if local_eq: sets the reference magnetic field to be the value of the toroidal field at R=r_geo on the flux surface of interest aSurf :: the real-space flux surface label of the surface with deltam and deltan shaping as well as a major radius of Rmaj (only used when geoType={1}) shift :: d rmaj/drho. (Typically < 0.) (only used when geoType={0,2,3}) -> Alternative definitions: Raxis :: major radial location of the magnetic axis. (only used when geoType={1}) shiftVert :: d zmaj/drho. (< 0 for a vertical shift) (only used when geoType={0,2,3}) -> Alternative definitions: Zaxis :: axial location of the magnetic axis. (only used when geoType={1}) mMode :: the first poloidal mode number (only used when geoType={2,3}) nMode :: the second poloidal mode number (only used when geoType={2,3}) deltam :: magnitude of the first poloidal shaping effect (only used when geoType={2,3}) -> Alternative definitions: delta2 :: identical to akappa (only used when geoType={1}) akappa :: flux surface elongation (only used when geoType={0}) deltan :: magnitude of the second poloidal shaping effect (only used when geoType={2,3}) -> Alternative definitions: delta3 :: magnitude of the triangulaity-like shaping effect (only used when geoType={1}) tri :: flux surface triangulaity (only used when geoType={0}) deltampri :: d(deltam)/d(rhoc) (only used when geoType={2,3}) -> Alternative definitions: akappapri :: d(akappa)/d(rhoc) (only used when geoType={0}) deltanpri :: d(deltan)/d(rhoc) (only used when geoType={2,3}) -> Alternative definitions: tripri :: d(tri)/d(rhoc) (only used when geoType={0}) thetam :: tilt angle of the first poloidal shaping effect (only used when geoType={2,3}) -> Alternative definitions: theta2 :: tilt angle of the flux surface elongation (only used when geoType={1}) thetak :: tilt angle of the flux surface elongation (only used when geoType={0}) thetan :: tilt angle of the first poloidal shaping effect (only used when geoType={2,3}) -> Alternative definitions: theta3 :: tilt angle of the triangularity-like shaping effect (only used when geoType={1}) thetad :: tilt angle of the flux surface triangularity (only used when geoType={0}) qinp :: safety factor q at surface of interest (local_eq = T only) shat :: d q/drho at surface of interest
delrho :: numerical parameter. Should be "small enough". Typically 0.001 ok.
force_sym :: true -> assume up-down symmetric equilibrium. writelots :: .true. for more output
Use: Include this module, set input variables appropriately and call eikcoefs.
I do not recommend calling other routines, but a few are available for diagnostic purposes.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
class(abstract_geo_type), | public, | allocatable | :: | geom | |||
type(flux_surface_type), | public | :: | surf | ||||
real, | public | :: | p_prime_input | ||||
real, | public | :: | invLp_input | ||||
real, | public | :: | beta_prime_input | ||||
real, | public | :: | alpha_input | ||||
real, | public | :: | dp_mult | ||||
integer, | public | :: | irho | ||||
integer, | public | :: | bishop | ||||
logical, | public | :: | gen_eq | ||||
logical, | public | :: | efit_eq | ||||
logical, | public | :: | ppl_eq | ||||
logical, | public | :: | local_eq | ||||
logical, | public | :: | chs_eq | ||||
logical, | public | :: | dfit_eq | ||||
logical, | public | :: | idfit_eq | ||||
logical, | public | :: | gs2d_eq | ||||
logical, | public | :: | transp_eq | ||||
logical, | public | :: | writelots | ||||
logical, | public | :: | equal_arc | ||||
logical, | public | :: | force_sym | ||||
logical, | public | :: | use_large_aspect | ||||
character(len=EQFILE_LENGTH), | public | :: | eqfile | ||||
character(len=EQFILE_LENGTH), | public | :: | eqnormfile | ||||
real, | private | :: | rpmin | ||||
real, | private | :: | rpmax | ||||
integer, | public | :: | big | = | 8 | ||
integer, | private | :: | nth |
The single-period theta domain "half-grid" size. That is, the number of points in
. |
|||
integer, | private | :: | ntgrid |
The full theta domain "half-grid" size. That is, the total size of the
geometry grid is |
|||
integer, | public | :: | verb | = | 2 |
Verbosity of print statements | Should we use runtime_tests:verbosity? |
|
logical, | private | :: | debug | = | .false. |
If debug then print out lots of debug info. debug is true if (verb >= 3) |
A type to hold a collection of geometrical quantities as calculated by eikcoefs
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | ntheta | ||||
real, | public | :: | drhodpsi | ||||
real, | public | :: | kxfac | ||||
real, | public | :: | qsf | ||||
real, | public | :: | surfarea | ||||
real, | public | :: | dvdrhon | ||||
real, | public | :: | shat | ||||
real, | public | :: | dbetadrho | ||||
real, | public | :: | rhoc | ||||
real, | public, | dimension(:), allocatable | :: | theta |
theta :: the theta grid of the results. The final grid has pi - 2pinperiod <= theta <= 2pinperiod - pi |
||
real, | public, | dimension(:), allocatable | :: | bmag |
bmag :: |B| |
||
real, | public, | dimension(:), allocatable | :: | gradpar |
gradpar :: coefficient of b_hat.grad operator |
||
real, | public, | dimension(:), allocatable | :: | grho |
grho :: grad(rho) |
||
real, | public, | dimension(:), allocatable | :: | cdrift |
cdrift :: part of coriolis drift operator independent of theta_0 |
||
real, | public, | dimension(:), allocatable | :: | cvdrift |
cvdrift :: coefficient of v_par**2/Omega b_hat x (b_hat.grad b_hat).grad |
||
real, | public, | dimension(:), allocatable | :: | gbdrift |
gbdrift :: coefficient of mu/Omega b_hat x (grad B).grad operator |
||
real, | public, | dimension(:), allocatable | :: | cdrift0 |
cdrift0 :: part of coriolis drift operator proportional to theta_0 |
||
real, | public, | dimension(:), allocatable | :: | cvdrift0 |
cvdrift0 :: part of total curvature drift operator proportional to theta_0 |
||
real, | public, | dimension(:), allocatable | :: | gbdrift0 |
gbdrift0 :: part of total grad B drift operator proportional to theta_0 |
||
real, | public, | dimension(:), allocatable | :: | gds2 |
gds2 :: part of grad_perp**2 operator independent of theta_0 |
||
real, | public, | dimension(:), allocatable | :: | gds21 |
gds21 :: part of grad_perp**2 operator proportional to theta_0 |
||
real, | public, | dimension(:), allocatable | :: | gds22 |
gds22 :: part of grad_perp2 operator proportional to theta_02 |
||
real, | public, | dimension(:), allocatable | :: | gds23 |
gds23 :: part of v_E . grad theta operator independent of theta_0 |
||
real, | public, | dimension(:), allocatable | :: | gds24 |
gds24 :: part of v_E . grad theta operator proportional to theta_0 |
||
real, | public, | dimension(:), allocatable | :: | gds24_noq |
gds24_noq :: gds24/dqdrp |
||
real, | public, | dimension(:), allocatable | :: | rplot |
rplot :: R(theta) - major radius on theta grid |
||
real, | public, | dimension(:), allocatable | :: | zplot |
zplot :: Z(theta) - height on theta grid |
||
real, | public, | dimension(:), allocatable | :: | aplot |
aplot :: a(theta) - alpha-phi on theta grid |
||
real, | public, | dimension(:), allocatable | :: | rprime |
rprime :: derivative wrt rho of R(theta) - major radius on theta grid |
||
real, | public, | dimension(:), allocatable | :: | zprime |
zprime :: derivative wrt rho Z(theta) - height on theta grid |
||
real, | public, | dimension(:), allocatable | :: | aprime |
aprime :: derivative wrt rho a(theta) - alpha-phi on theta grid |
||
real, | public, | dimension(:), allocatable | :: | bpol |
bpol :: Poloidal magnetic field strength |
Returns the diameter of the surface labelled by rp.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rp |
FIXME : Add documentation
Only used with ~undocumented irho = 4
option. Ideally would update to take pbar
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rp |
Calculates the normalised radial coordinate
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rp |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | r | |||
real, | intent(in) | :: | theta |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rho |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rho |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rp |
Calculate the derivative of the flux label rho with respect to the internal flux label rp (the normalised poloidal flux)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rp | |||
real, | intent(in) | :: | dr |
FIXME : Add documentation | Only used for writelots output and if irho=1
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rp | |||
real, | intent(in) | :: | drhodrp | |||
real, | intent(in) | :: | dr |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-nth:nth) | :: | f | ||
real, | intent(in), | dimension(-nth:nth) | :: | theta | ||
real, | intent(in), | dimension(-nth:nth) | :: | jacob |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-nth:nth) | :: | b_mag | ||
real, | intent(in), | dimension(-nth:nth) | :: | theta | ||
real, | intent(in), | dimension(-nth:nth) | :: | jacob |
Run eikcoefs with ntheta = ntheta_input
and then resample onto
a theta grid with resolution ntheta_target
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(inout) | :: | ntheta_input |
This is the number of theta grid points per two pi period used in the geometry calculation (actually one less than this). |
||
integer, | intent(in) | :: | ntheta_target |
This is the desired number of theta grid points per two pi period in the resulting geometrical coefficients used elsewhere. |
||
integer, | intent(in) | :: | nperiod | |||
type(eikcoefs_output_type), | intent(out) | :: | results | |||
integer, | intent(in), | optional | :: | job_id |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(inout) | :: | ntheta_input |
Sets the number of theta grid points. If odd will be updated to be ntheta_input - 1 |
||
integer, | intent(in) | :: | nperiod | |||
type(eikcoefs_output_type), | intent(out) | :: | outputs |
A type containing all of the outputs of eikcoefs |
||
integer, | intent(in), | optional | :: | job_id_in |
job_id is the job id in a Trinity or list mode run. Only for debug |
Allocate the geom instance to the correct type
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(abstract_geo_type), | intent(inout), | allocatable | :: | geom |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:, :) | :: | gradstot | ||
real, | intent(in), | dimension(-ntgrid:, :) | :: | gradrptot | ||
real, | intent(in), | dimension(-ntgrid:, :) | :: | gradthtot | ||
real, | intent(in) | :: | dpsidrho | |||
real, | intent(in) | :: | dqdrp | |||
real, | intent(in) | :: | qval | |||
real, | intent(in) | :: | rhoc | |||
real, | intent(in), | dimension(-ntgrid:) | :: | bmag | ||
logical, | intent(in) | :: | force_sym | |||
real, | intent(out), | dimension(-ntgrid:) | :: | gds2 | ||
real, | intent(out), | dimension(-ntgrid:) | :: | gds21 | ||
real, | intent(out), | dimension(-ntgrid:) | :: | gds22 | ||
real, | intent(out), | dimension(-ntgrid:) | :: | gds23 | ||
real, | intent(out), | dimension(-ntgrid:) | :: | gds24 | ||
real, | intent(out), | dimension(-ntgrid:) | :: | gds24_noq |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(in), | dimension(-ntgrid:) | :: | gradpar | ||
real, | intent(in), | dimension(-ntgrid:) | :: | bmag | ||
real, | intent(in), | dimension(-ntgrid:) | :: | grho | ||
real, | intent(in) | :: | dpsidrho | |||
real, | intent(in) | :: | rhoc | |||
integer, | intent(in) | :: | nth | |||
real, | intent(out) | :: | surfarea | |||
real, | intent(out) | :: | dvdrhon |
Computes something like x cross y dot grad zeta.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:, :) | :: | x | ||
real, | intent(in), | dimension(-ntgrid:, :) | :: | y | ||
real, | intent(in), | dimension(-ntgrid:) | :: | invR | ||
real, | intent(out), | dimension(-ntgrid:) | :: | val |
Calculate rgrid for the case of eeq and deq: in this case, rgrid is the distance to the magnetic axis as a function of theta.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(out), | dimension(-ntgrid:) | :: | rgrid | ||
real, | intent(in) | :: | rp | |||
real, | intent(in), | dimension(-ntgrid:) | :: | theta |
Calculates rmajor = R(theta), where R is the distance to the central axis (rmajor not to be confused with rmaj, r_geo etc) rgrid is the submodule radial grid (r_circ for EFIT & deq, rp otherwise)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | rgrid | ||
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(out), | dimension(-ntgrid:) | :: | rmajor |
Calculates invR = 1/R(theta), where R is the distance to the central axis (rmajor not to be confused with rmaj, r_geo etc) rgrid is the submodule radial grid (r_circ for EFIT & deq, rp otherwise)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | rgrid | ||
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(out), | dimension(-ntgrid:) | :: | invR |
Calculates the vertical height, Z, on the theta grid
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | rgrid | ||
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(out), | dimension(-ntgrid:) | :: | Zoftheta |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid: ) | :: | invR | ||
integer, | intent(in) | :: | nth | |||
real, | intent(in), | dimension(-ntgrid:,:) | :: | gradrptot | ||
real, | intent(in) | :: | dpsidrp | |||
real, | intent(in) | :: | bi | |||
real, | intent(out), | dimension(-ntgrid:,:) | :: | bvector |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | invR | ||
real, | intent(in), | dimension(-ntgrid:, :) | :: | grads | ||
real, | intent(out), | dimension(-ntgrid:, :) | :: | gradstot |
Calculate the cross product c = a x b on the [-nth, nth]
domain
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:, :) | :: | a | ||
real, | intent(in), | dimension(-ntgrid:, :) | :: | b | ||
real, | intent(out), | dimension(-ntgrid:, :) | :: | c |
Calculate the dot product c = a.b where a, b and c are on the [-ntgrid, ntgrid]
domain
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:,:) | :: | a | ||
real, | intent(in), | dimension(-ntgrid:,:) | :: | b | ||
real, | intent(out), | dimension(-ntgrid:) | :: | c |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension (-ntgrid:) | :: | rgrid | ||
real, | intent(in) | :: | rp | |||
real, | intent(in), | dimension (-ntgrid:, :) | :: | gradstot | ||
real, | intent(in), | dimension (-ntgrid:, :) | :: | gradrptot | ||
real, | intent(in), | dimension (-ntgrid:, :) | :: | gradztot | ||
real, | intent(in), | dimension (-ntgrid:, :) | :: | gradbtot | ||
real, | intent(in) | :: | dqdrp | |||
real, | intent(in) | :: | dpsidrho | |||
real, | intent(in) | :: | dpsidrp | |||
real, | intent(in), | dimension (-ntgrid:) | :: | bmag | ||
real, | intent(in), | dimension (-ntgrid:) | :: | invR | ||
real, | intent(in), | dimension (-ntgrid:) | :: | theta | ||
real, | intent(in) | :: | dbetadrho | |||
real, | intent(in) | :: | bi | |||
logical, | intent(in) | :: | force_sym | |||
integer, | intent(in) | :: | nperiod | |||
real, | intent(out), | dimension(-ntgrid:ntgrid) | :: | gbdrift | ||
real, | intent(out), | dimension(-ntgrid:ntgrid) | :: | gbdrift0 | ||
real, | intent(out), | dimension(-ntgrid:ntgrid) | :: | cvdrift | ||
real, | intent(out), | dimension(-ntgrid:ntgrid) | :: | cvdrift0 | ||
real, | intent(out), | dimension(-ntgrid:ntgrid) | :: | cdrift | ||
real, | intent(out), | dimension(-ntgrid:ntgrid) | :: | cdrift0 |
Forces symmetry in theta in the passed array. If isign == 0 then
a(+/- theta)
is replaced by the average of a(+/-theta)
, otherwise
a(+/- theta)
is replaced by +/- (a(theta)-a(-theta))
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(inout), | dimension(-ntgrid:) | :: | a | ||
integer, | intent(in) | :: | isign | |||
integer, | intent(in) | :: | ntgrid |
This is not right for cart types, but ok for those compatible with bishop=0.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in) | :: | rp | |||
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(in), | dimension(-ntgrid:) | :: | invR | ||
real, | intent(inout) | :: | qval | |||
real, | intent(in) | :: | bi | |||
integer, | intent(in) | :: | nperiod | |||
real, | intent(out), | dimension(-ntgrid:) | :: | seik |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(in), | dimension(-ntgrid:) | :: | invR | ||
real, | intent(in), | dimension(-ntgrid:) | :: | trip | ||
real, | intent(inout) | :: | qval | |||
real, | intent(in) | :: | bi | |||
real, | intent(out), | dimension(-ntgrid:) | :: | seik | ||
real, | intent(out), | dimension(-ntgrid:) | :: | dsdthet | ||
real, | intent(out) | :: | dpsidrp |
Integrates arg(grid) from i = 0 (generally theta = 0) to either
ends of the grid (usually +/- pi) storing the cumulative integral
in the output ans
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | arg | ||
real, | intent(in), | dimension(-ntgrid:) | :: | grid | ||
real, | intent(out), | dimension(-ntgrid:) | :: | ans | ||
integer, | intent(in) | :: | n |
Find soln where f(soln) ~= fval. Works best if a and b bracket the root. We can then first bisect the bracket to narrow in on the solution. Following this we use the secant method to find the solution. The input xerrbi influences the number of bisection iterations whilst xerrsec is used to set the stopping tolerance on the secant method.
Type | Intent | Optional | Attributes | Name | ||||||||||||||||
---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|
private function f(x_)Arguments
Return Value real |
||||||||||||||||||||
real, | intent(in) | :: | fval | |||||||||||||||||
real, | intent(in) | :: | a | |||||||||||||||||
real, | intent(in) | :: | b | |||||||||||||||||
real, | intent(in) | :: | xerrbi | |||||||||||||||||
real, | intent(in) | :: | xerrsec | |||||||||||||||||
integer, | intent(out) | :: | ier | |||||||||||||||||
real, | intent(out) | :: | soln |
Given an initial theta and gradpar (=dtheta/dl) compute a new theta grid, returned in arcl, which has gradpar = dg/dl = constant.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ntheta | |||
integer, | intent(in) | :: | nperiod | |||
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(inout), | dimension(-ntgrid:) | :: | gpar | ||
real, | intent(out), | dimension(-ntgrid:) | :: | arcl |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension (-ntgrid:) | :: | rgrid | ||
real, | intent(in), | dimension (-ntgrid:) | :: | theta | ||
real, | intent(out), | dimension (-ntgrid:) | :: | ltheta |
Compute df / dltheta using a periodic central difference. The
input ext
can be used to control the effective jump in f
between both ends of the "periodic" grid.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension (-ntgrid:) | :: | ltheta | ||
real, | intent(in), | dimension (-ntgrid:) | :: | f | ||
real, | intent(out), | dimension (-ntgrid:) | :: | dfdl | ||
real, | intent(in) | :: | ext | |||
integer, | intent(in) | :: | n |
Finds the angle between the surface tangent vector in the poloidal plane and the direction of increasing major radius.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension (-ntgrid:, :) | :: | rpgrad | ||
real, | intent(out), | dimension (-ntgrid:) | :: | th_bish | ||
integer, | intent(in) | :: | nth |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension (-ntgrid:) | :: | invR | ||
real, | intent(in), | dimension (-ntgrid:, :) | :: | rpgrad | ||
real, | intent(in) | :: | dpsidrp | |||
real, | intent(out), | dimension (-ntgrid:) | :: | bpol | ||
integer, | intent(in) | :: | nth |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension (-ntgrid:) | :: | invR | ||
real, | intent(in), | dimension (-ntgrid:) | :: | bpol | ||
real, | intent(out), | dimension(-ntgrid:) | :: | bmag | ||
real, | intent(in) | :: | bi | |||
integer, | intent(in) | :: | nth |
General calculation of the poloidal magnetic field radius of curvature (see section 2.3 of "GS2 analytic geometry specification") ! JRB
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | theta | ||
real, | intent(in), | dimension(-ntgrid:) | :: | th_bish | ||
real, | intent(in), | dimension(-ntgrid:) | :: | ltheta | ||
real, | intent(out), | dimension(-ntgrid:) | :: | invRpol | ||
integer, | intent(in) | :: | nth |
FIXME : Add documentation
Initialse geometry to a uniform grid with 2*((2*nperiod - 1)*(nt/2)) + 1
points. Also set nth, geometry and [[geometry::ntgrid]]
consistently; importantly nth is forced to be even. If ntheta
is odd, then
nth == ntheta - 1
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(inout) | :: | ntheta | |||
real, | intent(inout), | dimension(:), allocatable | :: | theta | ||
integer, | intent(in) | :: | nperiod | |||
integer, | intent(out) | :: | ntgrid | |||
integer, | intent(out) | :: | nth |
Copy the [-nth,nth] = [-pi,pi]
section of the passed array into
the theta > +/- pi
parts, adding on k * ext
(with k an integer
labelling the sections/periods). Here we expect ext == array(nth) - array(-nth)
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(inout), | dimension(-ntgrid:) | :: | a | ||
real, | intent(in) | :: | ext | |||
integer, | intent(in) | :: | nperiod |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(-ntgrid:) | :: | bmag | ||
real, | intent(in), | dimension(-ntgrid:) | :: | bpol | ||
real, | intent(in), | dimension(-ntgrid:) | :: | invRpol | ||
real, | intent(in), | dimension(-ntgrid:) | :: | th_bish | ||
real, | intent(in), | dimension(-ntgrid:) | :: | ltheta | ||
real, | intent(in), | dimension(-ntgrid:) | :: | rmajor | ||
real, | intent(in) | :: | dp | |||
real, | intent(in) | :: | bi | |||
real, | intent(out), | dimension(-ntgrid:,:) | :: | gradB |
Perform minor checks on consistency of flags
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
logical, | intent(in) | :: | geq | |||
logical, | intent(in) | :: | eeq | |||
logical, | intent(in) | :: | peq | |||
logical, | intent(in) | :: | leq | |||
logical, | intent(in) | :: | deq | |||
logical, | intent(in) | :: | ideq | |||
logical, | intent(in) | :: | ceq | |||
logical, | intent(in) | :: | gs2deq | |||
logical, | intent(in) | :: | teq | |||
integer, | intent(in) | :: | bishop |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension (-ntgrid:) | :: | rgrid | ||
real, | intent(in), | dimension (-ntgrid:) | :: | seik | ||
real, | intent(in), | dimension (-ntgrid:, :) | :: | grads | ||
real, | intent(in) | :: | dpsidrho | |||
real, | intent(in) | :: | dr | |||
real, | intent(in), | dimension (-ntgrid:) | :: | theta | ||
real, | intent(out), | dimension (-ntgrid:) | :: | rplot | ||
real, | intent(out), | dimension (-ntgrid:) | :: | zplot | ||
real, | intent(out), | dimension (-ntgrid:) | :: | aplot | ||
real, | intent(out), | dimension (-ntgrid:) | :: | rprime | ||
real, | intent(out), | dimension (-ntgrid:) | :: | zprime | ||
real, | intent(out), | dimension (-ntgrid:) | :: | aprime |