ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
p_j_iteration.f90
Go to the documentation of this file.
1 subroutine p_j_iteration(equil_in, equil_out, first, a, cx, cy, psaxis, &
2  xaxis, yaxis, nax, rax, sax)
3 !-----------------------------------------------------------------------------
4 ! subroutine to solve the Grad-Shafranov equation for either p and j_tor
5 ! or dp and j_tor input
6 ! The following input combinations are allowed so far:
7 ! 1.) input_type = "p and j_tor", radial_coordinate = "psi":
8 ! HELENA will solve the GS equation for dp/dpsi(psi) and j_tor(psi) after
9 ! deriving dp/dpsi(psi) from p(psi) and psi specified in equil_in;
10 ! whether psi is maintained as in equil_in or rather Ip depends on the
11 ! specification of either alfa or Ip
12 ! 2.) input_type = "p and j_tor", radial_coordinate = "rho_vol"
13 ! HELENA will solve the GS equation for dp/dpsi(psi) and j_tor(psi) by
14 ! updating them on every iteration step. I.e. for a current solution with
15 ! psi and rho_vol, psi is transferred to the rho_vol grid of equil_in
16 ! and the input p and j_tor are interpolated and derived with respect to
17 ! this transferred psi grid and then evaluated at the current psi knots.
18 ! dp/dpsi(psi) and j_tor(psi) therefore change with every iteration step.
19 ! 3.) input_type = "p' and j_tor", radial_coordinate = "psi"
20 ! HELENA will solve the GS equation for dp/dpsi(psi) and j_tor(psi) as
21 ! specified in equil_in.
22 ! 4.) input_type = "p' and FF'", radial_coordinate = "psi"
23 ! HELENA solves the GS equation by a simple call to dp_fdf_iteration
24 ! (no current iteration needed)
25 ! Attention: still assumes that equil_in%profiles_1d%psi be set on input for
26 ! all cases (will be lifted later)
27 ! Subroutine assumes that equil_out be allocated on entry!
28 !-----------------------------------------------------------------------------
29 
30  use itm_types
31  use mod_dat
32  use mod_equilibrium
34  use mod_helena_io, only : iu6, out_he
36  use mod_mesh
37  use mod_output
38  use mod_parameter, only : mu0
39  use mod_profiles
41  use math_functions, only : integrate
42  use euitm_schemas
43  use plot_data
44  use helena_spline
45 
46  implicit none
47 
48  type(type_equilibrium), intent(inout) :: equil_in, equil_out
49  type(spline_coefficients) :: p_spline, coeff
50 
51  logical, intent(inout) :: first
52 
53  real(r8), intent(inout) :: a
54  real(r8), intent(inout) :: xaxis, yaxis
55  real(r8), intent(inout) :: cx, cy
56  real(r8), intent(inout) :: rax, sax
57 
58  real(r8), intent(out) :: psaxis
59  integer(itm_i4), intent(out) :: nax
60 
61 
62  real(r8), parameter :: fdf_eps = 1.e-7_r8
63 
64  character(len = 2) :: file_no
65  real(r8), dimension(nr) :: psi_hat
66  real(r8), dimension(size(equil_in%profiles_1d%psi)) :: rho_sqr_in
67  real(r8), dimension(size(equil_out%profiles_1d%psi)) :: rho_sqr_out
68  real(r8), dimension(npts) :: psi_tmp
69  real(r8), dimension(npts) :: fdf_old
70  real(r8), dimension(3) :: abltg
71  real(r8) :: pressure
72  real(r8) :: diff_fdf
73  real(r8) :: b_ref
74  integer(itm_i4) :: i, j
75  integer(itm_i4) :: meshno
76  logical :: converged_current
77 
78 !-- initializations
79  converged_current = .false.
80  meshno = 0
81 
82 !-- prepare and check input for "p and j_tor" cases
83  if (trim(input_type) == "p and j_tor") then
84  if (.not. associated(equil_in%profiles_1d%pressure)) then
85  write(iu6, *) 'ERROR in p_j_iteration: pressure not provided ', &
86  'on input.'
87  stop
88  end if
89  if (.not. associated(equil_in%profiles_1d%jphi)) then
90  write(iu6, *) 'ERROR in p_j_iteration: jphi not provided on input.'
91  stop
92  end if
93 
94  if (trim(radial_coordinate) == "psi") then
95 
96  if (.not. associated(equil_in%profiles_1d%psi)) then
97  write(iu6, *) 'ERROR in p_j_iteration: psi not provided ', &
98  'on input.'
99  stop
100  end if
101 
102  else if (trim(radial_coordinate) == "rho_vol") then
103 
104  if (.not. associated(equil_in%profiles_1d%rho_vol)) then
105  write(iu6, *) 'ERROR in p_j_iteration: rho_vol not provided ', &
106  'on input.'
107  stop
108  end if
109 
110 !-- avoid singular derivative at axis by interpolating rho**2
111  rho_sqr_in = equil_in%profiles_1d%rho_vol**2
112  call transfer_1d_profile(equil_in%profiles_1d%psi, &
113  rho_sqr_in, equil_out%profiles_1d%psi, &
114  rho_sqr_out, spline_type = 1, monotonic = .true.)
115  equil_out%profiles_1d%rho_vol = rho_sqr_out**0.5
116 
117  else if (trim(radial_coordinate) == "rho_tor") then
118 
119  if (.not. associated(equil_in%profiles_1d%rho_tor)) then
120  write(iu6, *) 'ERROR in p_j_iteration: rho_tor not provided ', &
121  'on input.'
122  stop
123  end if
124 
125 !-- avoid singular derivative at axis by interpolating rho**2
126  rho_sqr_in = equil_in%profiles_1d%rho_tor**2
127  call transfer_1d_profile(equil_in%profiles_1d%psi, &
128  rho_sqr_in, equil_out%profiles_1d%psi, &
129  rho_sqr_out, spline_type = 1, monotonic = .true.)
130  equil_out%profiles_1d%rho_tor = rho_sqr_out**0.5
131 
132  end if
133 
134 !-- prepare and check input for "p' and j_tor" cases
135  else if (trim(input_type) == "p' and j_tor") then
136 
137  if (.not. associated(equil_in%profiles_1d%pprime)) then
138  write(iu6, *) 'ERROR in p_j_iteration: pprime not provided ', &
139  'on input.'
140  stop
141  end if
142  if (.not. associated(equil_in%profiles_1d%jphi)) then
143  write(iu6, *) 'ERROR in p_j_iteration: jphi not provided on input.'
144  stop
145  end if
146 
147  call transfer_1d_profile(equil_in%profiles_1d%psi, &
148  equil_in%profiles_1d%pprime, equil_out%profiles_1d%psi, &
149  equil_out%profiles_1d%pprime, spline_type = 1)
150 
151  call transfer_1d_profile(equil_in%profiles_1d%psi, &
152  equil_in%profiles_1d%jphi, equil_out%profiles_1d%psi, &
153  equil_out%profiles_1d%jphi, spline_type = 1)
154 
155  end if
156 
157 !-- beginning of current loop
158  do i = 1, nmesh
159  if (verbosity > 3) then
160  write(file_no, '(i2.2)') i
161  end if
162 
163  if (verbosity > 1) then
164  write(iu6, *) 'input_type ', input_type
165  write(iu6, *) 'radial_coordinate ', radial_coordinate
166  end if
167 
168  if (trim(input_type) == "p and j_tor") then
169  if (trim(radial_coordinate) == "rho_vol") then
170 !-- interpolate new psi (equil_out) on new rho_vol, to get new psi
171 ! on reference rho_vol (equil_in)
172 !-- avoid singular derivative at axis by interpolating rho**2
173  rho_sqr_in = equil_in%profiles_1d%rho_vol**2
174  rho_sqr_out = equil_out%profiles_1d%rho_vol**2
175  call transfer_1d_profile(rho_sqr_out, &
176  equil_out%profiles_1d%psi, rho_sqr_in, &
177  psi_tmp, spline_type = 1, monotonic = .true.)
178 
179  if (verbosity > 2) then
180  write(iu6, *) 'old cpsurfin = ', &
181  equil_in%profiles_1d%psi(size(equil_in%profiles_1d%psi))
182  write(iu6, *) 'new cpsurfin = ', psi_tmp(size(psi_tmp))
183  end if
184 
185  else if (trim(radial_coordinate) == "rho_tor") then
186 !-- interpolate new psi (equil_out) on new rho_tor, to get new psi
187 ! on reference rho_tor (equil_in)
188 !-- avoid singular derivative at axis by interpolating rho**2
189  rho_sqr_in = (equil_in%profiles_1d%rho_tor &
190  / equil_in%profiles_1d%rho_tor(nr))**2
191  rho_sqr_out = (equil_out%profiles_1d%rho_tor &
192  / equil_out%profiles_1d%rho_tor(nr))**2
193  call transfer_1d_profile(rho_sqr_out, &
194  equil_out%profiles_1d%psi, rho_sqr_in, &
195  psi_tmp, spline_type = 1, monotonic = .true.)
196 
197  if (verbosity > 2) then
198  write(iu6, *) 'old cpsurfin = ', &
199  equil_in%profiles_1d%psi(size(equil_in%profiles_1d%psi))
200  write(iu6, *) 'new cpsurfin = ', psi_tmp(size(psi_tmp))
201  write(iu6, *) 'old rho_tor(nr) = ', &
202  equil_in%profiles_1d%rho_tor(size(equil_in%profiles_1d%rho_tor))
203  write(iu6, *) 'new rho_tor(nr) = ', &
204  equil_out%profiles_1d%rho_tor(size(equil_out%profiles_1d%rho_tor))
205  end if
206 
207 
208  else if (trim(radial_coordinate) == "psi") then
209 !-- maintain p(psi_bar) and j_tor(psi_bar) with psi_bar being the normalized
210 ! psi specified on input
211  psi_tmp = (equil_in%profiles_1d%psi &
212  / equil_in%profiles_1d%psi(size(equil_in%profiles_1d%psi))) &
213  * equil_out%profiles_1d%psi(size(equil_out%profiles_1d%psi))
214  end if
215 
216  if (i <= 50) then
217 !-- calculate dp/dpsi for new psi (equil_out) from reference
218 ! p given on reference rho_vol
219  call transfer_1d_profile(psi_tmp, equil_in%profiles_1d%pressure, &
220  equil_out%profiles_1d%psi, equil_out%profiles_1d%pressure, &
221  equil_out%profiles_1d%pprime, spline_type = 3)
222  else
223 !-- calculate dp/dpsi for new psi (equil_out) from reference
224 ! pprime given on reference rho_vol
225  call transfer_1d_profile(equil_in%profiles_1d%psi, &
226  equil_in%profiles_1d%pprime, equil_out%profiles_1d%psi, &
227  equil_out%profiles_1d%pprime, spline_type = 3)
228  end if
229  if (xmgrace_output .and. verbosity > 3) then
230  call plot_data_1d('line', equil_in%profiles_1d%psi, &
231  equil_in%profiles_1d%pressure, size(equil_in%profiles_1d%psi), &
232  'xmgr_p_in_' // file_no)
233  call plot_data_1d('line', equil_out%profiles_1d%psi, &
234  equil_out%profiles_1d%pressure, size(equil_out%profiles_1d%psi), &
235  'xmgr_p_out_' // file_no)
236  call plot_data_1d('line', equil_out%profiles_1d%psi, &
237  equil_out%profiles_1d%pprime, size(equil_out%profiles_1d%psi), &
238  'xmgr_pprime_out_' // file_no)
239  end if
240 
241 !-- set new reference profiles after 5 iterations
242  if (i == 50) then
243  equil_in%profiles_1d%psi = psi_tmp
244  call transfer_1d_profile(equil_out%profiles_1d%psi, &
245  equil_out%profiles_1d%pressure, psi_tmp, &
246  equil_in%profiles_1d%pressure, equil_in%profiles_1d%pprime, &
247  spline_type = 3)
248  end if
249 !-- calculate jphi for new psi (equil_out) from reference
250 ! jphi given on reference grid (rho_vol or rho_tor)
251  call transfer_1d_profile(psi_tmp, equil_in%profiles_1d%jphi, &
252  equil_out%profiles_1d%psi, equil_out%profiles_1d%jphi, spline_type = 3)
253 
254 !TODO: The following should be one contingent block (no call to
255 ! fill_equilibrium allowed in this block!
256 !-- calculate dpres on non-equidistant psivec grid with npts points
257  call deallocate_spline_coefficients(dp_spline)
258  call allocate_spline_coefficients(dp_spline, nr)
259 
260 !-- splining has to on psi_hat to be able to transfer to psivec which
261 ! goes from 0 to 1 !
262  psi_hat = equil_out%profiles_1d%psi / equil_out%profiles_1d%psi(nr)
263  call spline(nr, psi_hat, equil_out%profiles_1d%pprime, 0._r8, &
264  0._r8, 2, dp_spline)
265 
266 !TODO: Is psi_hat really different from psivec? Isn't npts=nr?
267  if (verbosity > 3) then
268  write(iu6, *) 'npts = ', npts, ' nr = ', nr
269  end if
270 !-- transfer equil_out%profiles_1d%pprime to psivec
271  do j = 1, npts
272  dpres(j) = spwert(nr, psivec(j), dp_spline, psi_hat, abltg, 0)
273  end do
274  call deallocate_spline_coefficients(dp_spline)
275 
276 !-- integrate to get new pressure
277  call allocate_spline_coefficients(dp_spline, npts)
278 !TODO: Is this extra splining necessary (not if psivec=psi_hat)?
279  call spline(npts, psivec, dpres, 0._r8, 0._r8, 2, dp_spline)
280  p_int = integrate(dp_spline, psivec, npts, .true.) &
281  * sign(1._r8, cpsurfin)
282  if (verbosity >= 2) &
283  write(iu6, *) 'new pressure on axis ', p_int(1)
284  end if
285 
286 !-- update FF' for input_type /= "p' and FF'"
287  if (trim(input_type) == "p' and j_tor" &
288  .or. trim(input_type) == "p and j_tor" &
289  .or. trim(input_type) == "p' and q") then
290 
291  if (i > 1) then
292  fdf_old = fdf
293  else
294  fdf_old = 0._r8
295  end if
296 
297 !-- calculate FF' from new p' and (old) jtor
298  call update_fdf(equil_out, xaxis)
299 !-- resulting FF' (fdf) already normalized to 1 on axis
300 
301  diff_fdf = 0._r8
302  do j = 1, npts
303  diff_fdf = diff_fdf + abs(fdf(j) - fdf_old(j)) / (fdf_eps &
304  + abs(fdf(j) + fdf_old(j)) / 2._r8)
305  end do
306 
307  if (verbosity > 1) &
308  write(iu6, *) "residual for FF' = ", diff_fdf
309 
310 !-- exit current loop if fdf converged
311  if (diff_fdf < errcur) converged_current = .true.
312 
313 !-- re-calculate gam_spline and gam_int after change of fdf (already
314 ! normalized)
315  dgam = fdf
316 !-- spline profiles
317  call spline(npts, psivec, dgam, 0._r8, 0._r8, 2, gam_spline)
318 !-- integrate profiles
319  gam_int = integrate(gam_spline, psivec, npts, .true.) &
320  * sign(1._r8, cpsurfin)
321 
322 !-- adapt B for i > 1
323  b_ref = b
324  if (verbosity > 2) &
325  write(iu6, *) 'old B = ', b_ref
326  if (hbt) then
327  b = 2._r8 * eps * equil_out%profiles_1d%pprime(1) &
328  / (equil_out%profiles_1d%pprime(1) &
329  + equil_out%profiles_1d%ffprime(1) / (mu0 * rvac**2))
330  else
331  b = mu0 * rvac**2 * equil_out%profiles_1d%pprime(1) &
332  / equil_out%profiles_1d%ffprime(1)
333  end if
334  if (verbosity > 2) &
335  write(iu6, *) 'new B = ', b
336 
337 !-- re-normalize pressure profiles (gamma profiles normalized in update_fdf)
338  pscale = dpres(1)
339  if (verbosity > 3) &
340  write(iu6, *) 'new pscale = ', pscale
341  dpres = dpres / pscale
342  p_int = p_int / pscale
343  dp_spline%sp1 = dp_spline%sp1 / pscale
344  dp_spline%sp2 = dp_spline%sp2 / pscale
345  dp_spline%sp3 = dp_spline%sp3 / pscale
346  dp_spline%sp4 = dp_spline%sp4 / pscale
347 
348  else
349 
350  converged_current = .true.
351 
352  end if
353 
354 !-- re-install old mesh for inner iteration
355  xx = xxold
356  yy = yyold
357  psi = psiold
358 
359 !-- inner iteration with fixed dp and fdf
360  call dp_fdf_iteration(first, a, psaxis, xaxis, yaxis, nax, rax, sax)
361 
362 !-- remeshing needed for calculation of flux surface quantities
363  call remesh(a, nr, np, meshno, cx, cy, xaxis, yaxis, nax, rax, sax)
364 
365 !-- adapt alfa to new a
367  call adapt_alfa(a, xaxis, cx, cy, npts)
369 
370 !-- fill equilibrium with new converged solution
371  if (trim(input_type) == "p and j_tor") then
372  call fill_equilibrium(a, xaxis, yaxis, current_averaging, equil_out)
373  if (trim(radial_coordinate) == "rho_tor") then
374 !-- calculate q and rho_tor
375  call q_calculation(xaxis, cx, cy, a, equil_out%profiles_1d%q)
376 
377  call allocate_spline_coefficients(coeff, nr)
378  call spline(nr, equil_out%profiles_1d%psi, equil_out%profiles_1d%q, &
379  0._r8, 0._r8, 3, coeff)
380  equil_out%profiles_1d%phi = integrate(coeff, &
381  equil_out%profiles_1d%psi, nr, .false.)
383 
384  equil_out%profiles_1d%rho_tor = sqrt(equil_out%profiles_1d%phi / pi &
385  / equil_out%global_param%toroid_field%b0)
386 
387  end if
388  end if
389 
390  if (verbosity > 1) &
391  write(iu6, *) 'current iteration ', i, ' completed'
392 
393 !-- calculate W_MHD
394 ! rminor = equil_out%eqgeometry%a_minor
395 ! mhd_energy = w_mhd(equil_out, rminor, eps)
396 ! volume = plasma_volume(equil_out, rminor, eps)
397 ! if (verbosity >= 2) then
398 ! write(iu6, *) 'W_MHD = ', mhd_energy, ' J'
399 ! write(iu6, *) 'plasma volume = ', volume, ' m^3'
400 ! end if
401 
402  if (verbosity > 1) &
403  write(iu6, *) "FF' converged? :", converged_current
404 
405  if (converged_current) exit
406 
407  end do
408 !-- end of loop over current
409 
410  return
411 end subroutine p_j_iteration
subroutine adapt_alfa(a, xaxis, cx, cy, npts)
subroutine dp_fdf_iteration(first, a, psaxis, xaxis, yaxis, nax, rax, sax)
subroutine, public q_calculation(xaxis, cx, cy, a, q)
real(r8) function, dimension(n), public integrate(f_spline, x, n, inv)
subroutine allocate_spline_coefficients(spline, n)
subroutine, public fill_equilibrium(a, xaxis, yaxis, current_averaging, equilibrium)
subroutine remesh(a, nrnew, npnew, meshno, cx, cy, xaxis, yaxis, nax, rax, sax)
Definition: remesh.f90:1
subroutine spline(N, X, Y, ALFA, BETA, TYP, A, B, C, D)
Definition: solution6.f90:655
REAL *8 function spwert(N, XWERT, A, B, C, D, X, ABLTG)
Definition: solution6.f90:155
subroutine update_fdf(equilibrium, xaxis)
Definition: update_fdf.f90:1
subroutine plot_data_1d(type_plot, x, y, np, printname, z)
Definition: plot_data.f90:11
subroutine deallocate_spline_coefficients(spline)
subroutine transfer_1d_profile(x_in, f_in, x_out, f_out, df_out, spline_type, monotonic)
subroutine, public current_calculations_destructor
subroutine p_j_iteration(equil_in, equil_out, first, a, cx, cy, psaxis, xaxis, yaxis, nax, rax, sax)
real(r8) function pressure(flux)
subroutine, public current_calculations_constructor(nr, np)