PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Go to the documentation of this file.
5 #include <PB3D_macros.h>
20 logical :: plot_info = .false.
47 integer function run_driver_eq(grid_eq_out,grid_eq_B_out,eq_1_out,&
48 &eq_2_out,vac)
result(ierr)
69 character(*),
parameter :: rout_name =
'run_driver_eq'
72 type(
grid_type),
intent(inout),
target :: grid_eq_out
73 type(
grid_type),
intent(inout),
pointer :: grid_eq_b_out
74 type(
eq_1_type),
intent(inout) :: eq_1_out
75 type(
eq_2_type),
intent(inout) :: eq_2_out
83 integer :: eq_limits(2)
84 integer :: rich_lvl_name
85 logical :: do_eq_1_ops
86 logical :: do_eq_2_ops
87 logical :: only_half_grid
88 logical :: dealloc_vars
89 character(len=8) :: flux_name
90 character(len=max_str_ln) :: grid_eq_b_name
110 call writo(
'Prepare to jump to solution')
119 &
rich_lvl=rich_lvl_name,grid_limits=eq_limits)
125 &
rich_lvl=rich_lvl_name,grid_limits=eq_limits)
129 call grid_eq%dealloc()
143 call writo(
'Skipping rest to jump to solution')
148 dealloc_vars = .true.
150 if (plot_info) dealloc_vars = .false.
151 if (
ltest) dealloc_vars = .false.
156 call writo(
'The equilibrium variables are processed')
161 call writo(
'With a single field-line')
163 call writo(
'with label alpha = '//&
172 flux_name =
'poloidal'
174 flux_name =
'toroidal'
176 call writo(
'using the '//trim(flux_name)//
' flux as the normal &
183 only_half_grid = .false.
185 only_half_grid = .true.
192 do_eq_1_ops = .false.
201 do_eq_2_ops = .false.
206 call writo(
'Determine the equilibrium grid')
229 call writo(
'Flux quantities plot not requested')
242 call writo(
'Calculate angular equilibrium grid')
248 call writo(
'Determine the field-aligned equilibrium grid')
250 &only_half_grid=only_half_grid)
260 if (do_eq_2_ops)
then
267 ierr =
calc_eq(grid_eq,eq_1,eq_2,dealloc_vars=dealloc_vars)
273 call plot_info_for_vmec_hel_comparision()
291 call writo(
'Copy the equilibrium grids and variables to &
294 call writo(
'Redistribute the equilibrium grids and variables &
300 if (do_eq_2_ops)
then
301 if (
associated(grid_eq_out%r_F))
call grid_eq_out%dealloc()
304 ierr = grid_eq%copy(grid_eq_out)
313 if (do_eq_1_ops)
then
314 if (
allocated(eq_1_out%pres_FD))
call eq_1_out%dealloc()
317 call eq_1%copy(grid_eq,eq_1_out)
326 if (do_eq_2_ops)
then
327 if (
allocated(eq_2_out%jac_FD))
call eq_2_out%dealloc()
330 call eq_2%copy(grid_eq,eq_2_out)
341 grid_eq_b_out => grid_eq_out
343 if (
associated(grid_eq_b_out))
then
344 call grid_eq_b_out%dealloc()
345 deallocate(grid_eq_b_out)
346 nullify(grid_eq_b_out)
348 allocate(grid_eq_b_out)
351 ierr = grid_eq_b%copy(grid_eq_b_out)
360 call grid_eq%dealloc()
362 if (do_eq_2_ops)
then
366 call grid_eq_b%dealloc()
376 ierr =
b_plot(grid_eq_out,eq_1_out,eq_2_out,&
381 ierr =
b_plot(grid_eq_out,eq_1_out,eq_2_out)
386 if (
eq_job_nr.eq.1)
call writo(
'Magnetic field plot not requested')
394 ierr =
j_plot(grid_eq_out,eq_1_out,eq_2_out,&
399 ierr =
j_plot(grid_eq_out,eq_1_out,eq_2_out)
412 ierr =
kappa_plot(grid_eq_out,eq_1_out,eq_2_out,&
417 ierr =
kappa_plot(grid_eq_out,eq_1_out,eq_2_out)
433 grid_eq_b_name =
'eq'
435 grid_eq_b_name =
'eq_B'
452 call grid_eq_b%dealloc()
455 if (
eq_job_nr.eq.1)
call writo(
'Magnetic grid plot not requested')
461 subroutine plot_info_for_vmec_hel_comparision()
466 real(dp),
allocatable :: x_plot(:,:,:), y_plot(:,:,:), z_plot(:,:,:)
471 call writo(
'Plotting information for comparison between VMEC and &
472 &HELENA',alert=.true.)
475 call writo(
'derivative in dim 1?')
477 call writo(
'derivative in dim 2?')
479 call writo(
'derivative in dim 3?')
483 allocate(x_plot(grid_eq%n(1),grid_eq%n(2),grid_eq%loc_n_r))
484 allocate(y_plot(grid_eq%n(1),grid_eq%n(2),grid_eq%loc_n_r))
485 allocate(z_plot(grid_eq%n(1),grid_eq%n(2),grid_eq%loc_n_r))
486 x_plot = grid_eq%theta_F
488 do id = 1,grid_eq%loc_n_r
489 z_plot(:,:,id) = grid_eq%loc_r_F(id)/maxval(grid_eq%r_F)
496 call print_ex_2d(
'pres_FD',
'pres_FD',eq_1%pres_FD(:,d(2)),&
497 &x=grid_eq%loc_r_F,draw=.false.)
498 call draw_ex([
'pres_FD'],
'pres_FD',1,1,.false.)
499 call print_ex_2d(
'q_saf_FD',
'q_saf_FD',eq_1%q_saf_FD(:,d(2)),&
500 &x=grid_eq%loc_r_F,draw=.false.)
501 call draw_ex([
'q_saf_FD'],
'q_saf_FD',1,1,.false.)
503 &eq_1%flux_p_FD(:,d(2)),x=grid_eq%loc_r_F,draw=.false.)
504 call draw_ex([
'flux_p_FD'],
'flux_p_FD',1,1,.false.)
506 &eq_1%flux_t_FD(:,d(2)),x=grid_eq%loc_r_F,draw=.false.)
507 call draw_ex([
'flux_t_FD'],
'flux_t_FD',1,1,.false.)
513 &eq_2%R_E(:,:,:,d(1),d(2),d(3)),x=x_plot,y=y_plot,&
516 &eq_2%Z_E(:,:,:,d(1),d(2),d(3)),x=x_plot,y=y_plot,&
519 &eq_2%L_E(:,:,:,d(1),d(2),d(3)),x=x_plot,y=y_plot,&
522 if (sum(d).eq.0)
then
524 &reshape(
r_h(:,grid_eq%i_min:grid_eq%i_max),&
525 &[grid_eq%n(1:2),grid_eq%loc_n_r]),&
526 &x=x_plot,y=y_plot,z=z_plot)
528 &reshape(
z_h(:,grid_eq%i_min:grid_eq%i_max),&
529 &[grid_eq%n(1:2),grid_eq%loc_n_r]),&
530 &x=x_plot,y=y_plot,z=z_plot)
532 call writo(
'R and Z can only be plot for d = 0')
538 &eq_2%jac_FD(:,:,:,d(1),d(2),d(3)),&
539 &x=x_plot,y=y_plot,z=z_plot)
542 &
'g_FD_'//trim(
i2str(id)),&
543 &eq_2%g_FD(:,:,:,id,d(1),d(2),d(3)),&
544 &x=x_plot,y=y_plot,z=z_plot)
546 &
'h_FD_'//trim(
i2str(id)),&
547 &eq_2%h_FD(:,:,:,id,d(1),d(2),d(3)),&
548 &x=x_plot,y=y_plot,z=z_plot)
552 deallocate(x_plot,y_plot,z_plot)
554 call writo(
'Want to redo the plotting?')
561 call writo(
'Done, paused',alert=.true.)
563 end subroutine plot_info_for_vmec_hel_comparision
integer, public x_grid_style
style for normal component of X grid (1: eq, 2: sol, 3: enriched)
integer function, public setup_grid_eq_b(grid_eq, grid_eq_B, eq, only_half_grid)
Sets up the field-aligned equilibrium grid.
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.
integer, parameter, public max_str_ln
maximum length of strings
Variables concerning Richardson extrapolation.
elemental character(len=max_str_ln) function, public i2str(k)
Convert an integer to string.
Operations on HDF5 and XDMF variables.
real(dp), public max_par_x
max. of parallel coordinate [ ] in field-aligned grid
Calculate the equilibrium quantities on a grid determined by straight field lines.
integer, public eq_job_nr
nr. of eq job
integer function, public b_plot(grid_eq, eq_1, eq_2, rich_lvl, plot_fluxes, XYZ)
Plots the magnetic fields.
Print 2-D output on a file.
integer function, public store_vac(grid, eq_1, eq_2, vac)
Stores calculation of the vacuum response by storing the necessary variables.
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.
logical, public plot_magn_grid
whether to plot the grid in real coordinates
Redistribute the equilibrium variables, but only the Flux variables are saved.
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 j_plot(grid_eq, eq_1, eq_2, rich_lvl, plot_fluxes, XYZ)
Plots the current.
integer function, public create_output_hdf5(HDF5_name)
Creates an HDF5 output file.
logical, public ltest
whether or not to call the testing routines
Variables pertaining to the vacuum quantities.
real(dp), dimension(:), allocatable, public alpha
field line label alpha
real(dp), dimension(:,:), allocatable, public r_h
major radius (xout)
Operations and variables pertaining to the vacuum response.
integer function, public reconstruct_pb3d_eq_2(grid_eq, eq, data_name, rich_lvl, tot_rich, lim_pos)
Reconstructs the equilibrium variables from PB3D HDF5 output.
real(dp), public min_par_x
min. of parallel coordinate [ ] in field-aligned grid
integer, public rich_restart_lvl
starting Richardson level (0: none [default])
real(dp), dimension(:,:), allocatable, public z_h
height (yout)
Print equilibrium quantities to an output file:
The last line will be deleted in html
subroutine, public draw_ex(var_names, draw_name, nplt, draw_dim, plot_on_screen, ex_plot_style, data_name, draw_ops, extra_ops, is_animated, ranges, delay, persistent)
Use external program to draw a plot.
integer function, dimension(:,:), allocatable, public derivs(order, dims)
Returns derivatives of certain order.
integer function, public magn_grid_plot(grid)
Plots the grid in real 3-D space.
Prints variables vars with names var_names in an HDF5 file with name c file_name and accompanying XDM...
Driver of the equilibrium part of PB3D.
logical, public plot_b
whether to plot the magnetic field in real coordinates
integer, public eq_style
either 1 (VMEC) or 2 (HELENA)
logical, public plot_j
whether to plot the current in real coordinates
logical, public plot_kappa
whether to plot curvature
integer function, public wait_mpi()
Wait for all processes, wrapper to MPI barrier.
integer, public alpha_style
style for alpha (1: one field line, many turns, 2: many field lines, one turn)
subroutine, public writo(input_str, persistent, error, warning, alert)
Write output to file identified by output_i.
Numerical utilities related to giving output.
Variables that have to do with HELENA quantities.
Operations on PB3D output.
real(dp), parameter, public pi
integer function, public reconstruct_pb3d_eq_1(grid_eq, eq, data_name, lim_pos)
Reconstructs the equilibrium variables from PB3D HDF5 output.
integer function, public kappa_plot(grid_eq, eq_1, eq_2, rich_lvl, XYZ)
Plots the curvature.
logical, public jump_to_sol
jump to solution
Variables pertaining to the different grids used.
integer function, public calc_ang_grid_eq_b(grid_eq, eq, only_half_grid)
Calculate equilibrium grid that follows magnetic field lines.
integer function, public reconstruct_pb3d_grid(grid, data_name, rich_lvl, tot_rich, lim_pos, grid_limits)
Reconstructs grid variables from PB3D HDF5 output.
integer function, public setup_grid_eq(grid_eq, eq_limits)
Sets up the equilibrium grid.
subroutine, public lvl_ud(inc)
Increases/decreases lvl of output.
integer function, public run_driver_eq(grid_eq_out, grid_eq_B_out, eq_1_out, eq_2_out, vac)
Main driver of PB3D equilibrium part.
logical, public use_pol_flux_f
whether poloidal flux is used in F coords.
integer, public rich_lvl
current level of Richardson extrapolation
integer, dimension(:,:), allocatable, public eq_jobs_lims
data about eq jobs: [ , ] for all jobs
integer, public n_alpha
nr. of field-lines
Operations concerning giving output, on the screen as well as in output files.
logical, public plot_flux_q
whether to plot flux quantities in real coordinates
integer function, public flux_q_plot(grid_eq, eq)
Plots the flux quantities in the solution grid.
integer function, public redistribute_output_grid(grid, grid_out, no_outer_trim)
Redistribute a grid to match the normal distribution of solution grid.
Operations on the equilibrium variables.
Operations that have to do with the grids and different coordinate systems.