PB3D  [2.45]
Ideal linear high-n MHD stability in 3-D
driver_X.f90
Go to the documentation of this file.
1 !------------------------------------------------------------------------------!
3 !------------------------------------------------------------------------------!
4 module driver_x
5 #include <PB3D_macros.h>
6 #include <wrappers.h>
7  use str_utilities
8  use output_ops
9  use messages
10  use num_vars, only: dp, pi, max_str_ln
11  use grid_vars, only: grid_type
12  use eq_vars, only: eq_1_type, eq_2_type
13  use x_vars, only: x_1_type, x_2_type, modes_type, &
14  &mds_x
15 
16  implicit none
17  private
18  public run_driver_x
19 
20  ! global variables
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
26 #if ldebug
27  logical :: debug_run_driver_X_1 = .false.
28  logical :: debug_run_driver_X_2 = .false.
29  logical :: plot_info= .false.
30 #endif
31 
32 contains
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)
61  use rich_vars, only: n_par_x, rich_lvl
65  &n_mod_x
66  use grid_ops, only: calc_norm_range
70  use files_utilities, only: delete_file
71  use grid_utilities, only: trim_grid
72  use mpi_utilities, only: get_ser_var
73 
74  character(*), parameter :: rout_name = 'run_driver_X'
75 
76  ! input / output
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
83  type(x_1_type), intent(inout) :: x_1
84  type(x_2_type), intent(inout) :: x_2
85 
86  ! local variables
87  type(grid_type) :: grid_eq_trim ! trimmed equilibrium grid
88  type(eq_2_type), pointer :: eq_2_b ! field-aligned metric equilibrium variables
89  integer :: eq_limits(2) ! min. and max. index of eq. grid for this process
90  integer :: norm_id(2) ! untrimmed normal indices for trimmed grids
91  real(dp), pointer :: r_f_eq(:) ! pointer to r_F of grid_eq
92  real(dp), pointer :: jq(:) ! q_saf (pol. flux) or rot_t (tor. flux) in Flux variables
93  real(dp), allocatable :: jq_ser(:) ! serial jq
94 
95  ! initialize ierr
96  ierr = 0
97 
98  ! set up whether Richardson level has to be appended to the name
99  select case (eq_style)
100  case (1) ! VMEC
101  rich_lvl_name = rich_lvl ! append richardson level
102  case (2) ! HELENA
103  rich_lvl_name = 0 ! do not append
104  end select
105 
106  ! Divide perturbation grid under group processes, calculating the limits
107  ! and the normal coordinate.
108  if (rich_lvl.eq.rich_restart_lvl .and. eq_job_nr.eq.1) then
109  ! initialize helper variables
110  ierr = trim_grid(grid_eq,grid_eq_trim,norm_id)
111  chckerr('')
112  if (use_pol_flux_f) then
113  jq => eq_1%q_saf_FD(:,0)
114  else
115  jq => eq_1%rot_t_FD(:,0)
116  end if
117  ierr = get_ser_var(jq(norm_id(1):norm_id(2)),jq_ser,scatter=.true.)
118  chckerr('')
119  call grid_eq_trim%dealloc()
120  eq_limits = [grid_eq%i_min,grid_eq%i_max]
121 
122  ! calculate normal range
123  call writo('Set up perturbation normal range')
124  call lvl_ud(1)
125  r_f_eq => grid_eq%r_F ! needed for compiler bug on Intel 12.0.2
126  ierr = calc_norm_range('PB3D_X',eq_limits=eq_limits,&
127  &x_limits=x_limits,r_f_eq=r_f_eq,r_f_x=r_f_x,jq=jq_ser)
128  chckerr('')
129  nullify(r_f_eq)
130  call lvl_ud(-1)
131 
132  ! clean up
133  nullify(jq)
134  end if
135 
136  ! jump to solution if requested
137  if (rich_lvl.eq.rich_restart_lvl .and. jump_to_sol) then
138  call writo('Prepare to jump to solution')
139  call lvl_ud(1)
140 
141  ! grid_X
142  ierr = reconstruct_pb3d_grid(grid_x,'X',rich_lvl=rich_lvl_name,&
143  &grid_limits=x_limits)
144  chckerr('')
145 
146  ! initialize modes and set up
147  ierr = init_modes(grid_eq,eq_1)
148  chckerr('')
149  ierr = setup_modes(mds_x,grid_eq,grid_x,plot_name='X')
150  chckerr('')
151 
152  ! tests
153  ierr = check_x_modes(grid_eq,eq_1)
154  chckerr('')
155 
156  ! X_1
157  if (eq_style.eq.2) then ! HELENA
158  ierr = reconstruct_pb3d_x_1(mds_x,grid_x,x_1,'X_1')
159  chckerr('')
160  end if
161 
162  ! X_2
163  ierr = reconstruct_pb3d_x_2(mds_x,grid_x,x_2,'X_2_int',&
164  &rich_lvl=rich_lvl,is_field_averaged=.true.)
165  chckerr('')
166 
167  call lvl_ud(-1)
168  call writo('Skipping rest to jump to solution')
169  return
170  end if
171 
172  ! user output
173  call writo('Perturbations are analyzed')
174  call lvl_ud(1)
175 
176  call writo('for '//trim(i2str(size(r_f_x)))//&
177  &' values on normal range '//trim(r2strt(min_r_sol))//'..'//&
178  &trim(r2strt(max_r_sol)))
179  select case (x_grid_style)
180  case (1,3) ! equilibrium or enriched
181  call lvl_ud(1)
182  call writo('interpolation to the solution grid with '//&
183  &trim(i2str(n_r_sol))//' points will happen in the &
184  &solution driver')
185  call lvl_ud(-1)
186  case (2) ! solution
187  ! do nothing
188  end select
189  call writo('for '//trim(i2str(n_par_x))//' values on parallel &
190  &range '//trim(r2strt(min_par_x*pi))//'..'//&
191  &trim(r2strt(max_par_x*pi)))
192 
193  ! set flux and mode names
194  if (use_pol_flux_f) then
195  flux_name = ['poloidal','toroidal']
196  mode_name = ['m','n']
197  else
198  flux_name = ['toroidal','poloidal']
199  mode_name = ['n','m']
200  end if
201 
202  ! user output
203  call writo('with '//flux_name(2)//' mode number '//mode_name(2)//&
204  &' = '//trim(i2str(prim_x)))
205  select case(x_style)
206  case (1) ! prescribed
207  call writo('and '//flux_name(1)//' mode number '//&
208  &mode_name(1)//' = '//trim(i2str(min_sec_x))//'..'//&
209  &trim(i2str(max_sec_x)))
210  case (2) ! fast
211  call writo('and the '//trim(i2str(n_mod_x))//&
212  &' '//flux_name(1)//' modes '//mode_name(1)//&
213  &' that resonate most')
214  end select
215 
216  call lvl_ud(-1)
217 
218  ! Run grid part of driver
219  ierr = run_driver_x_0(grid_eq,grid_eq_b,grid_x,grid_x_b,eq_1,&
220  &eq_2,eq_2_b)
221  chckerr('')
222 
223  ! setup nm and check the X modes if restart Richardson level
224  if (rich_lvl.eq.rich_restart_lvl .and. eq_job_nr.eq.1) then
225  ! initialize modes and set up
226  ierr = init_modes(grid_eq,eq_1)
227  chckerr('')
228  ierr = setup_modes(mds_x,grid_eq,grid_x,plot_name='X')
229  chckerr('')
230 
231  ! tests
232  ierr = check_x_modes(grid_eq,eq_1)
233  chckerr('')
234  end if
235 
236  ! plot resonances if requested
237  if (plot_resonance .and. rich_lvl.eq.1 .and. eq_job_nr.eq.1) then
238  ierr = resonance_plot(mds_x,grid_eq,eq_1)
239  chckerr('')
240  else
241  call writo('Resonance plot not requested')
242  end if
243 
244  ! Run vectorial part of driver
245  ierr = run_driver_x_1(grid_eq,grid_x,eq_1,eq_2,x_1)
246  chckerr('')
247 
248  ! Run Tensorial part of driver
249  ierr = run_driver_x_2(grid_eq_b,grid_x,grid_x_b,eq_1,eq_2_b,x_1,x_2)
250  chckerr('')
251 
252  ! clean up
253  call writo('Clean up')
254  call lvl_ud(1)
255  if (eq_style.eq.2) then
256  call eq_2_b%dealloc()
257  deallocate(eq_2_b)
258  end if
259  nullify(eq_2_b)
260  call lvl_ud(-1)
261  end function run_driver_x
262 
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)
269  use rich_vars, only: rich_lvl
270  use helena_ops, only: interp_hel_on_grid
271 
272  character(*), parameter :: rout_name = 'run_driver_X_0'
273 
274  ! input / output
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
279  type(eq_1_type), intent(in) :: eq_1
280  type(eq_2_type), intent(inout), target :: eq_2
281  type(eq_2_type), intent(inout), pointer :: eq_2_b
282 
283  ! local variables
284  logical :: do_x_2_ops ! whether specific calculations for X_2 are necessary
285  integer :: rich_lvl_name ! either the Richardson level or zero, to append to names
286 
287  ! initialize ierr
288  ierr = 0
289 
290  ! user output
291  call writo('Setting up perturbation grid')
292  call lvl_ud(1)
293 
294  ! decide whether to do certain perturbation calculations and whether to
295  ! append Richardson level to name
296  select case (eq_style)
297  case (1) ! VMEC
298  do_x_2_ops = .true.
299  rich_lvl_name = rich_lvl ! append richardson level
300  case (2) ! HELENA
301  if (rich_lvl.eq.rich_restart_lvl .and. eq_job_nr.eq.1) then
302  do_x_2_ops = .true.
303  else
304  do_x_2_ops = .false.
305  end if
306  rich_lvl_name = 0 ! do not append name
307  end select
308 
309  call writo('Reconstruct equilibrium variables')
310  call lvl_ud(1)
311  select case (eq_style)
312  case (1) ! VMEC
313  ! point: grid is already field-aligned
314  eq_2_b => eq_2
315  case (2) ! HELENA
316  ! allocate field-aligned equilibrium variables
317  allocate(eq_2_b)
318 
319  ! interpolate field-aligned metric equilibrium variables
320  call eq_2_b%init(grid_eq_b,setup_e=.false.,setup_f=.true.)
321  ierr = interp_hel_on_grid(grid_eq,grid_eq_b,eq_2=eq_2,&
322  &eq_2_out=eq_2_b,eq_1=eq_1,grid_name='field-aligned &
323  &equilibrium grid')
324  chckerr('')
325  end select
326  call lvl_ud(-1)
327 
328  ! create and setup perturbation grid with division limits
329  if (do_x_2_ops) then
330  call writo('Determine the perturbation grid')
331  call lvl_ud(1)
332  if (associated(grid_x%r_F)) call grid_x%dealloc() ! deallocate if necessary
333  ierr = setup_grid_x(grid_eq,grid_x,r_f_x,x_limits)
334  chckerr('')
335  call lvl_ud(-1)
336 
337  ! print output
338  ierr = print_output_grid(grid_x,'perturbation','X',&
339  &rich_lvl=rich_lvl_name,par_div=size(eq_jobs_lims,2).gt.1)
340  chckerr('')
341  end if
342 
343  ! grid_X_B
344  select case (eq_style)
345  case (1) ! VMEC
346  ! do nothing: grid is already field-aligned
347  grid_x_b => grid_x
348  case (2) ! HELENA
349  ! set up field-aligned perturbation grid
350  if (associated(grid_x_b)) then
351  call grid_x_b%dealloc()
352  deallocate(grid_x_b)
353  nullify(grid_x_b)
354  end if
355  allocate(grid_x_b)
356  ierr = setup_grid_x(grid_eq_b,grid_x_b,r_f_x,x_limits)
357  chckerr('')
358 
359  ! print field-aligned output
360  ierr = print_output_grid(grid_x_b,'field-aligned &
361  &perturbation','X_B',rich_lvl=rich_lvl,par_div=.true.)
362  chckerr('')
363  end select
364 
365  call lvl_ud(-1)
366  call writo('Perturbation grid set up')
367  end function run_driver_x_0
368 
375  integer function run_driver_x_1(grid_eq,grid_X,eq_1,eq_2,X_1) result(ierr)
377  use x_ops, only: calc_x, print_output_x
378  use rich_vars, only: rich_lvl
379 #if ldebug
380  use x_ops, only: print_debug_x_1
381 #endif
382 
383  character(*), parameter :: rout_name = 'run_driver_X_1'
384 
385  ! input / output
386  type(grid_type), intent(in) :: grid_eq
387  type(grid_type), intent(in) :: grid_x
388  type(eq_1_type), intent(in) :: eq_1
389  type(eq_2_type), intent(in) :: eq_2
390  type(x_1_type), intent(inout) :: x_1
391 
392  ! local variables
393  logical :: do_x_1_ops ! whether specific calculations for X_1 are necessary
394 
395  ! initialize ierr
396  ierr = 0
397 
398  ! user output
399  call writo('Calculating U up to order '//trim(i2str(u_style)))
400 
401  ! decide whether to do certain perturbation calculations and whether to
402  ! append Richardson level to name
403  select case (eq_style)
404  case (1) ! VMEC
405  do_x_1_ops = .true.
406  case (2) ! HELENA
407  if (rich_lvl.eq.rich_restart_lvl .and. eq_job_nr.eq.1) then
408  do_x_1_ops = .true.
409  else
410  do_x_1_ops = .false.
411  end if
412  end select
413 
414  if (do_x_1_ops) then
415  ! possibly deallocate
416  if (allocated(x_1%U_0)) call x_1%dealloc()
417 
418  ! calc variables
419  ierr = calc_x(mds_x,grid_eq,grid_x,eq_1,eq_2,x_1)
420  chckerr('')
421  else
422  ! For HELENA, everything has been done
423  return
424  end if
425 
426  ! write vectorial perturbation variables to output if HELENA
427  ! Note: There should be only 1 parallel job
428  if (eq_style.eq.2) then ! only if first Richardson level
429  ierr = print_output_x(grid_x,x_1,'X_1',par_div=.false.)
430  chckerr('')
431  end if
432 
433 #if ldebug
434  ! write vectorial perturbation quantities to output
435  if (debug_run_driver_x_1) then
436  ierr = print_debug_x_1(mds_x,grid_x,x_1)
437  chckerr('')
438  end if
439 
440  if (plot_info) then
441  ! plot information for comparison between VMEC and HELENA
442  call plot_info_for_vmec_hel_comparison(grid_x,x_1)
443  end if
444 #endif
445 #if ldebug
446  contains
448  subroutine plot_info_for_vmec_hel_comparison(grid_X,X)
450  use num_vars, only: use_pol_flux_f
451 
452  ! input / output
453  type(grid_type), intent(in) :: grid_x
454  type(x_1_type), intent(in) :: x
455 
456  ! local variables
457  real(dp), allocatable :: x_plot(:,:,:), y_plot(:,:,:), z_plot(:,:,:)
458  integer :: id, ld
459  logical :: not_ready
460 
461  call writo('Plotting information for comparison between VMEC and &
462  &HELENA',alert=.true.)
463  not_ready = .true.
464  do while (not_ready)
465  call writo('Plots for which mode number?')
466  ld = get_int(lim_lo=1,lim_hi=x%n_mod)
467 
468  ! x, y, z
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))
472  if (use_pol_flux_f) then
473  x_plot = grid_x%theta_F
474  else
475  x_plot = grid_x%zeta_F
476  end if
477  y_plot = 0._dp
478  do id = 1,grid_x%loc_n_r
479  z_plot(:,:,id) = grid_x%loc_r_F(id)/maxval(grid_x%r_F)
480  end do
481  call plot_hdf5('x_plot','x_plot',x_plot)
482  call plot_hdf5('y_plot','y_plot',y_plot)
483  call plot_hdf5('z_plot','z_plot',z_plot)
484 
485  ! U_0 and U_1
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)
494 
495  ! DU_0 and DU_1
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)
504 
505  ! clean up
506  deallocate(x_plot,y_plot,z_plot)
507 
508  call writo('Want to redo the plotting?')
509  not_ready = get_log(.true.)
510  end do
511 
512  call writo('Done, paused',alert=.true.)
513  call pause_prog()
514  end subroutine plot_info_for_vmec_hel_comparison
515 #endif
516  end function run_driver_x_1
517 
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)
529  use rich_vars, only: rich_lvl
531  use helena_ops, only: interp_hel_on_grid
532  use x_utilities, only: do_x
533 #if ldebug
534  use x_ops, only: print_debug_x_2
535 #endif
536 
537  character(*), parameter :: rout_name = 'run_driver_X_2'
538 
539  ! input / output
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
543  type(eq_1_type), intent(in) :: eq_1
544  type(eq_2_type), intent(in), pointer :: eq_2_b
545  type(x_1_type), intent(in) :: x_1
546  type(x_2_type), intent(inout) :: x_2_int
547 
548  ! local variables
549  logical :: no_output_loc ! local copy of no_output
550  type(x_1_type) :: x_1_loc(2) ! local X_1
551  type(x_2_type) :: x_2 ! tensorial X variables
552  integer :: id ! counter
553  integer :: var_size_without_par ! size of array for jobs division
554  integer :: prev_style ! style to treat prev_X
555  integer :: lims_loc(2,2) ! local X_jobs_lims
556 
557  ! initialize ierr
558  ierr = 0
559 
560  call writo('Setting up Tensorial perturbation jobs')
561  call lvl_ud(1)
562 
563  ! divide perturbation jobs, tensor phase
564  var_size_without_par = ceiling(product(grid_x_b%n(1:3))*1._dp)
565  ierr = divide_x_jobs(var_size_without_par)
566  chckerr('')
567 
568  ! if first equilibrium job of first Richardson level, set up integrated
569  ! tensorial perturbation variables
570  if (rich_lvl.eq.rich_restart_lvl .and. eq_job_nr.eq.1) then
571  if (rich_lvl.eq.1) then
572  call x_2_int%init(mds_x,grid_x_b,is_field_averaged=.true.)
573  x_2_int%PV_0 = 0._dp
574  x_2_int%PV_1 = 0._dp
575  x_2_int%PV_2 = 0._dp
576  x_2_int%KV_0 = 0._dp
577  x_2_int%KV_1 = 0._dp
578  x_2_int%KV_2 = 0._dp
579  else
580  ierr = reconstruct_pb3d_x_2(mds_x,grid_x,x_2_int,'X_2_int',&
581  &rich_lvl=rich_lvl-1,is_field_averaged=.true.)
582  chckerr('')
583  end if
584  end if
585 
586  ! user output
587  call writo('The interpolation method used for the magnetic integrals:')
588  call lvl_ud(1)
589  select case (magn_int_style)
590  case (1)
591  call writo('magnetic interpolation style: 1 (trapezoidal rule)')
592  case (2)
593  call writo('magnetic interpolation style: 2 (Simpson 3/8 rule)')
594  end select
595  call lvl_ud(-1)
596 
597  call lvl_ud(-1)
598  call writo('Tensorial perturbation jobs set up')
599 
600  ! main loop over tensorial jobs
601  x_job_nr = 0
602  call writo('Starting tensorial perturbation jobs')
603  call lvl_ud(1)
604  x_jobs: do while(do_x())
605  ! user output
606  call print_info_x_2()
607  call lvl_ud(1)
608  no_output_loc = no_output
609  no_output = .true.
610 
611  ! set local X_jobs limits
612  lims_loc = reshape(x_jobs_lims(:,x_job_nr),[2,2])
613 
614  ! set local X_1
615  do id = 1,2
616  ! create X_1_loc
617  call x_1_loc(id)%init(mds_x,grid_x,lim_sec_x=lims_loc(:,id))
618 
619  ! fill X_1_loc
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))
624 
625  ! possibly interpolate
626  select case (eq_style)
627  case (1) ! VMEC
628  ! do nothing
629  case (2) ! HELENA
630  ierr = interp_hel_on_grid(grid_x,grid_x_b,&
631  &x_1=x_1_loc(id),grid_name=&
632  &'field-aligned perturbation grid')
633  chckerr('')
634  end select
635  end do
636 
637  ! calculate X variables, tensor phase
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)
640  chckerr('')
641 
642  ! integrate magnetic integrals of tensorial perturbation variables
643  ! over field-aligned grid. If not first Richardson level or
644  ! equilibrium job, previous integrated X_2 is used:
645  ! - rich_lvl = 1, eq_job_nr = 1: just calculate integral
646  ! - rich_lvl = 1, eq_job_nr > 1: sum to previous integral
647  ! - rich_lvl > 1, eq_job_nr = 1: sum to previous integral / 2
648  ! - rich_lvl > 1, eq_job_nr > 1: sum to previous integral
649  if (rich_lvl.eq.1) then
650  if (eq_job_nr.eq.1) then
651  prev_style = 0 ! don't add
652  else
653  prev_style = 1 ! add
654  end if
655  else
656  if (eq_job_nr.eq.1) then
657  prev_style = 2 ! divide by 2 and add, change indices
658  else
659  prev_style = 3 ! add, change indices
660  end if
661  end if
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)
664  chckerr('')
665 
666  ! clean up
667  do id = 1,2
668  call x_1_loc(id)%dealloc()
669  end do
670  call x_2%dealloc()
671 
672  ! user output
673  no_output = no_output_loc
674  call lvl_ud(-1)
675  end do x_jobs
676 
677 #if ldebug
678  ! write integrated field-aligned tensorial perturbation quantities to
679  ! output and plot
680  if (debug_run_driver_x_2) then
681  ierr = print_debug_x_2(mds_x,grid_x,x_2_int)
682  chckerr('')
683  end if
684 #endif
685 
686  call lvl_ud(-1)
687  call writo('Tensorial perturbation jobs finished')
688 
689  ! write field-averaged tensorial perturbation variables to output
690  if (eq_job_nr.eq.size(eq_jobs_lims,2)) then
691  ierr = print_output_x(grid_x_b,x_2_int,'X_2_int',rich_lvl=rich_lvl,&
692  &is_field_averaged=.true.)
693  chckerr('')
694  end if
695  end function run_driver_x_2
696 
698  subroutine print_info_x_2()
700 
701  ! user output
702  call writo('Job '//trim(i2str(x_job_nr))//'/'//&
703  &trim(i2str(size(x_jobs_lims,2)))//': mode numbers ('//&
704  &trim(i2str(x_jobs_lims(1,x_job_nr)))//'..'//&
705  &trim(i2str(x_jobs_lims(2,x_job_nr)))//')x('//&
706  &trim(i2str(x_jobs_lims(3,x_job_nr)))//'..'//&
707  &trim(i2str(x_jobs_lims(4,x_job_nr)))//')')
708  end subroutine print_info_x_2
709 end module driver_x
mpi_utilities::get_ser_var
Gather parallel variable in serial version on group master.
Definition: MPI_utilities.f90:55
num_vars::x_grid_style
integer, public x_grid_style
style for normal component of X grid (1: eq, 2: sol, 3: enriched)
Definition: num_vars.f90:99
num_vars::u_style
integer, public u_style
style for calculation of U (1: ord.2, 2: ord.1, 1: ord.0)
Definition: num_vars.f90:91
driver_x::run_driver_x_1
integer function run_driver_x_1(grid_eq, grid_X, eq_1, eq_2, X_1)
Part 1 of driver_x: Vectorial jobs.
Definition: driver_X.f90:376
num_vars::plot_resonance
logical, public plot_resonance
whether to plot the q-profile or iota-profile with resonances
Definition: num_vars.f90:102
files_utilities::delete_file
integer function, public delete_file(file_name)
Removes a file.
Definition: files_utilities.f90:144
x_ops
Operations considering perturbation quantities.
Definition: X_ops.f90:4
num_vars::dp
integer, parameter, public dp
double precision
Definition: num_vars.f90:46
eq_vars
Variables that have to do with equilibrium quantities and the grid used in the calculations:
Definition: eq_vars.f90:27
mpi_utilities
Numerical utilities related to MPI.
Definition: MPI_utilities.f90:20
num_vars
Numerical variables used by most other modules.
Definition: num_vars.f90:4
input_utilities::get_log
logical function, public get_log(yes, ind)
Queries for a logical value yes or no, where the default answer is also to be provided.
Definition: input_utilities.f90:22
x_vars::mds_x
type(modes_type), public mds_x
modes variables for perturbation grid
Definition: X_vars.f90:124
num_vars::max_str_ln
integer, parameter, public max_str_ln
maximum length of strings
Definition: num_vars.f90:50
pb3d_ops::reconstruct_pb3d_x_2
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.
Definition: PB3D_ops.f90:992
rich_vars
Variables concerning Richardson extrapolation.
Definition: rich_vars.f90:4
str_utilities::i2str
elemental character(len=max_str_ln) function, public i2str(k)
Convert an integer to string.
Definition: str_utilities.f90:18
grid_vars::max_par_x
real(dp), public max_par_x
max. of parallel coordinate [ ] in field-aligned grid
Definition: grid_vars.f90:25
x_ops::print_debug_x_2
integer function, public print_debug_x_2(mds, grid_X, X_2_int)
Prints debug information for X_2 driver.
Definition: X_ops.f90:3633
x_vars::min_r_sol
real(dp), public min_r_sol
min. normal range for pert.
Definition: X_vars.f90:135
x_ops::init_modes
integer function, public init_modes(grid_eq, eq)
Initializes some variables concerning the mode numbers.
Definition: X_ops.f90:919
grid_ops::setup_grid_x
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.
Definition: grid_ops.f90:655
num_vars::eq_job_nr
integer, public eq_job_nr
nr. of eq job
Definition: num_vars.f90:79
num_vars::x_style
integer, public x_style
style for secondary mode numbers (1: prescribed, 2: fast)
Definition: num_vars.f90:95
num_vars::x_job_nr
integer, public x_job_nr
nr. of X job
Definition: num_vars.f90:78
grid_ops::calc_norm_range
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.
Definition: grid_ops.f90:46
str_utilities
Operations on strings.
Definition: str_utilities.f90:4
x_vars::min_sec_x
integer, public min_sec_x
m_X (pol. flux) or n_X (tor. flux) (only for X style 1)
Definition: X_vars.f90:127
grid_ops::print_output_grid
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.
Definition: grid_ops.f90:1489
grid_vars::grid_type
Type for grids.
Definition: grid_vars.f90:59
str_utilities::r2strt
elemental character(len=max_str_ln) function, public r2strt(k)
Convert a real (double) to string.
Definition: str_utilities.f90:54
eq_vars::eq_1_type
flux equilibrium type
Definition: eq_vars.f90:63
driver_x::run_driver_x
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.
Definition: driver_X.f90:59
grid_vars::n_r_sol
integer, public n_r_sol
nr. of normal points in solution grid
Definition: grid_vars.f90:22
grid_vars::min_par_x
real(dp), public min_par_x
min. of parallel coordinate [ ] in field-aligned grid
Definition: grid_vars.f90:24
x_vars::prim_x
integer, public prim_x
n_X (pol. flux) or m_X (tor. flux)
Definition: X_vars.f90:126
x_vars::max_r_sol
real(dp), public max_r_sol
max. normal range for pert.
Definition: X_vars.f90:136
num_vars::rich_restart_lvl
integer, public rich_restart_lvl
starting Richardson level (0: none [default])
Definition: num_vars.f90:173
driver_x::print_info_x_2
subroutine print_info_x_2()
Prints information for tensorial perturbation job.
Definition: driver_X.f90:699
x_vars::x_2_type
tensorial perturbation type
Definition: X_vars.f90:81
output_ops::plot_hdf5
Prints variables vars with names var_names in an HDF5 file with name c file_name and accompanying XDM...
Definition: output_ops.f90:137
num_vars::eq_style
integer, public eq_style
either 1 (VMEC) or 2 (HELENA)
Definition: num_vars.f90:89
x_ops::setup_modes
integer function, public setup_modes(mds, grid_eq, grid, plot_name)
Sets up some variables concerning the mode numbers.
Definition: X_ops.f90:1060
rich_vars::n_par_x
integer, public n_par_x
nr. of parallel points in field-aligned grid
Definition: rich_vars.f90:20
x_vars
Variables pertaining to the perturbation quantities.
Definition: X_vars.f90:4
helena_ops::interp_hel_on_grid
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).
Definition: HELENA_ops.f90:550
x_ops::check_x_modes
integer function, public check_x_modes(grid_eq, eq)
Checks whether the high-n approximation is valid:
Definition: X_ops.f90:1350
driver_x::run_driver_x_0
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.
Definition: driver_X.f90:267
x_utilities::do_x
logical function, public do_x()
Tests whether this perturbatino job should be done.
Definition: X_utilities.f90:179
helena_ops
Operations on HELENA variables.
Definition: HELENA_ops.f90:4
x_utilities
Numerical utilities related to perturbation operations.
Definition: X_utilities.f90:4
x_vars::modes_type
mode number type
Definition: X_vars.f90:36
x_ops::print_output_x
Print either vectorial or tensorial perturbation quantities of a certain order to an output file.
Definition: X_ops.f90:85
x_ops::divide_x_jobs
integer function, public divide_x_jobs(arr_size)
Divides the perturbation jobs.
Definition: X_ops.f90:3364
num_vars::no_output
logical, public no_output
no output shown
Definition: num_vars.f90:145
num_vars::magn_int_style
integer, public magn_int_style
style for magnetic integrals (1: trapezoidal, 2: Simpson 3/8)
Definition: num_vars.f90:124
messages::writo
subroutine, public writo(input_str, persistent, error, warning, alert)
Write output to file identified by output_i.
Definition: messages.f90:275
pb3d_ops::reconstruct_pb3d_x_1
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.
Definition: PB3D_ops.f90:833
messages
Numerical utilities related to giving output.
Definition: messages.f90:4
pb3d_ops
Operations on PB3D output.
Definition: PB3D_ops.f90:8
x_ops::resonance_plot
integer function, public resonance_plot(mds, grid_eq, eq)
plot -profile or -profile in flux coordinates with or indicated if requested.
Definition: X_ops.f90:1814
num_vars::pi
real(dp), parameter, public pi
Definition: num_vars.f90:83
grid_utilities
Numerical utilities related to the grids and different coordinate systems.
Definition: grid_utilities.f90:4
x_ops::calc_magn_ints
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...
Definition: X_ops.f90:3081
num_vars::jump_to_sol
logical, public jump_to_sol
jump to solution
Definition: num_vars.f90:141
grid_vars
Variables pertaining to the different grids used.
Definition: grid_vars.f90:4
x_vars::max_sec_x
integer, public max_sec_x
m_X (pol. flux) or n_X (tor. flux) (only for\ c X style 1)
Definition: X_vars.f90:128
pb3d_ops::reconstruct_pb3d_grid
integer function, public reconstruct_pb3d_grid(grid, data_name, rich_lvl, tot_rich, lim_pos, grid_limits)
Reconstructs grid variables from PB3D HDF5 output.
Definition: PB3D_ops.f90:442
messages::lvl_ud
subroutine, public lvl_ud(inc)
Increases/decreases lvl of output.
Definition: messages.f90:254
input_utilities::get_int
integer function, public get_int(lim_lo, lim_hi, ind)
Queries for user input for an integer value, where allowable range can be provided as well.
Definition: input_utilities.f90:152
input_utilities::pause_prog
subroutine, public pause_prog(ind)
Pauses the running of the program.
Definition: input_utilities.f90:226
x_ops::print_debug_x_1
integer function, public print_debug_x_1(mds, grid_X, X_1)
Prints debug information for X_1 driver.
Definition: X_ops.f90:3490
x_vars::n_mod_x
integer, public n_mod_x
size of m_X (pol. flux) or n_X (tor. flux)
Definition: X_vars.f90:129
driver_x::run_driver_x_2
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.
Definition: driver_X.f90:526
num_vars::use_pol_flux_f
logical, public use_pol_flux_f
whether poloidal flux is used in F coords.
Definition: num_vars.f90:114
x_ops::calc_x
Calculates either vectorial or tensorial perturbation variables.
Definition: X_ops.f90:39
rich_vars::rich_lvl
integer, public rich_lvl
current level of Richardson extrapolation
Definition: rich_vars.f90:19
x_vars::x_1_type
vectorial perturbation type
Definition: X_vars.f90:51
num_vars::eq_jobs_lims
integer, dimension(:,:), allocatable, public eq_jobs_lims
data about eq jobs: [ , ] for all jobs
Definition: num_vars.f90:77
output_ops
Operations concerning giving output, on the screen as well as in output files.
Definition: output_ops.f90:5
input_utilities
Numerical utilities related to input.
Definition: input_utilities.f90:4
files_utilities
Numerical utilities related to files.
Definition: files_utilities.f90:4
eq_vars::eq_2_type
metric equilibrium type
Definition: eq_vars.f90:114
num_vars::x_jobs_lims
integer, dimension(:,:), allocatable, public x_jobs_lims
data about X jobs: [ , , , ] for all jobs
Definition: num_vars.f90:76
driver_x
Driver of the perturbation part of PB3D.
Definition: driver_X.f90:4
grid_utilities::trim_grid
integer function, public trim_grid(grid_in, grid_out, norm_id)
Trim a grid, removing any overlap between the different regions.
Definition: grid_utilities.f90:1636
grid_ops
Operations that have to do with the grids and different coordinate systems.
Definition: grid_ops.f90:4