1 subroutine diagnostics(a, xaxis, yaxis, zvol, zvolp, current)
26 real(r8),
intent(in) :: a, xaxis, yaxis
27 real(r8),
dimension(nr),
intent(out) :: zvol, zvolp
28 real(r8),
intent(out) ::
current
30 real(r8),
dimension(nr) :: avc
31 real(r8),
dimension(nr - 1) :: zjar, zpar, zar, zps, xl
32 real(r8),
dimension(3) :: ablt
34 real(r8) :: area, volume
35 real(r8) :: parea, carea, parea2
36 real(r8) :: grps2,
bp2, bp2vol
37 real(r8) :: r, wr, s, ws
38 real(r8) :: x, xr, xs, y, yr, ys, ps, psr, pss
39 real(r8) :: xjac, arhs
41 real(r8) :: beta, betapl, betastar
42 real(r8) :: xlength, avcur, sumc3, sumc5
43 real(r8) :: bigr, zjdchi, dl, dlv
44 real(r8) :: prs, bussac, err
45 integer(itm_i4) :: i, j, ngr, ngs
46 integer(itm_i4) :: nelm
47 integer(itm_i4) :: n1, n2, n3, n4
50 if (ias == 1) factas = 1._r8
61 call
spline(npts, psivec, p_int, dpres(1) * sign(1._r8, cpsurfin), &
62 dpres(npts) * sign(1._r8, cpsurfin), 1, p_spline)
67 nelm = (i - 1) * (np - 1) + j
77 call
interpolation(1, xx(1, n1), xx(1, n2), xx(1, n3), xx(1, n4), &
79 call
interpolation(1, yy(1, n1), yy(1, n2), yy(1, n3), yy(1, n4), &
81 call
interpolation(1, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) &
82 + 1), psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, ps, &
84 xjac = xr * ys - xs * yr
85 area = area + wr * ws * xjac
86 volume = volume + (1. + eps * x) * wr * ws * xjac
87 grps2 = psr**2 * (xs**2 + ys**2) / xjac**2
88 bp2 = grps2 / (1._r8 + eps * x)**2
89 bp2vol = bp2vol + (1._r8 + eps * x) *
bp2 *wr * ws * xjac
90 parea = parea - wr * ws * xjac *
spwert(npts, ps, p_spline, &
92 parea2 = parea2 + wr * ws * xjac *
spwert(npts, ps, p_spline, &
94 arhs =
rh_side(1._r8, x, 0._r8, 0._r8, 0._r8, 0._r8, ps, 0._r8, &
96 carea = carea + wr * ws * xjac * arhs
101 zjar(i) = -factas * a * carea
102 zpar(i) = -factas * 0.5_r8 * a * b * parea
104 zjar(i) = -factas * a * carea
105 zpar(i) = -factas * eps * a * b * parea
107 zar(i) = factas * abs(area)
108 zvol(nr - i + 1) = factas * abs(volume)
110 area = factas * abs(area)
111 volume = factas * abs(volume)
114 parea = factas * 0.5_r8 * a * b * parea
115 parea2 = factas * (0.5_r8*a*b)**2 * parea2
116 carea = factas * a * carea
118 parea = factas * eps * a * b * parea
119 parea2 = factas * (eps * a * b)**2 * parea2
120 carea = factas * a * carea
123 volume = 2._r8 * pi * volume
124 xli = abs(4._r8 * factas * pi * bp2vol / (rav * carea**2))
125 beta = 2._r8 * (eps / alfa**2) * abs(parea) / area
126 betastar = 2._r8 * (eps / alfa**2) * sqrt(abs(parea2) / area)
127 betapl = (8._r8 * pi / eps) * abs(parea) / carea**2
128 current = abs(eps * carea / alfa)
133 if (standard_output)
then
135 write(out_he, *)
'***************************************'
136 write(out_he, 11) xaxis, yaxis
137 write(out_he, 2) betapl
138 write(out_he, 3) beta
139 write(out_he, 8) betastar
140 write(out_he, 12) 1.256637_r8 * beta /
current * 1.e+6_r8
142 write(out_he, 5) area
143 write(out_he, 6) volume
145 write(out_he, 9) alfa
146 write(out_he, 13) a, b
147 write(out_he, *)
'***************************************'
151 write(txtout(20), 2) betapl
152 write(txtout(21), 3) beta
153 write(txtout(22), 8) betastar
155 write(txtout(24), 5) area
156 write(txtout(25), 6) volume
157 write(txtout(26), 7) xli
158 write(txtout(27), 9) alfa
168 nelm = (i - 1) * (np - 1) + j
174 call
interpolation(1, xx(1, n1), xx(1, n2), xx(1, n3), xx(1, n4), &
176 call
interpolation(1, yy(1, n1), yy(1, n2), yy(1, n3), yy(1, n4), &
178 call
interpolation(1, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) + 1), &
179 psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, ps, &
182 xjac = xr * ys - xs * yr
183 bigr = (1._r8 + eps * x)
184 zjdchi = bigr * xjac / abs(psr)
185 dl = sqrt(xs**2 + ys**2)
186 xlength = xlength + dl * ws
187 sumc5 = sumc5 + ws * zjdchi
188 if (trim(current_averaging) ==
'line')
then
193 sumc3 = sumc3 + ws * dlv
194 arhs =
rh_side(1._r8, x, 0._r8, 0._r8, 0._r8, 0._r8, ps, 0._r8, 0._r8)
195 avcur = avcur + ws * arhs * dlv
199 xl(i) = factas * xlength
200 avc(i) = factas * a * eps * avcur / (alfa * sumc3)
201 zvolp(nr - i + 1) = factas * abs(sumc5)
204 if (standard_output)
then
206 write(out_he, *)
'***********************************************' &
208 write(out_he, *)
'* i psi s <j> error length ' &
210 write(out_he, *)
'***********************************************' &
212 avc(nr) = avc(nr - 1) - (avc(nr - 2) - avc(nr - 1)) &
213 / (zps(nr - 2) - zps(nr - 1)) * zps(nr - 1)
218 prs = 0.5 * a * b *
spwert(npts, ps, p_spline, psivec, &
221 prs = eps * a * b *
spwert(npts, ps, p_spline, psivec, &
224 bussac = -(2._r8 / eps) * (prs - zpar(i) / zar(i)) / (zjar(i) &
226 err = avc(i) / avc(nr) -
j_phi(zps(i))
227 write(out_he, 72) i, zps(i), sqrt(zps(i)), avc(i) / avc(nr), &
230 write(out_he, *)
'***********************************************' &
236 2
format(
' poloidal beta : ', e12.4)
237 3
format(
' toroidal beta : ', e12.4)
238 4
format(
' total current : ', e12.4)
239 5
format(
' total area : ', e12.4)
240 6
format(
' total volume : ', e12.4)
241 7
format(
' int. inductance : ', e12.4)
242 8
format(
' beta star : ', e12.4)
243 9
format(
' pol. flux : ', e12.4)
244 11
format(
' magnetic axis : ', 2f9.5)
245 12
format(
' norm. beta : ', f9.5)
246 13
format(
' a, b : ', 2e12.4)
247 72
format(i3, 3f12.4, es10.2, 2f12.4)
real(r8) function bp2(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine interpolation(type_interp, xn1, xn2, xn3, xn4, r, s, x, xr, xs, xrs, xrr, xss, yn1, yn2, yn3, yn4, pn1, pn2, pn3, pn4, yr, ys, ps)
subroutine allocate_spline_coefficients(spline, n)
subroutine spline(N, X, Y, ALFA, BETA, TYP, A, B, C, D)
REAL *8 function spwert(N, XWERT, A, B, C, D, X, ABLTG)
real(r8) function rh_side(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine current(GEOMETRY, PROFILES, TRANSPORT, SOURCES, EVOLUTION, CONTROL, j_boun, ifail, failstring)
CURRENT TRANSPORT EQUATION.
subroutine deallocate_spline_coefficients(spline)
subroutine set_node_number(n_node, n1, n2, n3, n4)
real(r8) function j_phi(flux)
subroutine diagnostics(a, xaxis, yaxis, zvol, zvolp, current)