1 subroutine helena(equilibrium_in, equilibrium_out, in_path, code_parameters)
139 type (type_equilibrium
),
pointer :: equilibrium_in(:)
140 type (type_equilibrium
),
pointer :: equilibrium_out(:)
141 character(len = 132),
optional :: in_path
142 type (type_param
) :: code_parameters
144 type (type_equilibrium
) :: equilibrium
148 real(r8),
parameter :: rel_change_w_mhd = 1.e-2_r8
149 integer(itm_i4),
parameter :: max_points_boundary = 1025
151 real(r8),
allocatable :: zvol(:), zvolp(:)
152 real(r8),
allocatable :: drmerc(:), dimerc(:), hh(:)
153 real(r8),
allocatable :: qprof(:), dqprof(:), geonc(:), zjpar(:)
154 real(r8),
allocatable :: fcirc(:), b02av(:), b0max(:), rav(:)
155 real(r8),
allocatable :: x_work(:), y_work(:), work(:), theta(:)
156 real(r8),
allocatable :: fr(:)
157 real(r8),
dimension(:),
pointer :: frmajor, fz
158 real(r8),
allocatable :: dummy1(:)
159 real(r8),
allocatable :: pskn1(:)
160 real(r8),
allocatable :: xx_center(:), yy_center(:)
162 real(r8),
dimension(3) :: abltg
164 real(r8) :: a, a_input
166 real(r8) :: rgeo, amin
170 real(r8) :: xaxis, yaxis, psaxis
173 real(r8) :: q95out, q1out
174 real(r8) ::
current, mhd_energy, volume
175 real(r8) :: w_mhd_ref = 0._r8, b_ref
177 integer(itm_i4) :: i, j, m
178 integer(itm_i4) :: nax
179 integer(itm_i4) :: n_bnd
180 integer(itm_i4) :: meshno = 0
181 integer(itm_i4) :: nelm, no, n1, n2, n3, n4
182 integer(itm_i4) :: end_of_path
184 logical :: converged_current = .false.
185 logical :: npts_set = .false.
188 integer(itm_i4) :: return_status
190 character(len = 132) :: codename(1) =
'HELENA'
191 character(len = 132) :: codeversion(1) =
'325'
193 character(len = 2) :: file_no
195 if (.not. present(in_path))
then
200 end_of_path = index(path, achar(0)) - 1
201 if (end_of_path <= 0)
then
202 end_of_path = len_trim(adjustl(path))
203 if (end_of_path == 0)
then
206 path = trim(adjustl(path))
209 path = trim(adjustl(path(1 : end_of_path)))
213 write(*, *)
'pre: helena_initialization'
222 write(*, *)
'pre: gaussian_points'
228 write(*, *)
'post: gaussian_points'
232 if (
associated(equilibrium_in(1)%eqgeometry%boundary(1)%r) &
233 .and.
associated(equilibrium_in(1)%eqgeometry%boundary(1)%r))
then
234 frmajor => equilibrium_in(1)%eqgeometry%boundary(1)%r
235 fz => equilibrium_in(1)%eqgeometry%boundary(1)%z
237 do n_bnd =
size(frmajor), 1, -1
238 if (frmajor(n_bnd) /= 0._r8 .or. fz(n_bnd) /= 0._r8)
exit
243 write(iu6, *)
'ERROR: boundary has only 1 point'
246 if (frmajor(1) == frmajor(n_bnd) .and. fz(1) == fz(n_bnd))
then
248 write(iu6, *)
'removed duplicated point from boundary'
253 write(iu6, *)
'n_bnd = ', n_bnd
254 rvac = (maxval(frmajor(1 : n_bnd)) + minval(frmajor(1 : n_bnd))) &
256 rminor = (maxval(frmajor(1 : n_bnd)) - minval(frmajor(1 : n_bnd))) &
258 write(*,*)
'Guessed rvac, rminor = ', rvac, rminor
260 write(iu6, *)
'ERROR: no plasma boundary specified'
266 if (
associated(equilibrium_in(1)%profiles_1d%psi))
then
267 cpsurfin = equilibrium_in(1)%profiles_1d%psi( &
268 size(equilibrium_in(1)%profiles_1d%psi)) &
269 - equilibrium_in(1)%profiles_1d%psi(1)
273 if (
exists_entry(equilibrium_in(1)%global_param%psi_ax))
then
274 psi_axis = equilibrium_in(1)%global_param%psi_ax
280 if (
exists_entry(equilibrium_in(1)%global_param%psi_bound) &
281 .and.
exists_entry(equilibrium_in(1)%global_param%psi_ax))
then
282 cpsurfin = (equilibrium_in(1)%global_param%psi_bound &
283 - equilibrium_in(1)%global_param%psi_ax)
285 if (
associated(equilibrium_in(1)%profiles_1d%psi))
then
286 equilibrium_in(1)%global_param%psi_ax &
287 = equilibrium_in(1)%profiles_1d%psi(1)
288 equilibrium_in(1)%global_param%psi_bound &
289 = equilibrium_in(1)%profiles_1d%psi(
size( &
290 equilibrium_in(1)%profiles_1d%psi))
293 if (
exists_entry(equilibrium_in(1)%eqgeometry%geom_axis%r)) &
294 rvac = equilibrium_in(1)%eqgeometry%geom_axis%r
295 if (
exists_entry(equilibrium_in(1)%eqgeometry%a_minor)) &
296 rminor = equilibrium_in(1)%eqgeometry%a_minor
298 if (
exists_entry(equilibrium_in(1)%global_param%mag_axis%bphi)) &
299 bvac = equilibrium_in(1)%global_param%mag_axis%bphi
300 if (
exists_entry(equilibrium_in(1)%global_param%toroid_field%r0) &
301 .and.
exists_entry(equilibrium_in(1)%global_param%toroid_field%b0)) &
302 bvac = equilibrium_in(1)%global_param%toroid_field%r0 &
303 * equilibrium_in(1)%global_param%toroid_field%b0 / rvac
304 if (
exists_entry(equilibrium_in(1)%global_param%i_plasma)) &
305 ip = equilibrium_in(1)%global_param%i_plasma
307 allocate(equilibrium_out(1))
308 allocate(equilibrium_out(1)%codeparam%codename(1))
309 allocate(equilibrium_out(1)%codeparam%codeversion(1))
310 if (.not.
associated(code_parameters%parameters))
then
311 write(iu6, *)
'ERROR: code parameters not associated!'
314 allocate(equilibrium_out(1)%codeparam%parameters(
size( &
315 code_parameters%parameters)))
319 equilibrium_out(1)%datainfo%cocos = 13
320 equilibrium_out(1)%codeparam%codename = codename
321 equilibrium_out(1)%codeparam%codeversion = codeversion
322 equilibrium_out(1)%codeparam%parameters = code_parameters%parameters
324 write(*, *)
'pre: helena_assign_code_parameters'
330 if (return_status /= 0)
then
331 write(iu6, *)
'ERROR: Could not assign code parameters.'
336 if (
exists_entry(equilibrium_in(1)%global_param%toroid_field%r0))
then
337 r_vessel = equilibrium_in(1)%global_param%toroid_field%r0
342 if (trim(match) ==
"Ip" .and. ip == 0._r8)
then
343 write(iu6, *)
"Error: for match == 'Ip' the total current must be" &
348 if (verbosity > 0)
then
349 write(iu6, *)
'******************************************************'
350 write(iu6, *)
'* HELENA *'
351 write(iu6, *)
'******************************************************'
354 if (verbosity > 0)
write(iu6, *)
'done assigning code parameters'
357 if (
output /=
'none' .and. (standard_output .or. elite_output &
358 .or. profiles_output .or. additional_output .or. xmgrace_output))
then
359 call system(
'mkdir -p ' // trim(adjustl(path)))
360 write(iu6, *)
'created run directory'
364 if (xmgrace_output)
then
365 call system(
'mkdir -p ' // trim(adjustl(path)) //
'plots')
366 write(iu6, *)
'created plots directory'
371 if (ias == 1 .and. mod(nchi, 2) /= 0)
then
372 write(iu6, *)
'Error: nchi must be 2^n'
375 if (ias == 0 .and. mod(nchi, 2) /= 1)
then
376 write(iu6, *)
'Error: nchi must be 2^n + 1'
380 if (n_acc_points < 0 .or. n_acc_points > 10)
then
381 write(iu6, *)
'Error: condition 0 <= n_acc_points <= 10 not fulfilled!'
384 if (equidistant < 0 .or. equidistant > 1)
then
385 write(iu6, *)
'Error: condition 0 <= equidistant <= 1 not fulfilled!'
389 if (te%shape > 2)
then
390 write(iu6, *)
'Error: 0 <= te%shape <= 2 !'
393 if (ne%shape > 2)
then
394 write(iu6, *)
'Error: 0 <= ne%shape <= 2 !'
398 if (wmhd < 0._r8)
then
399 write(iu6, *)
'Error: wmhd must be positive !'
403 if (verbosity > 1)
then
404 write(iu6, *)
'cpsurfin ', cpsurfin
405 write(iu6, *)
'rvac ', rvac
406 write(iu6, *)
'rminor ', rminor
407 write(iu6, *)
'eps ', eps
408 write(iu6, *)
'bvac ', bvac
411 if (verbosity > 2)
then
412 write(iu6, *)
'W_MHD = ', wmhd
413 write(iu6, *)
'mfm = ', mfm
419 write(iu6, *)
'ERROR: mfm not a power of 2'
426 write(*, *)
'pre: shape'
434 allocate(fr(mfm + 2))
437 if (verbosity > 1)
write(iu6, *)
'shape from discrete points'
447 theta(m) = 2._r8 * pi * (m - 1) / dble(mfm - 1)
450 if (xmgrace_output)
then
451 call
plot_data_1d(
'line', theta, fr, mfm,
'xmgr_boundary_orig')
458 call
rft2(fr, mfm, 1)
461 fr(m) = 2._r8 * fr(m) / dble(mfm)
469 write(*, *)
'post: shape'
474 if (verbosity > 2)
write(iu6, *)
'minor radius a = ', rminor
476 if (standard_output) &
477 open (unit = out_he, file = trim(adjustl(path)) // file_he_out, &
478 status =
'replace', form =
'formatted', action =
'write', iostat = i_error)
482 a = 2._r8 * (1._r8 + (1._r8 - eps**2 / 4._r8) / ellip**2)
483 b = 2._r8 * eps + tria / ellip**2 / (1._r8 + (1._r8 - eps**2 &
485 cpsurf = 0.5_r8 * eps**2 / (1._r8 + eps**2)**2 * ellip &
486 / sqrt(1._r8 - eps**2 / 4._r8)
487 radius = sqrt(eps**2 / (1._r8 + eps**2))
488 b0 = sqrt(1._r8 + eps**2)
489 alfa = radius**2 * b0 / cpsurf
498 if (
associated(equilibrium_in(1)%profiles_1d%pprime))
then
499 if (verbosity > 1)
write(iu6, *)
'pprime specified on input'
501 npts =
size(equilibrium_in(1)%profiles_1d%pprime)
503 dp = equilibrium_in(1)%profiles_1d%pprime
505 if (
associated(equilibrium_in(1)%profiles_1d%pressure))
then
506 if (npts_set .and.
size(equilibrium_in(1)%profiles_1d%pressure) &
508 write(iu6, *)
'ERROR: input profiles have inconsistent sizes'
511 if (verbosity > 1)
write(iu6, *)
'pressure specified on input'
513 npts =
size(equilibrium_in(1)%profiles_1d%pressure)
514 p_sep = equilibrium_in(1)%profiles_1d%pressure(npts)
519 if (
associated(equilibrium_in(1)%profiles_1d%ffprime))
then
520 if (npts_set .and.
size(equilibrium_in(1)%profiles_1d%ffprime) &
522 write(iu6, *)
'ERROR: input profiles have inconsistent sizes'
525 if (verbosity > 1)
write(iu6, *)
'ffprime specified on input'
527 npts =
size(equilibrium_in(1)%profiles_1d%ffprime)
529 fdf = equilibrium_in(1)%profiles_1d%ffprime
532 if (
associated(equilibrium_in(1)%profiles_1d%jphi))
then
533 if (npts_set .and.
size(equilibrium_in(1)%profiles_1d%jphi) &
535 write(iu6, *)
'ERROR: input profiles have inconsistent sizes'
538 if (verbosity > 1)
write(iu6, *)
'j_tor specified on input'
540 npts =
size(equilibrium_in(1)%profiles_1d%jphi)
541 allocate(j_tor(npts))
542 j_tor = equilibrium_in(1)%profiles_1d%jphi
545 if (
associated(equilibrium_in(1)%profiles_1d%q))
then
546 if (npts_set .and.
size(equilibrium_in(1)%profiles_1d%q) &
548 write(iu6, *)
'ERROR: input profiles have inconsistent sizes'
551 if (verbosity > 1)
write(iu6, *)
'q specified on input'
553 npts =
size(equilibrium_in(1)%profiles_1d%q)
555 q_in = equilibrium_in(1)%profiles_1d%q
559 if (.not. npts_set)
then
560 write(iu6, *)
'fatal error: no input profiles specified'
564 if (
associated(equilibrium_in(1)%profiles_1d%psi))
then
565 if (npts_set .and.
size(equilibrium_in(1)%profiles_1d%psi) &
567 write(iu6, *)
'ERROR: input profiles have inconsistent sizes'
570 if (verbosity > 1)
write(iu6, *)
'psi specified on input'
572 npts =
size(equilibrium_in(1)%profiles_1d%psi)
573 allocate(psi_in(npts))
574 if (equilibrium_in(1)%profiles_1d%psi(1) /= 0._r8)
then
575 if (verbosity > 2)
write(iu6, *)
'psi shifted to be 0 on axis'
577 psi_in = (equilibrium_in(1)%profiles_1d%psi &
578 - equilibrium_in(1)%profiles_1d%psi(1))
583 if (verbosity >1)
then
584 write(iu6, *)
'no psi specified on input'
585 write(iu6, *)
'assuming equidistant psi'
587 allocate(psi_in(npts))
589 psi_in(i) = dble(i - 1) / dble(npts - 1) * cpsurfin
593 if (.not.
associated(equilibrium_in(1)%profiles_1d%psi)) &
594 allocate(equilibrium_in(1)%profiles_1d%psi(npts))
597 equilibrium_in(1)%profiles_1d%psi = psi_in
601 if (trim(input_type) ==
"p and j_tor")
then
603 call
spline(npts, psi_in, equilibrium_in(1)%profiles_1d%pressure, &
604 0._r8, 0._r8, 2, p_spline)
605 if (.not.
associated(equilibrium_in(1)%profiles_1d%pprime)) &
606 allocate(equilibrium_in(1)%profiles_1d%pprime(npts))
609 equilibrium_in(1)%profiles_1d%pprime(i) = abltg(1)
613 write(iu6, *)
"interpolated p to derive p'"
617 if (xmgrace_output)
then
618 if (
associated(equilibrium_in(1)%profiles_1d%pprime)) &
619 call
plot_data_1d(
'line', psi_in, equilibrium_in(1)%profiles_1d%pprime, npts, &
621 if (
associated(equilibrium_in(1)%profiles_1d%pressure)) &
622 call
plot_data_1d(
'line', psi_in, equilibrium_in(1)%profiles_1d%pressure, npts, &
623 'xmgr_input_pressure')
624 if (
associated(equilibrium_in(1)%profiles_1d%ffprime)) &
625 call
plot_data_1d(
'line', psi_in, equilibrium_in(1)%profiles_1d%ffprime, npts, &
626 'xmgr_input_ffprime')
627 if (
associated(equilibrium_in(1)%profiles_1d%jphi)) &
628 call
plot_data_1d(
'line', psi_in, equilibrium_in(1)%profiles_1d%jphi, npts, &
630 if (
associated(equilibrium_in(1)%profiles_1d%q)) &
631 call
plot_data_1d(
'line', psi_in, equilibrium_in(1)%profiles_1d%q, npts, &
635 write(*, *)
'pre: profiles'
639 if (.not.
allocated(fdf))
allocate(fdf(npts))
642 allocate(zvol(nrmap))
643 allocate(zvolp(nrmap))
645 allocate(psivec(npts))
646 allocate(dpres(npts))
647 allocate(p_int(npts))
649 allocate(gam_int(npts))
656 if (psi_in(1) /= 0._r8)
then
657 write(iu6, *)
'fatal error: psi on axis must be zero'
660 psivec = psi_in / psi_in(npts)
663 dpres = equilibrium_in(1)%profiles_1d%pprime
664 if (trim(input_type) ==
"p' and FF'")
then
665 dgam = equilibrium_in(1)%profiles_1d%ffprime
669 call
spline(npts, psivec, dpres, 0._r8, 0._r8, 2, dp_spline)
670 if (trim(input_type) ==
"p' and FF'")
then
671 call
spline(npts, psivec, dgam, 0._r8, 0._r8, 2, gam_spline)
673 if (
allocated(j_tor)) &
674 call
spline(npts, psivec, j_tor, 0._r8, 0._r8, 2, j_spline)
677 p_int =
integrate(dp_spline, psivec, npts, .true.) &
678 * sign(1.0_r8, cpsurfin)
679 if (trim(input_type) ==
"p' and FF'")
then
680 gam_int =
integrate(gam_spline, psivec, npts, .true.) &
681 * sign(1.0_r8, cpsurfin)
685 if (trim(input_type) ==
"p' and FF'")
then
686 a = -dgam(1) * (twopi * eps * rvac)**2 / cpsurfin
690 if (verbosity > 2)
write(iu6, *)
'initial a = ', a
695 alfa = twopi * (eps * rvac)**2 * bvac / cpsurfin
696 if (verbosity > 2)
write(iu6, *)
'initial alfa = ', alfa
699 if (xmgrace_output)
then
700 call
plot_data_1d(
'line', psivec, dpres, npts,
'xmgr_initial_dp_psipol')
702 'xmgr_initial_p_psipol')
703 if (trim(input_type) ==
"p' and FF'")
then
706 'xmgr_initial_ffprime_psipol')
707 work = 2._r8 * gam_int
708 call
plot_data_1d(
'line', psivec, work, npts,
'xmgr_initial_f2_psipol')
711 work = sqrt(abs((rvac * bvac)**2 + work)) * sign(1.0_r8, bvac)
713 'xmgr_initial_rbphi_psipol')
718 write(*, *)
'pre: equilibrium constructor'
727 allocate(kkbig(kklda, 4 * (nr + 2) * np))
729 if (trim(input_type) ==
"p' and j_tor" &
730 .or. trim(input_type) ==
"p and j_tor" &
731 .or. trim(input_type) ==
"p' and q")
then
738 if (trim(input_type) ==
"p' and FF'")
then
741 gam_int = gam_int / fscale
742 gam_spline%sp1 = gam_spline%sp1 / fscale
743 gam_spline%sp2 = gam_spline%sp2 / fscale
744 gam_spline%sp3 = gam_spline%sp3 / fscale
745 gam_spline%sp4 = gam_spline%sp4 / fscale
749 dpres = dpres / pscale
750 p_int = p_int / pscale
751 dp_spline%sp1 = dp_spline%sp1 / pscale
752 dp_spline%sp2 = dp_spline%sp2 / pscale
753 dp_spline%sp3 = dp_spline%sp3 / pscale
754 dp_spline%sp4 = dp_spline%sp4 / pscale
757 if (trim(input_type) ==
"p' and j_tor" &
758 .or. trim(input_type) ==
"p and j_tor" &
759 .or. trim(input_type) ==
"p' and q")
then
761 j_tor = j_tor / cscale
762 j_spline%sp1 = j_spline%sp1 / cscale
763 j_spline%sp2 = j_spline%sp2 / cscale
764 j_spline%sp3 = j_spline%sp3 / cscale
765 j_spline%sp4 = j_spline%sp4 / cscale
770 if (b == 0._r8 .and. trim(input_type) ==
"p' and FF'")
then
771 if (verbosity > 1)
then
772 write(iu6, *)
'no B specified on input'
773 write(iu6, *)
'automatic calculation of B'
776 b = 2._r8 * eps * mu0 * rvac**2 * pscale / fscale
778 b = mu0 * rvac**2 * pscale / fscale
782 if (trim(input_type) ==
"p' and FF'")
then
783 if (verbosity > 2)
then
784 write(iu6, *)
"================================================="
785 write(iu6, *)
"consistency check for input parameter B: "
786 write(iu6, *)
"input p'(psi) = dp/dpsi and"
787 write(iu6, *)
"F(psi)F'(psi) = 0.5 * dF^2/dpsi"
788 write(iu6, *)
"should be in SI units"
789 write(iu6, *)
"p'(0) (on input) = ", pscale
790 write(iu6, *)
"F(0)F'(0) (on input) = ", fscale
792 write(iu6, *)
"B = ", b,
" ?= ", 2._r8 * eps * mu0 * rvac**2 &
795 write(iu6, *)
"B = ", b,
" ?= ", mu0 * rvac**2 * pscale / fscale
797 write(iu6, *)
"================================================="
803 if (standard_output)
then
805 write(out_he, *)
' soloviev equilibrium '
808 write(out_he, 4) alfa
809 write(out_he, 5) radius
811 write(out_he, 7) cpsurf
813 txtout(1) =
'soloviev equilibrium'
816 if (standard_output)
then
817 write(out_he, *)
'******************************************'
818 write(out_he, *)
'* input profiles : *'
819 write(out_he, *)
'******************************************'
820 if (trim(input_type) ==
"p' and FF'" .and. .not. hbt)
then
822 '* psi, dp/dpsi, FdF/dpsi *'
823 else if (trim(input_type) ==
"p' and FF'" .and. hbt)
then
825 '* psi, dp/dpsi, gamma *'
826 else if (trim(input_type) ==
"p' and q" .and. hbt)
then
828 '* psi, dp/dpsi, gamma, q *'
829 else if (trim(input_type) ==
"p' and q" .and. .not. hbt)
then
831 '* psi, dp/dpsi, FdF/dpsi, q *'
832 else if (trim(input_type) ==
"p' and j_tor" .and. .not. hbt)
then
834 '* psi, dp/dpsi, FdF/dpsi, <j_tor> *'
835 else if (trim(input_type) ==
"p' and j_tor" .and. hbt)
then
837 '* psi, dp/dpsi, gamma, <j_tor> *'
838 else if (trim(input_type) ==
"p and j_tor" .and. .not. hbt)
then
840 '* psi, p, FdF/dpsi, <j_tor> *'
841 else if (trim(input_type) ==
"p and j_tor" .and. hbt)
then
843 '* psi, p, gamma, <j_tor> *'
846 '******************************************************'
849 write(*, *)
'pre: allocate'
853 allocate(psi(4 * nr * np))
854 allocate(xx(4, nr * np))
855 allocate(yy(4, nr * np))
856 allocate(xxold(4, nr * np))
857 allocate(yyold(4, nr * np))
858 allocate(psiold(4 * nr * np))
859 allocate(nodeno((nr - 1) * (np - 1), 4))
868 allocate(x_work(nr * np))
869 allocate(y_work(nr * np))
874 write(*, *)
'pre: initial grid'
884 if (xmgrace_output)
then
885 x_work = xx(1, 1 : nr * np) * rminor + rvac
886 y_work = yy(1, 1 : nr * np) * rminor
888 'xmgr_initial_flux_surfaces_real')
890 'xmgr_initial_radial_real')
892 'xmgr_initial_grid_x_y_real')
895 write(*, *)
'pre: node_numbers'
901 write(*, *)
'pre: flxint'
905 if (trim(input_type) ==
"p' and j_tor" &
906 .or. trim(input_type) ==
"p and j_tor" &
907 .or. trim(input_type) ==
"p' and q")
then
909 if (trim(input_type) ==
"p' and j_tor" &
910 .or. trim(input_type) ==
"p and j_tor")
then
911 call
flxint(xaxis, a, 0, fscale)
914 call
flxint2(xaxis, a, cx, cy, sumdq)
921 call
spline(npts, psivec, dgam, 0._r8, 0._r8, 2, gam_spline)
923 gam_int =
integrate(gam_spline, psivec, npts, .true.)
929 allocate(ddpsikn(nr))
932 psikn(nr - i + 1) = sg(i)**2
933 dpsikn(nr - i + 1) = 2._r8 * sg(i) * dsg(i)
934 ddpsikn(nr - i + 1) = 2._r8 * sg(i) * ddsg(i) + 2._r8 * dsg(i)**2
935 radpsi(nr - i + 1) = dble(i - 1) / dble(nr - 1)
943 equilibrium%eqgeometry%geom_axis%r = rvac
947 equilibrium%profiles_1d%psi(i) = cpsurfin * psikn(nr + 1 - i)
949 equilibrium%profiles_1d%psi(1) = 0._r8
951 write(*, *)
'pre: transfer'
955 if (
associated(equilibrium_in(1)%profiles_1d%jphi))
then
957 equilibrium_in(1)%profiles_1d%jphi, &
958 equilibrium%profiles_1d%psi, equilibrium%profiles_1d%jphi, &
961 if (
associated(equilibrium_in(1)%profiles_1d%pprime))
then
963 equilibrium_in(1)%profiles_1d%pprime, &
964 equilibrium%profiles_1d%psi, equilibrium%profiles_1d%pprime, &
967 if (
associated(equilibrium_in(1)%profiles_1d%pressure))
then
969 equilibrium_in(1)%profiles_1d%pressure, &
970 equilibrium%profiles_1d%psi, equilibrium%profiles_1d%pressure, &
978 write(*, *)
'pre: w_mhd'
981 if (verbosity > 3)
then
982 mhd_energy =
w_mhd(equilibrium, rminor, eps)
984 write(iu6, *)
'initial W_MHD = ', mhd_energy,
' J'
985 write(iu6, *)
'initial plasma volume = ', volume,
' m^3'
1007 write(*, *)
'pre: p_j'
1011 call
p_j_iteration(equilibrium_in(1), equilibrium, first, a, cx, cy, &
1012 psaxis, xaxis, yaxis, nax, rax, sax)
1089 deallocate(x_work, y_work)
1090 deallocate(psikn, dpsikn, ddpsikn, radpsi)
1092 if (verbosity > 2)
then
1093 write(iu6, *)
"================================================"
1094 write(iu6, *)
"consistency check for constrained parameter A: "
1095 write(iu6, *)
"A on input = ", a_input
1096 write(iu6, *)
"A after final iteration = ", a
1097 write(iu6, *)
"================================================"
1101 if (verbosity > 1)
write(iu6, *)
' final remesh'
1103 allocate(psikn(nrmap))
1104 allocate(dpsikn(nrmap))
1105 allocate(ddpsikn(nrmap))
1106 allocate(radpsi(nrmap))
1108 write(*, *)
'pre: remesh'
1111 call
remesh(a, nrmap, npmap, meshno, cx, cy, xaxis, yaxis, nax, rax, sax)
1112 if (verbosity > 1)
write(iu6, *)
' done final remesh'
1116 write(*, *)
'pre: fill_equilibrium'
1121 write(*, *)
'pre: w_mhd'
1124 if (verbosity > 2)
then
1125 mhd_energy =
w_mhd(equilibrium, rminor, eps)
1127 write(iu6, *)
'W_MHD = ', mhd_energy,
' J'
1128 write(iu6, *)
'plasma volume = ', volume,
' m^3'
1131 if (xmgrace_output)
then
1132 call
plot_data_1d(
'line', equilibrium%profiles_1d%volume, &
1133 equilibrium%profiles_1d%pressure, nr,
'xmgr_p_final')
1134 call
plot_data_1d(
'line', equilibrium%profiles_1d%volume, &
1135 equilibrium%profiles_1d%pprime, nr,
'xmgr_pprime_final')
1136 call
plot_data_1d(
'line', equilibrium%profiles_1d%volume, &
1137 equilibrium%profiles_1d%jphi, nr,
'xmgr_jphi_final')
1140 deallocate(xxold, yyold, psiold)
1141 deallocate(sg, dsg, ddsg)
1152 if (.not.
allocated(j_tor))
allocate(j_tor(npts))
1154 write(*, *)
'pre: rescale alfa'
1158 if (trim(match) ==
"q95")
then
1159 if (verbosity > 1)
write(iu6, *)
' safety factor', q95
1160 allocate(dummy1(npts))
1163 if (verbosity > 1)
write(iu6, *) q1out, q95, alfa
1164 alfa = alfa * q95 / q1out
1167 else if (trim(match) ==
"Ip")
then
1168 if (verbosity > 1)
write(iu6, *)
' current Ip = ', ip
1170 if (verbosity > 1)
write(iu6, *)
' calculated current = ',
current
1172 if (verbosity > 1)
write(iu6, *)
' alfa : ', alfa
1177 write(*, *)
'pre: si_currents'
1184 pskn1(i) = psikn(nr + 1 - i)
1188 allocate(xx_center((nr - 1) * (np - 1)))
1189 allocate(yy_center((nr - 1) * (np - 1)))
1197 xx_center(no) = xaxis
1198 yy_center(no) = yaxis
1202 nelm = (nr - i) * (np - 1) + j
1204 xx_center(no) = (xx(1, n1) + xx(1, n2) + xx(1, n3) + xx(1, n4)) &
1206 yy_center(no) = (yy(1, n1) + yy(1, n2) + yy(1, n3) + yy(1, n4)) &
1213 if (xmgrace_output)
then
1214 call
plot_data_1d(
'line', pskn1, j_tor_av, nr,
'xmgr_j_tor_av_psi')
1215 call
plot_data_1d(
'line', pskn1, i_tor, nr,
'xmgr_I_tor_psi')
1216 call
plot_data_1d(
'3d', xx_center, yy_center, nr - 1,
'idl_j_tor_loc', &
1220 deallocate(xx_center, yy_center)
1226 if (xmgrace_output)
then
1227 call
plot_data_1d(
'line', abs(equilibrium%profiles_1d%psi), equilibrium%profiles_1d%ffprime, nr, &
1228 'xmgr_resulting_ffprime')
1235 allocate(j_tor(npts))
1237 psivec, j_tor, spline_type = 1)
1238 call
spline(npts, psivec, j_tor, 0._r8, 0._r8, 2, j_spline)
1242 j_tor = j_tor / cscale
1243 j_spline%sp1 = j_spline%sp1 / cscale
1244 j_spline%sp2 = j_spline%sp2 / cscale
1245 j_spline%sp3 = j_spline%sp3 / cscale
1246 j_spline%sp4 = j_spline%sp4 / cscale
1248 write(*, *)
'pre: diagnostics'
1251 if (verbosity > 3)
write(txtout(13), 13) xaxis, yaxis
1255 if (verbosity > 2)
then
1257 '-----------------------------------------------------'
1259 'consistency check for Ip and computed alfa:'
1261 'current = ',
current,
', unnormalized current (from Ip) = ', &
1262 ip / 1.e+3_r8,
' kA, alfa = ', alfa
1264 'current = ',
current,
', unnormalized current (calculated) = ', &
1265 current / 1.e+3_r8,
' kA, alfa = ',alfa
1267 '-----------------------------------------------------'
1276 deallocate(p0, rbphi, dp0, drbphi)
1279 allocate(g11_hel(nr, nchi))
1280 allocate(g12_hel(nr, nchi))
1281 allocate(g33_hel(nr, nchi))
1283 allocate(chikn(nchi))
1290 allocate(drbphi(nr))
1292 write(*, *)
'pre: mapping'
1295 call
mapping(cx, cy, xaxis, yaxis, a, equilibrium_out(1))
1296 write(*, *)
'post: mapping'
1299 if (xmgrace_output)
then
1300 call
plot_data_2d(
'3d', equilibrium_out(1)%coord_sys%position%r, &
1301 equilibrium_out(1)%coord_sys%position%z,
'idl_jphi_x_y', &
1302 equilibrium_out(1)%profiles_2d(1)%jphi)
1303 call
plot_data_2d(
'3d', equilibrium_out(1)%coord_sys%position%r, &
1304 equilibrium_out(1)%coord_sys%position%z,
'idl_jpar_x_y', &
1305 equilibrium_out(1)%profiles_2d(1)%jpar)
1311 write(*, *)
'pre: w_mhd'
1314 if (verbosity > 2)
then
1315 mhd_energy =
w_mhd(equilibrium_out(1), rminor, eps)
1317 write(iu6, *)
'W_MHD = ', mhd_energy,
' J'
1318 write(iu6, *)
'plasma volume = ', volume,
' m^3'
1323 allocate(dimerc(nr))
1324 allocate(drmerc(nr))
1327 allocate(dqprof(nr))
1333 if (diagnostics_on)
then
1334 call
mercier(a, dimerc, drmerc, hh, qprof, dqprof, geonc, zjpar)
1338 deallocate(g11_hel, g12_hel, g33_hel)
1339 deallocate(chi, chikn)
1342 deallocate(p0, dp0, rbphi, drbphi)
1343 deallocate(zvol, zvolp)
1347 if (diagnostics_on)
then
1348 call
si_output(a, rvac, bvac, fcirc, qprof, dqprof, rav, geonc, &
1349 zjpar, dimerc, drmerc, hh,
current)
1352 deallocate(psikn, dpsikn, ddpsikn, radpsi)
1353 deallocate(dimerc, drmerc)
1364 if ((trim(
output) ==
'equilibrium' .or. trim(
output) ==
'full') &
1365 .and. eqdsk_file)
then
1378 deallocate(dpres, p_int, dgam, gam_int)
1384 deallocate(psi, xx, yy)
1388 if (
allocated(dp))
deallocate(dp)
1389 if (
allocated(fdf))
deallocate(fdf)
1390 if (
allocated(j_tor))
deallocate(j_tor)
1391 if (
allocated(q_in))
deallocate(q_in)
1392 if (
allocated(psi_in))
deallocate(psi_in)
1393 deallocate(s_acc, sig, weights)
1396 2
format(
'a : ', e12.4)
1397 3
format(
'b : ', e14.6)
1398 4
format(
'alfa : ', e14.6)
1399 5
format(
'radius : ', e14.6)
1400 6
format(
'b0 : ', e14.6)
1401 7
format(
'cpsurf : ', e14.6)
1402 13
format(
'magnetic axis : x = ', e14.6,
' y = ', e14.6)
1404 write(*, *)
'helena end'
1407 if (standard_output)
close(out_he)
subroutine suydam_ballooning(equilibrium_out, zvol, zvolp, xaxis)
subroutine write_eqdsk(a, curr, r0, b0, qprof, path)
subroutine helena_initialization
subroutine helena(equilibrium_in, equilibrium_out, in_path, code_parameters)
subroutine mapping(cx, cy, xaxis, yaxis, a, equilibrium_out)
subroutine, public equilibrium_constructor(equilibrium, nr, nchi)
subroutine plot_data_2d(type_plot, x, y, printname, z)
subroutine si_output(a, r0, b0, fcirc, qprof, dqprof, rav, geonc, zjpar, dimerc, drmerc, hh, current)
real(r8) function, dimension(n), public integrate(f_spline, x, n, inv)
subroutine, public equidistant_profile(f, x, f_equi, x_equi)
subroutine helena_assign_code_parameters(code_parameters, return_status)
subroutine output(NGRID, betpol, zli3)
subroutine allocate_spline_coefficients(spline, n)
subroutine flxint(xaxis, a, nmg, fscale)
subroutine passing_particles(a, fcirc, b02av, b0max, rav, equilibrium)
subroutine, public fill_equilibrium(a, xaxis, yaxis, current_averaging, equilibrium)
subroutine remesh(a, nrnew, npnew, meshno, cx, cy, xaxis, yaxis, nax, rax, sax)
subroutine spline(N, X, Y, ALFA, BETA, TYP, A, B, C, D)
subroutine deallocate_solver
REAL *8 function spwert(N, XWERT, A, B, C, D, X, ABLTG)
subroutine rft2(data, nr, kr)
subroutine plot_data_1d(type_plot, x, y, np, printname, z)
subroutine write_equidistant_profiles(equilibrium, pskn1, nr)
subroutine mercier(a, dimerc, drmerc, hh, qq, dq, geonc, zjpar)
subroutine gaussian_points
subroutine, public si_currents
subroutine current(GEOMETRY, PROFILES, TRANSPORT, SOURCES, EVOLUTION, CONTROL, j_boun, ifail, failstring)
CURRENT TRANSPORT EQUATION.
subroutine initial_grid(fr)
real(r8) function, public plasma_volume(equilibrium, rminor, eps)
subroutine deallocate_spline_coefficients(spline)
subroutine flxint2(xaxis, a, cx, cy, sumdq)
subroutine, public total_current(a, xaxis, toroidal_current)
real(r8) function, public w_mhd(equilibrium, rminor, eps)
subroutine, public equilibrium_destructor(equilibrium)
subroutine transfer_1d_profile(x_in, f_in, x_out, f_out, df_out, spline_type, monotonic)
subroutine, public current_calculations_destructor
subroutine set_node_number(n_node, n1, n2, n3, n4)
real(r8) function integrate_1d(f, x, x_low, x_high, sum)
subroutine allocate_solver(nr, np)
subroutine p_j_iteration(equil_in, equil_out, first, a, cx, cy, psaxis, xaxis, yaxis, nax, rax, sax)
subroutine write_parameters
subroutine current_at_nodes(a, equilibrium)
subroutine safety_factor(xaxis, a, cx, cy, qprf, q95out, q1out)
real(r8) function pressure(flux)
subroutine, public current_calculations_constructor(nr, np)
subroutine shape_from_points(R_bnd, Z_bnd, n_bnd, mfm, Rgeo, amin, fr)
subroutine diagnostics(a, xaxis, yaxis, zvol, zvolp, current)