ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
suydam_ballooning.f90
Go to the documentation of this file.
1 subroutine suydam_ballooning(equilibrium_out, zvol, zvolp, xaxis)
2 !-----------------------------------------------------------------------
3 ! program to determine the ballooning stability of helena equilibria
4 ! - reads the mapping file as used by castor
5 ! - calculates stability index using suydam method
6 ! - for symmetric and asymmetric plasma boundary shapes
7 !
8 ! version : 1 date : 28-09-95
9 !-----------------------------------------------------------------------
10 ! status :
11 !
12 ! 28/9/95 - tested for symmetric soloviev equilibrium with e=1.41
13 ! eps=0.382, compares well with turnbull results.
14 ! - also tested for sym. soloviev asymmetric helena
15 !
16 !-----------------------------------------------------------------------
17 
18  use itm_types
19  use mod_parameter
20  use mod_dat
21  use mod_map
22  use mod_mesh
23  use mod_profiles
24  use mod_spline
25  use mod_suydam
26  use mod_helena_io, only : out_he
27  use mod_output
28 
29  use euitm_schemas
30 
31  implicit none
32 
33  type(type_equilibrium), intent(in) :: equilibrium_out
34 
35  real(r8), dimension(nr), intent(in) :: zvol, zvolp
36  real(r8), intent(in) :: xaxis
37 
38  integer(itm_i4) :: i, j
39  real(r8) :: delf, fact
40  real(r8) :: funst, fstab, fmarg
41  real(r8) :: shear, alpha, shear1, alpha2
42 ! real :: alpha1, alpha11
43  real(r8) :: rho, drho_ds
44  real(r8) :: tbb, tbf
45  character(len = 25) :: bal,inibal
46 
47  call suydam_initialization(ias, tbb, tbf)
48 
49 !---------------------- read helena mapping file -----------------------
50  call suydam_read_helena(equilibrium_out)
51 
52 !---------------------- calculate b-field on grid points ---------------
53  call suydam_bfield
54 
55 !---------------------- calculate p and q coeff. on all gridpoints -----
56  call suydam_pq(ias)
57 
58 !---------------------- calculate stability index ----------------------
59  if (standard_output) then
60  write(out_he,*) '****************************************************', &
61  '************************************'
62  write(out_he,*) '* i, flux, rho, q, shear, shear1,', &
63  ' alpha, alpha1, fmarg, ballooning *'
64  write(out_he,*) '****************************************************', &
65  '************************************'
66  end if
67  do i = 2, npsi
68  fact = 1._r8
69  call suydam(i, 0._r8, tbb, tbf, 1._r8, bal)
70 !---------------------- find distance from stability boundary ----------
71  inibal = bal
72  if (bal == ' stable') then
73  delf = 2._r8
74  else
75  delf = 0.5_r8
76  end if
77  fact = delf * fact
78  do j = 1, 10
79  call suydam(i, 0._r8, tbb, tbf, fact, bal)
80  if (bal == inibal) then
81  fact = delf * fact
82  else
83  goto 10
84  end if
85  end do
86 !--------------------- upper and lower limits established --------------
87 10 if (inibal == ' stable') then
88  funst = fact
89  fstab = fact / 2._r8
90  else
91  fstab = fact
92  funst = fact * 2._r8
93  end if
94 !------------------------------------- bisection to find factor
95  do j = 1, 10
96  fact = (funst + fstab) / 2._r8
97  call suydam(i, 0._r8, tbb, tbf, fact, bal)
98  if (bal == ' stable') then
99  fstab = fact
100  else
101  funst = fact
102  end if
103  end do
104  shear = cs(i) / qs(i) * dqs(i)
105  alpha = -2._r8 * qs(i)**2 * p_spline%sp2(i) / eps
106 
107  shear1 = 2._r8* zvol(i) / qs(i) * dqs(i) / (2._r8 * cs(i) *zvolp(i))
108 !-------------------------------------------- lao's definition
109 ! alpha1 = - p_spline%sp2(i) / (2._r8 * cs(i)) * zvolp(i) / cpsurf**2
110 ! > * sqrt(zvol(i) / (4._r8 * pi * (1._r8 + eps *xaxis))) * eps**3.
111 !
112 !--------------------------------------------- use rho as radius
113  rho = sqrt(zvol(i) / zvol(nr))
114  drho_ds = cs(i) / (rho * zvol(nr)) * zvolp(i)
115  alpha2 = alpha / drho_ds
116 !----------------------------------------- lao corrected? needs check
117 ! alpha11 = - p_spline%sp2(i) / (2_r8 * cs(i)) * zvolp(i) / cpsurf**2
118 ! > / (pi**2 * eps) * rho * zvol(nr) * eps**4
119 !----------------------------------------------------------------------
120 
121  fmarg = (fstab + funst) / 2._r8
122 
123 !-----------------------------temporary fix for negative shear and ias=1
124 
125  if (shear < 0._r8 .and. ias == 1) then
126  fmarg = 0._r8
127  inibal = ' stable'
128  end if
129 
130  if (standard_output) &
131  write(out_he,11) i, cs(i)**2, rho, qs(i), shear, shear1, &
132  alpha, alpha2, fmarg, inibal
133  end do
134 
135  if (standard_output) &
136  write(out_he, *) '****************************************************', &
137  '*****************'
138 
139  deallocate(b02)
140  deallocate(dt_b02)
141  deallocate(ds_b02)
142  deallocate(cp0)
143  deallocate(cp1)
144  deallocate(cp2)
145  deallocate(cq0)
146  deallocate(cq1)
147 
148  call deallocate_spline_coefficients(q_spline)
149  call deallocate_spline_coefficients(p_spline)
150  call deallocate_spline_coefficients(rbphi_spline)
151 
152  11 format(' ',i3,' ',f9.6,' ',f9.6,' ',f7.4,' ',f8.4,' ',f8.4,' ', &
153  e10.4,' ',e10.4,' ',f8.3,' ',a25)
154 
155  return
156 end subroutine suydam_ballooning
subroutine suydam_ballooning(equilibrium_out, zvol, zvolp, xaxis)
subroutine suydam_pq(ias)
Definition: suydam_pq.f90:1
subroutine suydam_bfield
subroutine suydam_initialization(ias, tbb, tbf)
subroutine deallocate_spline_coefficients(spline)
subroutine suydam(ipsi, t0, tbb, tbf, fact, bal)
Definition: suydam.f90:1
subroutine suydam_read_helena(equilibrium)