PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Go to the documentation of this file.
5 #include <PB3D_macros.h>
21 integer :: X_limits(2)
22 real(dp),
allocatable :: r_F_X(:)
23 character(len=8) :: flux_name(2)
24 character(len=1) :: mode_name(2)
25 integer :: rich_lvl_name
27 logical :: debug_run_driver_X_1 = .false.
28 logical :: debug_run_driver_X_2 = .false.
29 logical :: plot_info= .false.
57 integer function run_driver_x(grid_eq,grid_eq_B,grid_X,grid_X_B,eq_1,eq_2,&
58 &X_1,X_2)
result(ierr)
74 character(*),
parameter :: rout_name =
'run_driver_X'
77 type(
grid_type),
intent(in),
target :: grid_eq
78 type(
grid_type),
intent(in),
pointer :: grid_eq_b
79 type(
grid_type),
intent(inout),
target :: grid_x
80 type(
grid_type),
intent(inout),
pointer :: grid_x_b
81 type(
eq_1_type),
intent(in),
target :: eq_1
82 type(
eq_2_type),
intent(inout),
target :: eq_2
89 integer :: eq_limits(2)
91 real(dp),
pointer :: r_f_eq(:)
92 real(dp),
pointer :: jq(:)
93 real(dp),
allocatable :: jq_ser(:)
110 ierr =
trim_grid(grid_eq,grid_eq_trim,norm_id)
113 jq => eq_1%q_saf_FD(:,0)
115 jq => eq_1%rot_t_FD(:,0)
117 ierr =
get_ser_var(jq(norm_id(1):norm_id(2)),jq_ser,scatter=.true.)
119 call grid_eq_trim%dealloc()
120 eq_limits = [grid_eq%i_min,grid_eq%i_max]
123 call writo(
'Set up perturbation normal range')
125 r_f_eq => grid_eq%r_F
127 &x_limits=x_limits,r_f_eq=r_f_eq,r_f_x=r_f_x,jq=jq_ser)
138 call writo(
'Prepare to jump to solution')
143 &grid_limits=x_limits)
149 ierr =
setup_modes(mds_x,grid_eq,grid_x,plot_name=
'X')
168 call writo(
'Skipping rest to jump to solution')
173 call writo(
'Perturbations are analyzed')
176 call writo(
'for '//trim(
i2str(
size(r_f_x)))//&
182 call writo(
'interpolation to the solution grid with '//&
195 flux_name = [
'poloidal',
'toroidal']
196 mode_name = [
'm',
'n']
198 flux_name = [
'toroidal',
'poloidal']
199 mode_name = [
'n',
'm']
203 call writo(
'with '//flux_name(2)//
' mode number '//mode_name(2)//&
207 call writo(
'and '//flux_name(1)//
' mode number '//&
212 &
' '//flux_name(1)//
' modes '//mode_name(1)//&
213 &
' that resonate most')
228 ierr =
setup_modes(mds_x,grid_eq,grid_x,plot_name=
'X')
241 call writo(
'Resonance plot not requested')
249 ierr =
run_driver_x_2(grid_eq_b,grid_x,grid_x_b,eq_1,eq_2_b,x_1,x_2)
253 call writo(
'Clean up')
256 call eq_2_b%dealloc()
265 integer function run_driver_x_0(grid_eq,grid_eq_B,grid_X,grid_X_B,eq_1,&
266 &eq_2,eq_2_B)
result(ierr)
272 character(*),
parameter :: rout_name =
'run_driver_X_0'
275 type(
grid_type),
intent(in),
target :: grid_eq
276 type(
grid_type),
intent(in),
pointer :: grid_eq_b
277 type(
grid_type),
intent(inout),
target :: grid_x
278 type(
grid_type),
intent(inout),
pointer :: grid_x_b
280 type(
eq_2_type),
intent(inout),
target :: eq_2
281 type(
eq_2_type),
intent(inout),
pointer :: eq_2_b
284 logical :: do_x_2_ops
285 integer :: rich_lvl_name
291 call writo(
'Setting up perturbation grid')
309 call writo(
'Reconstruct equilibrium variables')
320 call eq_2_b%init(grid_eq_b,setup_e=.false.,setup_f=.true.)
322 &eq_2_out=eq_2_b,eq_1=eq_1,grid_name=
'field-aligned &
330 call writo(
'Determine the perturbation grid')
332 if (
associated(grid_x%r_F))
call grid_x%dealloc()
350 if (
associated(grid_x_b))
then
351 call grid_x_b%dealloc()
366 call writo(
'Perturbation grid set up')
375 integer function run_driver_x_1(grid_eq,grid_X,eq_1,eq_2,X_1)
result(ierr)
383 character(*),
parameter :: rout_name =
'run_driver_X_1'
390 type(
x_1_type),
intent(inout) :: x_1
393 logical :: do_x_1_ops
416 if (
allocated(x_1%U_0))
call x_1%dealloc()
435 if (debug_run_driver_x_1)
then
442 call plot_info_for_vmec_hel_comparison(grid_x,x_1)
448 subroutine plot_info_for_vmec_hel_comparison(grid_X,X)
457 real(dp),
allocatable :: x_plot(:,:,:), y_plot(:,:,:), z_plot(:,:,:)
461 call writo(
'Plotting information for comparison between VMEC and &
462 &HELENA',alert=.true.)
465 call writo(
'Plots for which mode number?')
466 ld =
get_int(lim_lo=1,lim_hi=x%n_mod)
469 allocate(x_plot(grid_x%n(1),grid_x%n(2),grid_x%loc_n_r))
470 allocate(y_plot(grid_x%n(1),grid_x%n(2),grid_x%loc_n_r))
471 allocate(z_plot(grid_x%n(1),grid_x%n(2),grid_x%loc_n_r))
473 x_plot = grid_x%theta_F
475 x_plot = grid_x%zeta_F
478 do id = 1,grid_x%loc_n_r
479 z_plot(:,:,id) = grid_x%loc_r_F(id)/maxval(grid_x%r_F)
486 call plot_hdf5(
'RE_U_0',
'RE_U_0',rp(x%U_0(:,:,:,ld)),&
487 &x=x_plot,y=y_plot,z=z_plot)
488 call plot_hdf5(
'IM_U_0',
'IM_U_0',ip(x%U_0(:,:,:,ld)),&
489 &x=x_plot,y=y_plot,z=z_plot)
490 call plot_hdf5(
'RE_U_1',
'RE_U_1',rp(x%U_1(:,:,:,ld)),&
491 &x=x_plot,y=y_plot,z=z_plot)
492 call plot_hdf5(
'IM_U_1',
'IM_U_1',ip(x%U_1(:,:,:,ld)),&
493 &x=x_plot,y=y_plot,z=z_plot)
496 call plot_hdf5(
'RE_DU_0',
'RE_DU_0',rp(x%DU_0(:,:,:,ld)),&
497 &x=x_plot,y=y_plot,z=z_plot)
498 call plot_hdf5(
'IM_DU_0',
'IM_DU_0',ip(x%DU_0(:,:,:,ld)),&
499 &x=x_plot,y=y_plot,z=z_plot)
500 call plot_hdf5(
'RE_DU_1',
'RE_DU_1',rp(x%DU_1(:,:,:,ld)),&
501 &x=x_plot,y=y_plot,z=z_plot)
502 call plot_hdf5(
'IM_DU_1',
'IM_DU_1',ip(x%DU_1(:,:,:,ld)),&
503 &x=x_plot,y=y_plot,z=z_plot)
506 deallocate(x_plot,y_plot,z_plot)
508 call writo(
'Want to redo the plotting?')
512 call writo(
'Done, paused',alert=.true.)
514 end subroutine plot_info_for_vmec_hel_comparison
524 integer function run_driver_x_2(grid_eq_B,grid_X,grid_X_B,eq_1,eq_2_B,X_1,&
525 &X_2_int)
result(ierr)
537 character(*),
parameter :: rout_name =
'run_driver_X_2'
540 type(
grid_type),
intent(in),
pointer :: grid_eq_b
541 type(
grid_type),
intent(in),
target :: grid_x
542 type(
grid_type),
intent(in),
pointer :: grid_x_b
544 type(
eq_2_type),
intent(in),
pointer :: eq_2_b
546 type(
x_2_type),
intent(inout) :: x_2_int
549 logical :: no_output_loc
553 integer :: var_size_without_par
554 integer :: prev_style
555 integer :: lims_loc(2,2)
560 call writo(
'Setting up Tensorial perturbation jobs')
564 var_size_without_par = ceiling(product(grid_x_b%n(1:3))*1._dp)
565 ierr = divide_x_jobs(var_size_without_par)
572 call x_2_int%init(
mds_x,grid_x_b,is_field_averaged=.true.)
587 call writo(
'The interpolation method used for the magnetic integrals:')
591 call writo(
'magnetic interpolation style: 1 (trapezoidal rule)')
593 call writo(
'magnetic interpolation style: 2 (Simpson 3/8 rule)')
598 call writo(
'Tensorial perturbation jobs set up')
602 call writo(
'Starting tensorial perturbation jobs')
604 x_jobs:
do while(
do_x())
617 call x_1_loc(id)%init(
mds_x,grid_x,lim_sec_x=lims_loc(:,id))
620 x_1_loc(id)%U_0 = x_1%U_0(:,:,:,lims_loc(1,id):lims_loc(2,id))
621 x_1_loc(id)%U_1 = x_1%U_1(:,:,:,lims_loc(1,id):lims_loc(2,id))
622 x_1_loc(id)%DU_0 = x_1%DU_0(:,:,:,lims_loc(1,id):lims_loc(2,id))
623 x_1_loc(id)%DU_1 = x_1%DU_1(:,:,:,lims_loc(1,id):lims_loc(2,id))
631 &x_1=x_1_loc(id),grid_name=&
632 &
'field-aligned perturbation grid')
638 ierr =
calc_x(
mds_x,grid_eq_b,grid_x_b,eq_1,eq_2_b,&
639 &x_1_loc(1),x_1_loc(2),x_2,lim_sec_x=lims_loc)
662 ierr = calc_magn_ints(grid_eq_b,grid_x_b,eq_2_b,x_2,x_2_int,&
663 &prev_style,lim_sec_x=lims_loc)
668 call x_1_loc(id)%dealloc()
680 if (debug_run_driver_x_2)
then
687 call writo(
'Tensorial perturbation jobs finished')
692 &is_field_averaged=.true.)
Gather parallel variable in serial version on group master.
integer, public x_grid_style
style for normal component of X grid (1: eq, 2: sol, 3: enriched)
integer, public u_style
style for calculation of U (1: ord.2, 2: ord.1, 1: ord.0)
integer function run_driver_x_1(grid_eq, grid_X, eq_1, eq_2, X_1)
Part 1 of driver_x: Vectorial jobs.
logical, public plot_resonance
whether to plot the q-profile or iota-profile with resonances
integer function, public delete_file(file_name)
Removes a file.
Operations considering perturbation quantities.
integer, parameter, public dp
double precision
Variables that have to do with equilibrium quantities and the grid used in the calculations:
Numerical utilities related to MPI.
Numerical variables used by most other modules.
type(modes_type), public mds_x
modes variables for perturbation grid
integer, parameter, public max_str_ln
maximum length of strings
integer function, public reconstruct_pb3d_x_2(mds, grid_X, X, data_name, rich_lvl, tot_rich, lim_sec_X, lim_pos, is_field_averaged)
Reconstructs the tensorial perturbation variables from PB3D HDF5 output.
Variables concerning Richardson extrapolation.
elemental character(len=max_str_ln) function, public i2str(k)
Convert an integer to string.
real(dp), public max_par_x
max. of parallel coordinate [ ] in field-aligned grid
integer function, public print_debug_x_2(mds, grid_X, X_2_int)
Prints debug information for X_2 driver.
real(dp), public min_r_sol
min. normal range for pert.
integer function, public init_modes(grid_eq, eq)
Initializes some variables concerning the mode numbers.
integer function, public setup_grid_x(grid_eq, grid_X, r_F_X, X_limits)
Sets up the general perturbation grid, in which the perturbation variables are calculated.
integer, public eq_job_nr
nr. of eq job
integer, public x_style
style for secondary mode numbers (1: prescribed, 2: fast)
integer, public x_job_nr
nr. of X job
integer function, public calc_norm_range(style, in_limits, eq_limits, X_limits, sol_limits, r_F_eq, r_F_X, r_F_sol, jq)
Calculates normal range for the input grid, the equilibrium grid and/or the solution grid.
integer, public min_sec_x
m_X (pol. flux) or n_X (tor. flux) (only for X style 1)
integer function, public print_output_grid(grid, grid_name, data_name, rich_lvl, par_div, remove_previous_arrs)
Print grid variables to an output file.
elemental character(len=max_str_ln) function, public r2strt(k)
Convert a real (double) to string.
integer function, public run_driver_x(grid_eq, grid_eq_B, grid_X, grid_X_B, eq_1, eq_2, X_1, X_2)
Main driver of PB3D perturbation part.
integer, public n_r_sol
nr. of normal points in solution grid
real(dp), public min_par_x
min. of parallel coordinate [ ] in field-aligned grid
integer, public prim_x
n_X (pol. flux) or m_X (tor. flux)
real(dp), public max_r_sol
max. normal range for pert.
integer, public rich_restart_lvl
starting Richardson level (0: none [default])
subroutine print_info_x_2()
Prints information for tensorial perturbation job.
tensorial perturbation type
Prints variables vars with names var_names in an HDF5 file with name c file_name and accompanying XDM...
integer, public eq_style
either 1 (VMEC) or 2 (HELENA)
integer function, public setup_modes(mds, grid_eq, grid, plot_name)
Sets up some variables concerning the mode numbers.
integer, public n_par_x
nr. of parallel points in field-aligned grid
Variables pertaining to the perturbation quantities.
integer function, public interp_hel_on_grid(grid_in, grid_out, eq_2, X_1, X_2, eq_2_out, X_1_out, X_2_out, eq_1, grid_name)
Interpolate variables resulting from HELENA equilibria to another grid (angularly).
integer function, public check_x_modes(grid_eq, eq)
Checks whether the high-n approximation is valid:
integer function run_driver_x_0(grid_eq, grid_eq_B, grid_X, grid_X_B, eq_1, eq_2, eq_2_B)
part 0 of driver_x: perturbation grid as well as reconstruction of variables.
logical function, public do_x()
Tests whether this perturbatino job should be done.
Operations on HELENA variables.
Numerical utilities related to perturbation operations.
Print either vectorial or tensorial perturbation quantities of a certain order to an output file.
integer function, public divide_x_jobs(arr_size)
Divides the perturbation jobs.
logical, public no_output
no output shown
integer, public magn_int_style
style for magnetic integrals (1: trapezoidal, 2: Simpson 3/8)
subroutine, public writo(input_str, persistent, error, warning, alert)
Write output to file identified by output_i.
integer function, public reconstruct_pb3d_x_1(mds, grid_X, X, data_name, rich_lvl, tot_rich, lim_sec_X, lim_pos)
Reconstructs the vectorial perturbation variables from PB3D HDF5 output.
Numerical utilities related to giving output.
Operations on PB3D output.
integer function, public resonance_plot(mds, grid_eq, eq)
plot -profile or -profile in flux coordinates with or indicated if requested.
real(dp), parameter, public pi
Numerical utilities related to the grids and different coordinate systems.
integer function, public calc_magn_ints(grid_eq, grid_X, eq, X, X_int, prev_style, lim_sec_X)
Calculate the magnetic integrals from and in an equidistant grid where the step size can vary depen...
logical, public jump_to_sol
jump to solution
Variables pertaining to the different grids used.
integer, public max_sec_x
m_X (pol. flux) or n_X (tor. flux) (only for\ c X style 1)
integer function, public reconstruct_pb3d_grid(grid, data_name, rich_lvl, tot_rich, lim_pos, grid_limits)
Reconstructs grid variables from PB3D HDF5 output.
subroutine, public lvl_ud(inc)
Increases/decreases lvl of output.
integer function, public print_debug_x_1(mds, grid_X, X_1)
Prints debug information for X_1 driver.
integer, public n_mod_x
size of m_X (pol. flux) or n_X (tor. flux)
integer function run_driver_x_2(grid_eq_B, grid_X, grid_X_B, eq_1, eq_2_B, X_1, X_2_int)
Part 2 of driver_X: Tensorial jobs.
logical, public use_pol_flux_f
whether poloidal flux is used in F coords.
Calculates either vectorial or tensorial perturbation variables.
integer, public rich_lvl
current level of Richardson extrapolation
vectorial perturbation type
integer, dimension(:,:), allocatable, public eq_jobs_lims
data about eq jobs: [ , ] for all jobs
Operations concerning giving output, on the screen as well as in output files.
Numerical utilities related to files.
integer, dimension(:,:), allocatable, public x_jobs_lims
data about X jobs: [ , , , ] for all jobs
Driver of the perturbation part of PB3D.
integer function, public trim_grid(grid_in, grid_out, norm_id)
Trim a grid, removing any overlap between the different regions.
Operations that have to do with the grids and different coordinate systems.