PB3D  [2.45]
Ideal linear high-n MHD stability in 3-D
PB3D_ops.f90
Go to the documentation of this file.
1 !------------------------------------------------------------------------------!
7 !------------------------------------------------------------------------------!
8 module pb3d_ops
9 #include <PB3D_macros.h>
10  use str_utilities
11  use messages
12  use output_ops
14  use grid_vars, only: grid_type
15  use eq_vars, only: eq_1_type, eq_2_type
16  use x_vars, only: x_1_type, x_2_type, modes_type
17  use vac_vars, only: vac_type
18  use sol_vars, only: sol_type
21 
22  implicit none
23  private
27 
28  ! global variables
29  real(dp), allocatable :: dum_1D(:) ! dummy variables
30  real(dp), allocatable :: dum_2D(:,:) ! dummy variables
31  real(dp), allocatable :: dum_3D(:,:,:) ! dummy variables
32  real(dp), allocatable :: dum_4D(:,:,:,:) ! dummy variables
33  real(dp), allocatable :: dum_6D(:,:,:,:,:,:) ! dummy variables
34  real(dp), allocatable :: dum_7D(:,:,:,:,:,:,:) ! dummy variables
35 
36 contains
40  integer function reconstruct_pb3d_in(data_name) result(ierr)
47  use hdf5_ops, only: read_hdf5_arr
48  use pb3d_utilities, only: conv_1d2nd
49  use eq_vars, only: r_0, pres_0, b_0, psi_0, rho_0, t_0, vac_perm, &
52  &n_alpha
54  &n_mod_x
55  use helena_vars, only: chi_h, flux_p_h, flux_t_h, r_h, z_h, nchi, ias, &
57  use vmec_vars, only: is_freeb_v, mnmax_v, mpol_v, ntor_v, is_asym_v, &
58  &gam_v, r_v_c, r_v_s, z_v_c, z_v_s, l_v_c, l_v_s, jac_v_c, &
60  &nfp_v, b_v_sub_c, b_v_sub_s
61  use helena_vars, only: h_h_11, h_h_12, h_h_33
62 #if ldebug
63  use vmec_vars, only: b_v_c, b_v_s, j_v_sup_int
64 #endif
65 
66  character(*), parameter :: rout_name = 'reconstruct_PB3D_in'
67 
68  ! input / output
69  character(len=*), intent(in) :: data_name
70 
71  ! local variables
72  character(len=max_str_ln) :: err_msg ! error message
73  type(var_1d_type) :: var_1d ! 1D variable
74  real(dp), parameter :: tol_version = 1.e-4_dp ! tolerance for version control
75  real(dp) :: pb3d_version ! version of PB3D variable read
76  logical :: debug_version_pb3d ! debug version of in
77  character(len=max_str_ln) :: use_debug_str(2) ! using debug or not
78 
79  ! initialize ierr
80  ierr = 0
81 
82  ! user output
83  call writo('Reconstructing input variables from PB3D output')
84  call lvl_ud(1)
85 
86  ! misc_in
87  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'misc_in')
88  chckerr('')
89  call conv_1d2nd(var_1d,dum_1d)
90  pb3d_version = dum_1d(1)
91  eq_style = nint(dum_1d(2))
92  rho_style = nint(dum_1d(3))
93  r_0 = dum_1d(4)
94  pres_0 = dum_1d(5)
95  b_0 = dum_1d(6)
96  psi_0 = dum_1d(7)
97  rho_0 = dum_1d(8)
98  t_0 = dum_1d(9)
99  vac_perm = dum_1d(10)
100  use_pol_flux_e = .false.
101  use_pol_flux_f = .false.
102  use_normalization = .false.
103  if (dum_1d(11).gt.0) use_pol_flux_e = .true.
104  if (dum_1d(12).gt.0) use_pol_flux_f = .true.
105  if (dum_1d(13).gt.0) use_normalization = .true.
106  norm_disc_prec_eq = nint(dum_1d(14))
107  n_r_in = nint(dum_1d(15))
108  n_r_eq = nint(dum_1d(16))
109  n_r_sol = nint(dum_1d(17))
110  max_flux_e = dum_1d(18)
111  max_flux_f = dum_1d(19)
112  debug_version_pb3d = .false.
113  if (dum_1d(20).gt.0) debug_version_pb3d = .true.
114  call dealloc_var_1d(var_1d)
115 
116  ! tests
117  select case (prog_style)
118  case (1) ! PB3D
119  ! do nothing
120  case (2) ! POST
121  call writo('Run tests')
122  call lvl_ud(1)
123 
124  call writo('PB3D version '//trim(r2strt(pb3d_version)))
125  if (pb3d_version.lt.min_pb3d_version*(1-tol_version)) then
126  ierr = 1
127  err_msg = 'Need at least PB3D version '//&
128  &trim(r2strt(min_pb3d_version))
129  chckerr(err_msg)
130  end if
131 
132  if (debug_version_pb3d) call writo('debug version')
133 
134  if (debug_version_pb3d.neqv.debug_version) then
135  ierr = 1
136  if (debug_version_pb3d) then
137  use_debug_str(1) = 'uses debug version'
138  else
139  use_debug_str(1) = 'uses release version'
140  end if
141  if (debug_version) then
142  use_debug_str(2) = 'uses debug version'
143  else
144  use_debug_str(2) = 'uses release version'
145  end if
146  call writo('The PB3D output '//trim(use_debug_str(1))//&
147  &' but POST '//trim(use_debug_str(2)))
148  err_msg = 'Need to use debug version for both, or not for &
149  &both'
150  chckerr(err_msg)
151  end if
152 
153  call lvl_ud(-1)
154  end select
155 
156  ! variables depending on equilibrium style
157  select case (eq_style)
158  case (1) ! VMEC
159  ! misc_in_V
160  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
161  &'misc_in_V')
162  chckerr('')
163  call conv_1d2nd(var_1d,dum_1d)
164  is_asym_v = .false.
165  is_freeb_v = .false.
166  if (dum_1d(1).gt.0) is_asym_v = .true.
167  if (dum_1d(2).gt.0) is_freeb_v = .true.
168  mnmax_v = nint(dum_1d(3))
169  mpol_v = nint(dum_1d(4))
170  ntor_v = nint(dum_1d(5))
171  nfp_v = nint(dum_1d(6))
172  gam_v = dum_1d(7)
173  call dealloc_var_1d(var_1d)
174 
175  ! flux_t_V
176  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
177  &'flux_t_V')
178  chckerr('')
179  call conv_1d2nd(var_1d,dum_2d)
180  allocate(flux_t_v(n_r_eq,0:max_deriv+1))
181  flux_t_v = dum_2d
182  call dealloc_var_1d(var_1d)
183 
184  ! flux_p_V
185  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
186  &'flux_p_V')
187  chckerr('')
188  call conv_1d2nd(var_1d,dum_2d)
189  allocate(flux_p_v(n_r_eq,0:max_deriv+1))
190  flux_p_v = dum_2d
191  call dealloc_var_1d(var_1d)
192 
193  ! pres_V
194  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'pres_V')
195  chckerr('')
196  call conv_1d2nd(var_1d,dum_2d)
197  allocate(pres_v(n_r_eq,0:max_deriv+1))
198  pres_v = dum_2d
199  call dealloc_var_1d(var_1d)
200 
201  ! rot_t_V
202  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'rot_t_V')
203  chckerr('')
204  call conv_1d2nd(var_1d,dum_2d)
205  allocate(rot_t_v(n_r_eq,0:max_deriv+1))
206  rot_t_v = dum_2d
207  call dealloc_var_1d(var_1d)
208 
209  ! q_saf_V
210  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'q_saf_V')
211  chckerr('')
212  call conv_1d2nd(var_1d,dum_2d)
213  allocate(q_saf_v(n_r_eq,0:max_deriv+1))
214  q_saf_v = dum_2d
215  call dealloc_var_1d(var_1d)
216 
217  ! mn_V
218  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'mn_V')
219  chckerr('')
220  call conv_1d2nd(var_1d,dum_2d)
221  allocate(mn_v(mnmax_v,2))
222  mn_v = nint(dum_2d)
223  call dealloc_var_1d(var_1d)
224 
225  ! RZL_V
226  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'RZL_V')
227  chckerr('')
228  call conv_1d2nd(var_1d,dum_4d)
229  allocate(r_v_c(mnmax_v,n_r_eq,0:max_deriv+1))
230  allocate(r_v_s(mnmax_v,n_r_eq,0:max_deriv+1))
231  allocate(z_v_c(mnmax_v,n_r_eq,0:max_deriv+1))
232  allocate(z_v_s(mnmax_v,n_r_eq,0:max_deriv+1))
233  allocate(l_v_c(mnmax_v,n_r_eq,0:max_deriv+1))
234  allocate(l_v_s(mnmax_v,n_r_eq,0:max_deriv+1))
235  r_v_c = dum_4d(:,:,:,1)
236  r_v_s = dum_4d(:,:,:,2)
237  z_v_c = dum_4d(:,:,:,3)
238  z_v_s = dum_4d(:,:,:,4)
239  l_v_c = dum_4d(:,:,:,5)
240  l_v_s = dum_4d(:,:,:,6)
241  call dealloc_var_1d(var_1d)
242 
243  ! jac_V
244  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'jac_V')
245  chckerr('')
246  call conv_1d2nd(var_1d,dum_4d)
247  allocate(jac_v_c(mnmax_v,n_r_eq,0:max_deriv))
248  allocate(jac_v_s(mnmax_v,n_r_eq,0:max_deriv))
249  jac_v_c = dum_4d(:,:,:,1)
250  jac_v_s = dum_4d(:,:,:,2)
251  call dealloc_var_1d(var_1d)
252 
253  ! B_V_sub
254  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'B_V_sub')
255  chckerr('')
256  call conv_1d2nd(var_1d,dum_4d)
257  allocate(b_v_sub_c(mnmax_v,n_r_eq,3))
258  allocate(b_v_sub_s(mnmax_v,n_r_eq,3))
259  b_v_sub_c = dum_4d(:,:,:,1)
260  b_v_sub_s = dum_4d(:,:,:,2)
261  call dealloc_var_1d(var_1d)
262 
263 #if ldebug
264  ! B_V
265  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'B_V')
266  chckerr('')
267  call conv_1d2nd(var_1d,dum_3d)
268  allocate(b_v_c(mnmax_v,n_r_eq))
269  allocate(b_v_s(mnmax_v,n_r_eq))
270  b_v_c = dum_3d(:,:,1)
271  b_v_s = dum_3d(:,:,2)
272  call dealloc_var_1d(var_1d)
273 
274  ! J_V_sup_int
275  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
276  &'J_V_sup_int')
277  chckerr('')
278  call conv_1d2nd(var_1d,dum_2d)
279  allocate(j_v_sup_int(n_r_eq,1:2))
280  j_v_sup_int = dum_2d
281  call dealloc_var_1d(var_1d)
282 #endif
283  case (2) ! HELENA
284  ! misc_in_H
285  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
286  &'misc_in_H')
287  chckerr('')
288  call conv_1d2nd(var_1d,dum_1d)
289  ias = nint(dum_1d(1))
290  nchi= nint(dum_1d(2))
291  call dealloc_var_1d(var_1d)
292 
293  ! RZ_H
294  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'RZ_H')
295  chckerr('')
296  call conv_1d2nd(var_1d,dum_3d)
297  allocate(r_h(nchi,n_r_eq))
298  allocate(z_h(nchi,n_r_eq))
299  r_h = dum_3d(:,:,1)
300  z_h = dum_3d(:,:,2)
301  call dealloc_var_1d(var_1d)
302 
303  ! chi_H
304  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'chi_H')
305  chckerr('')
306  call conv_1d2nd(var_1d,dum_1d)
307  allocate(chi_h(nchi))
308  chi_h = dum_1d
309  call dealloc_var_1d(var_1d)
310 
311  ! flux_p_H
312  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
313  &'flux_p_H')
314  chckerr('')
315  call conv_1d2nd(var_1d,dum_2d)
316  allocate(flux_p_h(n_r_eq,0:max_deriv+1))
317  flux_p_h = dum_2d
318  call dealloc_var_1d(var_1d)
319 
320  ! flux_t_H
321  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
322  &'flux_t_H')
323  chckerr('')
324  call conv_1d2nd(var_1d,dum_2d)
325  allocate(flux_t_h(n_r_eq,0:max_deriv+1))
326  flux_t_h = dum_2d
327  call dealloc_var_1d(var_1d)
328 
329  ! q_saf_H
330  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'q_saf_H')
331  chckerr('')
332  call conv_1d2nd(var_1d,dum_2d)
333  allocate(q_saf_h(n_r_eq,0:max_deriv+1))
334  q_saf_h = dum_2d
335  call dealloc_var_1d(var_1d)
336 
337  ! rot_t_H
338  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'rot_t_H')
339  chckerr('')
340  call conv_1d2nd(var_1d,dum_2d)
341  allocate(rot_t_h(n_r_eq,0:max_deriv+1))
342  rot_t_h = dum_2d
343  call dealloc_var_1d(var_1d)
344 
345  ! pres_H
346  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'pres_H')
347  chckerr('')
348  call conv_1d2nd(var_1d,dum_2d)
349  allocate(pres_h(n_r_eq,0:max_deriv+1))
350  pres_h = dum_2d
351  call dealloc_var_1d(var_1d)
352 
353  ! RBphi_H
354  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'RBphi_H')
355  chckerr('')
356  call conv_1d2nd(var_1d,dum_2d)
357  allocate(rbphi_h(n_r_eq,0:max_deriv+1))
358  rbphi_h = dum_2d
359  call dealloc_var_1d(var_1d)
360 
361  ! h_H
362  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'h_H')
363  chckerr('')
364  call conv_1d2nd(var_1d,dum_3d)
365  allocate(h_h_11(nchi,n_r_eq))
366  allocate(h_h_12(nchi,n_r_eq))
367  allocate(h_h_33(nchi,n_r_eq))
368  h_h_11 = dum_3d(:,:,1)
369  h_h_12 = dum_3d(:,:,2)
370  h_h_33 = dum_3d(:,:,3)
371  call dealloc_var_1d(var_1d)
372  end select
373 
374  ! misc_X
375  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'misc_X')
376  chckerr('')
377  call conv_1d2nd(var_1d,dum_1d)
378  prim_x = nint(dum_1d(1))
379  n_mod_x = nint(dum_1d(2))
380  min_sec_x = nint(dum_1d(3))
381  max_sec_x = nint(dum_1d(4))
382  norm_disc_prec_x = nint(dum_1d(5))
383  norm_style = nint(dum_1d(6))
384  u_style = nint(dum_1d(7))
385  x_style = nint(dum_1d(8))
386  matrix_slepc_style = nint(dum_1d(9))
387  select case(prog_style)
388  case (1) ! PB3D
389  magn_int_style = nint(dum_1d(10))
390  case (2) ! POST
391  magn_int_style = 1 ! integration is done in volume, currently only with trapezoidal rule
392  end select
393  k_style = nint(dum_1d(11))
394  alpha_style = nint(dum_1d(12))
395  x_grid_style = nint(dum_1d(13))
396  min_alpha = dum_1d(14)
397  max_alpha = dum_1d(15)
398  n_alpha = nint(dum_1d(16))
399  max_njq_change = dum_1d(17)
400  call dealloc_var_1d(var_1d)
401 
402  ! misc_sol
403  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'misc_sol')
404  chckerr('')
405  call conv_1d2nd(var_1d,dum_1d)
406  min_r_sol = dum_1d(1)
407  max_r_sol = dum_1d(2)
408  norm_disc_prec_sol = nint(dum_1d(3))
409  norm_disc_style_sol = nint(dum_1d(4))
410  bc_style(1) = nint(dum_1d(5))
411  bc_style(2) = nint(dum_1d(6))
412  ev_style = nint(dum_1d(7))
413  ev_bc = dum_1d(8)
414  call dealloc_var_1d(var_1d)
415 
416  ! user output
417  call lvl_ud(-1)
418  call writo('Input variables from PB3D output reconstructed')
419  end function reconstruct_pb3d_in
420 
440  integer function reconstruct_pb3d_grid(grid,data_name,rich_lvl,tot_rich,&
441  &lim_pos,grid_limits) result(ierr)
442 
443  use num_vars, only: pb3d_name
444  use hdf5_ops, only: read_hdf5_arr
445  use pb3d_utilities, only: conv_1d2nd
446 
447  character(*), parameter :: rout_name = 'reconstruct_PB3D_grid'
448 
449  ! input / output
450  type(grid_type), intent(inout) :: grid
451  character(len=*), intent(in) :: data_name
452  integer, intent(in), optional :: rich_lvl
453  logical, intent(in), optional :: tot_rich
454  integer, intent(in), optional :: lim_pos(3,2)
455  integer, intent(in), optional :: grid_limits(2)
456 
457  ! local variables
458  type(var_1d_type) :: var_1d ! 1D variable
459  integer :: lim_mem(3,2) ! memory limits for variables
460  integer :: n(3) ! total n
461  integer :: id ! counter
462  integer :: rich_lvl_loc ! local rich_lvl
463  integer :: par_id(3) ! parallel indices (start, end, stride)
464  integer :: par_id_mem(2) ! parallel indices (start, end) in memory (stride 1)
465  integer :: par_lim(2) ! parallel limits
466  integer :: rich_id(2) ! richardson level indices (start, end)
467 
468  ! initialize ierr
469  ierr = 0
470 
471  ! set up local rich_lvl
472  rich_lvl_loc = 0
473  if (present(rich_lvl)) rich_lvl_loc = rich_lvl
474 
475  ! setup rich_id
476  rich_id = setup_rich_id(rich_lvl_loc,tot_rich)
477 
478  ! get total grid size
479  ierr = get_pb3d_grid_size(n,data_name,rich_lvl,tot_rich)
480  chckerr('')
481 
482  ! possibly change n to user-specified
483  if (present(lim_pos)) then
484  if (lim_pos(1,1).ge.0 .and. lim_pos(1,2).ge.0) n(1) = &
485  &lim_pos(1,2)-lim_pos(1,1)+1
486  if (lim_pos(2,1).ge.0 .and. lim_pos(2,2).ge.0) n(2) = &
487  &lim_pos(2,2)-lim_pos(2,1)+1
488  if (lim_pos(3,1).ge.0 .and. lim_pos(3,2).ge.0) n(3) = &
489  &lim_pos(3,2)-lim_pos(3,1)+1
490  where (lim_pos(:,1).eq.lim_pos(:,2) .and. lim_pos(:,1).eq.[0,0,0]) &
491  &n = 0
492  end if
493 
494  ! set up parallel limits
495  par_lim = [1,n(1)]
496  if (present(lim_pos)) then
497  where (lim_pos(1,:).ge.0) par_lim = lim_pos(1,:)
498  end if
499 
500  ! create grid
501  ierr = grid%init(n,i_lim=grid_limits)
502  chckerr('')
503 
504  ! restore looping over richardson levels
505  do id = rich_id(2),rich_id(1),-1
506  ! setup par_id
507  par_id = setup_par_id(grid,rich_lvl_loc,id,tot_rich=tot_rich,&
508  &par_lim=par_lim,par_id_mem=par_id_mem)
509 
510  ! set up local limits for HDF5 reconstruction of full vars
511  lim_mem(1,:) = par_id_mem
512  lim_mem(2,:) = [1,grid%n(2)]
513  lim_mem(3,:) = [1,grid%n(3)]
514  if (present(lim_pos)) then
515  lim_mem(2,:) = lim_pos(2,:)
516  lim_mem(3,:) = lim_pos(3,:)
517  end if
518 
519  ! r_F
520  ierr = read_hdf5_arr(var_1d,pb3d_name,'grid_'//&
521  &trim(data_name),'r_F',rich_lvl=id,lim_loc=lim_mem(3:3,:))
522  chckerr('')
523  call conv_1d2nd(var_1d,dum_1d)
524  grid%r_F = dum_1d
525  call dealloc_var_1d(var_1d)
526 
527  ! r_E
528  ierr = read_hdf5_arr(var_1d,pb3d_name,'grid_'//&
529  &trim(data_name),'r_E',rich_lvl=id,lim_loc=lim_mem(3:3,:))
530  chckerr('')
531  call conv_1d2nd(var_1d,dum_1d)
532  grid%r_E = dum_1d
533  call dealloc_var_1d(var_1d)
534 
535  ! loc_r_F
536  grid%loc_r_F = grid%r_F(grid%i_min:grid%i_max)
537 
538  ! loc_r_E
539  grid%loc_r_E = grid%r_E(grid%i_min:grid%i_max)
540 
541  ! overwrite local limits for HDF5 reconstruction of divided vars
542  lim_mem(1,:) = par_id_mem
543  lim_mem(2,:) = [1,grid%n(2)]
544  lim_mem(3,:) = [grid%i_min,grid%i_max]
545  if (present(lim_pos)) then
546  lim_mem(2,:) = lim_pos(2,:)
547  lim_mem(3,:) = lim_mem(3,:) + lim_pos(3,1) - 1 ! take into account the grid limits (relative to position subset)
548  end if
549 
550  ! only for 3D grids
551  if (product(grid%n(1:2)).ne.0) then
552  ! theta_F
553  ierr = read_hdf5_arr(var_1d,pb3d_name,'grid_'//&
554  &trim(data_name),'theta_F',rich_lvl=id,lim_loc=lim_mem)
555  chckerr('')
556  call conv_1d2nd(var_1d,dum_3d)
557  grid%theta_F(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
558  call dealloc_var_1d(var_1d)
559 
560  ! theta_E
561  ierr = read_hdf5_arr(var_1d,pb3d_name,'grid_'//&
562  &trim(data_name),'theta_E',rich_lvl=id,lim_loc=lim_mem)
563  chckerr('')
564  call conv_1d2nd(var_1d,dum_3d)
565  grid%theta_E(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
566  call dealloc_var_1d(var_1d)
567 
568  ! zeta_F
569  ierr = read_hdf5_arr(var_1d,pb3d_name,'grid_'//&
570  &trim(data_name),'zeta_F',rich_lvl=id,lim_loc=lim_mem)
571  chckerr('')
572  call conv_1d2nd(var_1d,dum_3d)
573  grid%zeta_F(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
574  call dealloc_var_1d(var_1d)
575 
576  ! zeta_E
577  ierr = read_hdf5_arr(var_1d,pb3d_name,'grid_'//&
578  &trim(data_name),'zeta_E',rich_lvl=id,lim_loc=lim_mem)
579  chckerr('')
580  call conv_1d2nd(var_1d,dum_3d)
581  grid%zeta_E(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
582  call dealloc_var_1d(var_1d)
583  end if
584  end do
585  end function reconstruct_pb3d_grid
586 
595  integer function reconstruct_pb3d_eq_1(grid_eq,eq,data_name,lim_pos) &
596  &result(ierr) ! flux version
597  use num_vars, only: pb3d_name
598  use hdf5_ops, only: read_hdf5_arr
599  use pb3d_utilities, only: conv_1d2nd
600 
601  character(*), parameter :: rout_name = 'reconstruct_PB3D_eq_1'
602 
603  ! input / output
604  type(grid_type), intent(in) :: grid_eq
605  type(eq_1_type), intent(inout), optional :: eq
606  character(len=*), intent(in) :: data_name
607  integer, intent(in), optional :: lim_pos(1,2)
608 
609  ! local variables
610  type(var_1d_type) :: var_1d ! 1D variable
611  integer :: lim_mem(2,2) ! memory limits for variables
612 
613  ! initialize ierr
614  ierr = 0
615 
616  ! prepare
617 
618  ! create equilibrium
619  call eq%init(grid_eq,setup_e=.false.,setup_f=.true.)
620 
621  ! set up local limits for HDF5 reconstruction
622  lim_mem(1,:) = [grid_eq%i_min,grid_eq%i_max]
623  lim_mem(2,:) = [0,-1] ! all derivatives, starting from 0
624  if (present(lim_pos)) then
625  lim_mem(1,:) = lim_mem(1,:) + lim_pos(1,1) - 1 ! take into account the grid limits (relative to position subset)
626  end if
627 
628  ! restore variables
629 
630  ! pres_FD
631  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'pres_FD',&
632  &lim_loc=lim_mem)
633  chckerr('')
634  call conv_1d2nd(var_1d,dum_2d)
635  eq%pres_FD = dum_2d
636  call dealloc_var_1d(var_1d)
637 
638  ! q_saf_FD
639  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'q_saf_FD',&
640  &lim_loc=lim_mem)
641  chckerr('')
642  call conv_1d2nd(var_1d,dum_2d)
643  eq%q_saf_FD = dum_2d
644  call dealloc_var_1d(var_1d)
645 
646  ! rot_t_FD
647  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'rot_t_FD',&
648  &lim_loc=lim_mem)
649  chckerr('')
650  call conv_1d2nd(var_1d,dum_2d)
651  eq%rot_t_FD = dum_2d
652  call dealloc_var_1d(var_1d)
653 
654  ! flux_p_FD
655  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'flux_p_FD',&
656  &lim_loc=lim_mem)
657  chckerr('')
658  call conv_1d2nd(var_1d,dum_2d)
659  eq%flux_p_FD = dum_2d
660  call dealloc_var_1d(var_1d)
661 
662  ! flux_t_FD
663  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'flux_t_FD',&
664  &lim_loc=lim_mem)
665  chckerr('')
666  call conv_1d2nd(var_1d,dum_2d)
667  eq%flux_t_FD = dum_2d
668  call dealloc_var_1d(var_1d)
669 
670  ! rho
671  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'rho',&
672  &lim_loc=lim_mem(1:1,:))
673  chckerr('')
674  call conv_1d2nd(var_1d,dum_1d)
675  eq%rho = dum_1d
676  call dealloc_var_1d(var_1d)
677  end function reconstruct_pb3d_eq_1
678 
693  integer function reconstruct_pb3d_eq_2(grid_eq,eq,data_name,rich_lvl,&
694  &tot_rich,lim_pos) result(ierr) ! metric version
695  use num_vars, only: pb3d_name
696  use hdf5_ops, only: read_hdf5_arr
697  use pb3d_utilities, only: conv_1d2nd
698 
699  character(*), parameter :: rout_name = 'reconstruct_PB3D_eq_2'
700 
701  ! input / output
702  type(grid_type), intent(in) :: grid_eq
703  type(eq_2_type), intent(inout) :: eq
704  character(len=*), intent(in) :: data_name
705  integer, intent(in), optional :: rich_lvl
706  logical, intent(in), optional :: tot_rich
707  integer, intent(in), optional :: lim_pos(3,2)
708 
709  ! local variables
710  type(var_1d_type) :: var_1d ! 1D variable
711  integer :: lim_mem(7,2) ! memory limits for variables
712  integer :: id ! counter
713  integer :: rich_lvl_loc ! local rich_lvl
714  integer :: par_id(3) ! parallel indices (start, end, stride)
715  integer :: par_id_mem(2) ! parallel indices (start, end) in memory (stride 1)
716  integer :: par_lim(2) ! parallel limits
717  integer :: rich_id(2) ! richardson level indices (start, end)
718 
719  ! initialize ierr
720  ierr = 0
721 
722  ! set up local rich_lvl
723  rich_lvl_loc = 0
724  if (present(rich_lvl)) rich_lvl_loc = rich_lvl
725 
726  ! setup rich_id
727  rich_id = setup_rich_id(rich_lvl_loc,tot_rich)
728 
729  ! set up parallel limits
730  par_lim = [1,grid_eq%n(1)]
731  if (present(lim_pos)) then
732  where (lim_pos(1,:).ge.0) par_lim = lim_pos(1,:)
733  end if
734 
735  ! create equilibrium
736  call eq%init(grid_eq,setup_e=.false.,setup_f=.true.)
737 
738  ! restore looping over richardson levels
739  do id = rich_id(2),rich_id(1),-1
740  ! setup par_id
741  par_id = setup_par_id(grid_eq,rich_lvl_loc,id,tot_rich=tot_rich,&
742  &par_lim=par_lim,par_id_mem=par_id_mem)
743 
744  ! set up local limits for HDF5 reconstruction
745  lim_mem(1,:) = par_id_mem
746  lim_mem(2,:) = [1,grid_eq%n(2)]
747  lim_mem(3,:) = [grid_eq%i_min,grid_eq%i_max]
748  lim_mem(4,:) = [-1,-1]
749  lim_mem(5,:) = [-1,-1]
750  lim_mem(6,:) = [-1,-1]
751  lim_mem(7,:) = [-1,-1]
752  if (present(lim_pos)) then
753  lim_mem(2,:) = lim_pos(2,:)
754  lim_mem(3,:) = lim_mem(3,:) + lim_pos(3,1) - 1 ! take into account the grid limits (relative to position subset)
755  end if
756 
757  ! g_FD
758  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'g_FD',&
759  &rich_lvl=id,lim_loc=lim_mem)
760  chckerr('')
761  call conv_1d2nd(var_1d,dum_7d)
762  eq%g_FD(par_id(1):par_id(2):par_id(3),:,:,:,:,:,:) = dum_7d
763  call dealloc_var_1d(var_1d)
764 
765  ! h_FD
766  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'h_FD',&
767  &rich_lvl=id,lim_loc=lim_mem)
768  chckerr('')
769  call conv_1d2nd(var_1d,dum_7d)
770  eq%h_FD(par_id(1):par_id(2):par_id(3),:,:,:,:,:,:) = dum_7d
771  call dealloc_var_1d(var_1d)
772 
773  ! jac_FD
774  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'jac_FD',&
775  &rich_lvl=id,lim_loc=lim_mem(1:6,:))
776  chckerr('')
777  call conv_1d2nd(var_1d,dum_6d)
778  eq%jac_FD(par_id(1):par_id(2):par_id(3),:,:,:,:,:) = dum_6d
779  call dealloc_var_1d(var_1d)
780 
781  ! S
782  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'S',&
783  &rich_lvl=id,lim_loc=lim_mem(1:3,:))
784  chckerr('')
785  call conv_1d2nd(var_1d,dum_3d)
786  eq%S(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
787  call dealloc_var_1d(var_1d)
788 
789  ! kappa_n
790  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
791  &'kappa_n',rich_lvl=id,lim_loc=lim_mem(1:3,:))
792  chckerr('')
793  call conv_1d2nd(var_1d,dum_3d)
794  eq%kappa_n(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
795  call dealloc_var_1d(var_1d)
796 
797  ! kappa_g
798  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
799  &'kappa_g',rich_lvl=id,lim_loc=lim_mem(1:3,:))
800  chckerr('')
801  call conv_1d2nd(var_1d,dum_3d)
802  eq%kappa_g(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
803  call dealloc_var_1d(var_1d)
804 
805  ! sigma
806  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'sigma',&
807  &rich_lvl=id,lim_loc=lim_mem(1:3,:))
808  chckerr('')
809  call conv_1d2nd(var_1d,dum_3d)
810  eq%sigma(par_id(1):par_id(2):par_id(3),:,:) = dum_3d
811  call dealloc_var_1d(var_1d)
812  end do
813  end function reconstruct_pb3d_eq_2
814 
831  integer function reconstruct_pb3d_x_1(mds,grid_X,X,data_name,rich_lvl,&
832  &tot_rich,lim_sec_X,lim_pos) result(ierr)
833  use num_vars, only: pb3d_name
834  use x_vars, only: n_mod_x
835  use hdf5_ops, only: read_hdf5_arr
836  use pb3d_utilities, only: conv_1d2nd
837 
838  character(*), parameter :: rout_name = 'reconstruct_PB3D_X_1'
839 
840  ! input / output
841  type(modes_type), intent(in) :: mds
842  type(grid_type), intent(in) :: grid_x
843  type(x_1_type), intent(inout) :: x
844  character(len=*), intent(in) :: data_name
845  integer, intent(in), optional :: rich_lvl
846  logical, intent(in), optional :: tot_rich
847  integer, intent(in), optional :: lim_sec_x(2)
848  integer, intent(in), optional :: lim_pos(3,2)
849 
850  ! local variables
851  type(var_1d_type) :: var_1d ! 1D variable
852  integer :: lim_sec_x_loc(2) ! local version of lim_sec_X
853  integer :: lim_mem(4,2) ! memory limits for variables
854  integer :: id ! counter
855  integer :: rich_lvl_loc ! local rich_lvl
856  integer :: par_id(3) ! parallel indices (start, end, stride)
857  integer :: par_id_mem(2) ! parallel indices (start, end) in memory (stride 1)
858  integer :: par_lim(2) ! parallel limits
859  integer :: rich_id(2) ! richardson level indices (start, end)
860 
861  ! initialize ierr
862  ierr = 0
863 
864  ! set up local rich_lvl
865  rich_lvl_loc = 0
866  if (present(rich_lvl)) rich_lvl_loc = rich_lvl
867 
868  ! setup rich_id
869  rich_id = setup_rich_id(rich_lvl_loc,tot_rich)
870 
871  ! set up local lim_sec_X
872  lim_sec_x_loc = [1,n_mod_x]
873  if (present(lim_sec_x)) lim_sec_x_loc = lim_sec_x
874 
875  ! set up parallel limits
876  par_lim = [1,grid_x%n(1)]
877  if (present(lim_pos)) then
878  where (lim_pos(1,:).ge.0) par_lim = lim_pos(1,:)
879  end if
880 
881  ! create X
882  call x%init(mds,grid_x,lim_sec_x)
883 
884  ! restore looping over richardson levels
885  do id = rich_id(2),rich_id(1),-1
886  ! setup par_id
887  par_id = setup_par_id(grid_x,rich_lvl_loc,id,tot_rich=tot_rich,&
888  &par_lim=par_lim,par_id_mem=par_id_mem)
889 
890  ! set up local limits for HDF5 reconstruction
891  lim_mem(1,:) = par_id_mem
892  lim_mem(2,:) = [1,grid_x%n(2)]
893  lim_mem(3,:) = [grid_x%i_min,grid_x%i_max]
894  lim_mem(4,:) = lim_sec_x_loc
895  if (present(lim_pos)) then
896  lim_mem(2,:) = lim_pos(2,:)
897  lim_mem(3,:) = lim_mem(3,:) + lim_pos(3,1) - 1 ! take into account the grid limits (relative to position subset)
898  end if
899 
900  ! RE_U_0
901  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
902  &'RE_U_0',rich_lvl=id,lim_loc=lim_mem)
903  chckerr('')
904  call conv_1d2nd(var_1d,dum_4d)
905  x%U_0(par_id(1):par_id(2):par_id(3),:,:,:) = dum_4d
906  call dealloc_var_1d(var_1d)
907 
908  ! IM_U_0
909  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
910  &'IM_U_0',rich_lvl=id,lim_loc=lim_mem)
911  chckerr('')
912  call conv_1d2nd(var_1d,dum_4d)
913  x%U_0(par_id(1):par_id(2):par_id(3),:,:,:) = &
914  x%U_0(par_id(1):par_id(2):par_id(3),:,:,:) + iu*dum_4d
915  call dealloc_var_1d(var_1d)
916 
917  ! RE_U_1
918  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
919  &'RE_U_1',rich_lvl=id,lim_loc=lim_mem)
920  chckerr('')
921  call conv_1d2nd(var_1d,dum_4d)
922  x%U_1(par_id(1):par_id(2):par_id(3),:,:,:) = dum_4d
923  call dealloc_var_1d(var_1d)
924 
925  ! IM_U_1
926  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
927  &'IM_U_1',rich_lvl=id,lim_loc=lim_mem)
928  chckerr('')
929  call conv_1d2nd(var_1d,dum_4d)
930  x%U_1(par_id(1):par_id(2):par_id(3),:,:,:) = &
931  x%U_1(par_id(1):par_id(2):par_id(3),:,:,:) + iu*dum_4d
932  call dealloc_var_1d(var_1d)
933 
934  ! RE_DU_0
935  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
936  &'RE_DU_0',rich_lvl=id,lim_loc=lim_mem)
937  chckerr('')
938  call conv_1d2nd(var_1d,dum_4d)
939  x%DU_0(par_id(1):par_id(2):par_id(3),:,:,:) = dum_4d
940  call dealloc_var_1d(var_1d)
941 
942  ! IM_DU_0
943  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
944  &'IM_DU_0',rich_lvl=id,lim_loc=lim_mem)
945  chckerr('')
946  call conv_1d2nd(var_1d,dum_4d)
947  x%DU_0(par_id(1):par_id(2):par_id(3),:,:,:) = &
948  x%DU_0(par_id(1):par_id(2):par_id(3),:,:,:) + iu*dum_4d
949  call dealloc_var_1d(var_1d)
950 
951  ! RE_DU_1
952  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
953  &'RE_DU_1',rich_lvl=id,lim_loc=lim_mem)
954  chckerr('')
955  call conv_1d2nd(var_1d,dum_4d)
956  x%DU_1(par_id(1):par_id(2):par_id(3),:,:,:) = dum_4d
957  call dealloc_var_1d(var_1d)
958 
959  ! IM_DU_1
960  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
961  &'IM_DU_1',rich_lvl=id,lim_loc=lim_mem)
962  chckerr('')
963  call conv_1d2nd(var_1d,dum_4d)
964  x%DU_1(par_id(1):par_id(2):par_id(3),:,:,:) = &
965  x%DU_1(par_id(1):par_id(2):par_id(3),:,:,:) + iu*dum_4d
966  call dealloc_var_1d(var_1d)
967  end do
968  end function reconstruct_pb3d_x_1
969 
990  integer function reconstruct_pb3d_x_2(mds,grid_X,X,data_name,rich_lvl,&
991  &tot_rich,lim_sec_X,lim_pos,is_field_averaged) result(ierr)
992  use num_vars, only: pb3d_name
993  use hdf5_ops, only: read_hdf5_arr
994  use pb3d_utilities, only: conv_1d2nd
995  use x_utilities, only: get_sec_x_range
996 
997  character(*), parameter :: rout_name = 'reconstruct_PB3D_X_2'
998 
999  ! input / output
1000  type(modes_type), intent(in) :: mds
1001  type(grid_type), intent(in) :: grid_x
1002  type(x_2_type), intent(inout) :: x
1003  character(len=*), intent(in) :: data_name
1004  integer, intent(in), optional :: rich_lvl
1005  logical, intent(in), optional :: tot_rich
1006  integer, intent(in), optional :: lim_sec_x(2,2)
1007  integer, intent(in), optional :: lim_pos(3,2)
1008  logical, intent(in), optional :: is_field_averaged
1009 
1010  ! local variables
1011  type(var_1d_type) :: var_1d ! 1D var
1012  logical :: is_field_averaged_loc ! local is_field_averaged
1013  integer :: id ! counter
1014  integer :: m, k ! counters
1015  integer :: sxr_loc(2,2) ! local secondary X limits for symmetric and asymmetric vars
1016  integer :: sxr_tot(2,2) ! total secondary X limits for symmetric and asymmetric vars
1017  integer :: nn_mod_loc(2) ! local nr. of modes for symmetric and asymmetric vars
1018  integer :: lim_mem(4,2,2) ! memory limits for variables for symmetric and asymmetric vars
1019  integer :: rich_lvl_loc ! local rich_lvl
1020  integer :: par_id(3) ! parallel indices (start, end, stride)
1021  integer :: par_id_mem(2) ! parallel indices (start, end) in memory (stride 1)
1022  integer :: par_lim(2) ! parallel limits
1023  integer :: rich_id(2) ! richardson level indices (start, end)
1024  logical :: read_this(2) ! whether symmetric and asymmetric variables need to be read
1025 
1026  ! initialize ierr
1027  ierr = 0
1028 
1029  ! set up local rich_lvl
1030  rich_lvl_loc = 0
1031  if (present(rich_lvl)) rich_lvl_loc = rich_lvl
1032 
1033  ! set up local is_field_averaged
1034  is_field_averaged_loc = .false.
1035  if (present(is_field_averaged)) is_field_averaged_loc = &
1036  &is_field_averaged
1037 
1038  ! setup rich_id
1039  rich_id = setup_rich_id(rich_lvl_loc,tot_rich)
1040 
1041  ! set up parallel limits
1042  par_lim = [1,grid_x%n(1)]
1043  if (present(lim_pos)) then
1044  where (lim_pos(1,:).ge.0) par_lim = lim_pos(1,:)
1045  end if
1046 
1047  ! create X
1048  call x%init(mds,grid_x,lim_sec_x,is_field_averaged)
1049 
1050  ! restore looping over richardson levels
1051  do id = rich_id(2),rich_id(1),-1
1052  ! setup par_id
1053  if (is_field_averaged_loc) then
1054  par_id = [1,1,1] ! only first element
1055  par_id_mem = [1,1]
1056  else
1057  par_id = setup_par_id(grid_x,rich_lvl_loc,id,tot_rich=tot_rich,&
1058  &par_lim=par_lim,par_id_mem=par_id_mem)
1059  end if
1060 
1061  ! loop over second dimension (horizontal)
1062  do m = 1,x%n_mod(2)
1063  ! get contiguous range of modes of this m
1064  call get_sec_x_range(sxr_loc(:,1),sxr_tot(:,1),m,.true.,&
1065  &lim_sec_x)
1066  call get_sec_x_range(sxr_loc(:,2),sxr_tot(:,2),m,.false.,&
1067  &lim_sec_x)
1068  nn_mod_loc = sxr_loc(2,:)-sxr_loc(1,:)+1
1069  read_this = .false.
1070  do k = 1,2
1071  if (sxr_loc(1,k).le.sxr_loc(2,k)) read_this(k) = .true. ! a bound is found
1072  end do
1073 
1074  ! set up local limits for HDF5 reconstruction of this m
1075  ! Note: It are the indices in total matrix sXr_tot that
1076  ! correspond to the local limits. "tot" just refers the to fact
1077  ! that they are valid for a submatrix of the total matrix; They
1078  ! have been set up using the local grid_X limits as well.
1079  lim_mem(1,:,1) = par_id_mem
1080  lim_mem(2,:,1) = [1,grid_x%n(2)]
1081  lim_mem(3,:,1) = [grid_x%i_min,grid_x%i_max]
1082  lim_mem(4,:,1) = [sxr_tot(1,1),sxr_tot(2,1)]
1083  lim_mem(1,:,2) = par_id_mem
1084  lim_mem(2,:,2) = [1,grid_x%n(2)]
1085  lim_mem(3,:,2) = [grid_x%i_min,grid_x%i_max]
1086  lim_mem(4,:,2) = [sxr_tot(1,2),sxr_tot(2,2)]
1087  if (present(lim_pos)) then
1088  lim_mem(2,:,1) = lim_pos(2,:)
1089  lim_mem(3,:,1) = lim_mem(3,:,1) + lim_pos(3,1) - 1 ! take into account the grid limits (relative to position subset)
1090  lim_mem(2,:,2) = lim_pos(2,:)
1091  lim_mem(3,:,2) = lim_mem(3,:,2) + lim_pos(3,1) - 1 ! take into account the grid limits (relative to position subset)
1092  end if
1093 
1094  if (read_this(1)) then
1095  ! RE_PV_0
1096  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1097  &'RE_PV_0',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1098  chckerr('')
1099  call conv_1d2nd(var_1d,dum_4d)
1100  x%PV_0(par_id(1):par_id(2):par_id(3),:,:,&
1101  &sxr_loc(1,1):sxr_loc(2,1)) = dum_4d
1102  call dealloc_var_1d(var_1d)
1103 
1104  ! IM_PV_0
1105  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1106  &'IM_PV_0',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1107  chckerr('')
1108  call conv_1d2nd(var_1d,dum_4d)
1109  x%PV_0(par_id(1):par_id(2):par_id(3),:,:,&
1110  &sxr_loc(1,1):sxr_loc(2,1)) = &
1111  &x%PV_0(par_id(1):par_id(2):par_id(3),:,:,&
1112  &sxr_loc(1,1):sxr_loc(2,1)) + iu*dum_4d
1113  call dealloc_var_1d(var_1d)
1114 
1115  ! RE_PV_2
1116  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1117  &'RE_PV_2',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1118  chckerr('')
1119  call conv_1d2nd(var_1d,dum_4d)
1120  x%PV_2(par_id(1):par_id(2):par_id(3),:,:,&
1121  &sxr_loc(1,1):sxr_loc(2,1)) = dum_4d
1122  call dealloc_var_1d(var_1d)
1123 
1124  ! IM_PV_2
1125  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1126  &'IM_PV_2',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1127  chckerr('')
1128  call conv_1d2nd(var_1d,dum_4d)
1129  x%PV_2(par_id(1):par_id(2):par_id(3),:,:,&
1130  &sxr_loc(1,1):sxr_loc(2,1)) = &
1131  &x%PV_2(par_id(1):par_id(2):par_id(3),:,:,&
1132  &sxr_loc(1,1):sxr_loc(2,1)) + iu*dum_4d
1133  call dealloc_var_1d(var_1d)
1134 
1135  ! RE_KV_0
1136  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1137  &'RE_KV_0',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1138  chckerr('')
1139  call conv_1d2nd(var_1d,dum_4d)
1140  x%KV_0(par_id(1):par_id(2):par_id(3),:,:,&
1141  &sxr_loc(1,1):sxr_loc(2,1)) = dum_4d
1142  call dealloc_var_1d(var_1d)
1143 
1144  ! IM_KV_0
1145  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1146  &'IM_KV_0',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1147  chckerr('')
1148  call conv_1d2nd(var_1d,dum_4d)
1149  x%KV_0(par_id(1):par_id(2):par_id(3),:,:,&
1150  &sxr_loc(1,1):sxr_loc(2,1)) = &
1151  &x%KV_0(par_id(1):par_id(2):par_id(3),:,:,&
1152  &sxr_loc(1,1):sxr_loc(2,1)) + iu*dum_4d
1153  call dealloc_var_1d(var_1d)
1154 
1155  ! RE_KV_2
1156  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1157  &'RE_KV_2',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1158  chckerr('')
1159  call conv_1d2nd(var_1d,dum_4d)
1160  x%KV_2(par_id(1):par_id(2):par_id(3),:,:,&
1161  &sxr_loc(1,1):sxr_loc(2,1)) = dum_4d
1162  call dealloc_var_1d(var_1d)
1163 
1164  ! IM_KV_2
1165  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1166  &'IM_KV_2',rich_lvl=id,lim_loc=lim_mem(:,:,1))
1167  chckerr('')
1168  call conv_1d2nd(var_1d,dum_4d)
1169  x%KV_2(par_id(1):par_id(2):par_id(3),:,:,&
1170  &sxr_loc(1,1):sxr_loc(2,1)) = &
1171  &x%KV_2(par_id(1):par_id(2):par_id(3),:,:,&
1172  &sxr_loc(1,1):sxr_loc(2,1)) + iu*dum_4d
1173  call dealloc_var_1d(var_1d)
1174  end if
1175 
1176  if (read_this(2)) then
1177  ! RE_PV_1
1178  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1179  &'RE_PV_1',rich_lvl=id,lim_loc=lim_mem(:,:,2))
1180  chckerr('')
1181  call conv_1d2nd(var_1d,dum_4d)
1182  x%PV_1(par_id(1):par_id(2):par_id(3),:,:,&
1183  &sxr_loc(1,2):sxr_loc(2,2)) = dum_4d
1184  call dealloc_var_1d(var_1d)
1185 
1186  ! IM_PV_1
1187  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1188  &'IM_PV_1',rich_lvl=id,lim_loc=lim_mem(:,:,2))
1189  chckerr('')
1190  call conv_1d2nd(var_1d,dum_4d)
1191  x%PV_1(par_id(1):par_id(2):par_id(3),:,:,&
1192  &sxr_loc(1,2):sxr_loc(2,2)) = &
1193  &x%PV_1(par_id(1):par_id(2):par_id(3),:,:,&
1194  &sxr_loc(1,2):sxr_loc(2,2)) + iu*dum_4d
1195  call dealloc_var_1d(var_1d)
1196 
1197  ! RE_KV_1
1198  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1199  &'RE_KV_1',rich_lvl=id,lim_loc=lim_mem(:,:,2))
1200  chckerr('')
1201  call conv_1d2nd(var_1d,dum_4d)
1202  x%KV_1(par_id(1):par_id(2):par_id(3),:,:,&
1203  &sxr_loc(1,2):sxr_loc(2,2)) = dum_4d
1204  call dealloc_var_1d(var_1d)
1205 
1206  ! IM_KV_1
1207  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1208  &'IM_KV_1',rich_lvl=id,lim_loc=lim_mem(:,:,2))
1209  chckerr('')
1210  call conv_1d2nd(var_1d,dum_4d)
1211  x%KV_1(par_id(1):par_id(2):par_id(3),:,:,&
1212  &sxr_loc(1,2):sxr_loc(2,2)) = &
1213  &x%KV_1(par_id(1):par_id(2):par_id(3),:,:,&
1214  &sxr_loc(1,2):sxr_loc(2,2)) + iu*dum_4d
1215  call dealloc_var_1d(var_1d)
1216  end if
1217  end do
1218  end do
1219  end function reconstruct_pb3d_x_2
1220 
1227  integer function reconstruct_pb3d_vac(vac,data_name,rich_lvl) result(ierr)
1229  use hdf5_ops, only: read_hdf5_arr
1230  use pb3d_utilities, only: conv_1d2nd
1231  use x_utilities, only: get_sec_x_range
1232 
1233  character(*), parameter :: rout_name = 'reconstruct_PB3D_vac'
1234 
1235  ! input / output
1236  type(vac_type), intent(inout) :: vac
1237  character(len=*), intent(in) :: data_name
1238  integer, intent(in), optional :: rich_lvl
1239 
1240  ! local variables
1241  type(var_1d_type) :: var_1d ! 1D var
1242  integer :: rich_lvl_loc ! local rich_lvl
1243  integer :: style ! style of vacuum
1244  integer :: n_bnd ! number of terms in boundary
1245  integer :: prim_x ! primary mode number
1246  integer :: n_ang(2) ! number of angles (1) and number of field lines (2)
1247  real(dp) :: jq ! iota or q
1248 
1249  ! initialize ierr
1250  ierr = 0
1251 
1252  ! set up local rich_lvl
1253  rich_lvl_loc = 0
1254  if (present(rich_lvl)) rich_lvl_loc = rich_lvl
1255 
1256  ! misc_vac
1257  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'misc_vac',&
1258  &rich_lvl=rich_lvl_loc)
1259  chckerr('')
1260  call conv_1d2nd(var_1d,dum_1d)
1261  style = nint(dum_1d(1))
1262  n_bnd = nint(dum_1d(2))
1263  prim_x = nint(dum_1d(3))
1264  n_ang = nint(dum_1d(4:5))
1265  jq = dum_1d(6)
1266  call dealloc_var_1d(var_1d)
1267 
1268  ! create vac
1269  ierr = vac%init(style,n_bnd,prim_x,n_ang,jq)
1270  chckerr('')
1271 
1272  ! secondary mode numbers
1273  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1274  &'sec_X',rich_lvl=rich_lvl_loc)
1275  chckerr('')
1276  call conv_1d2nd(var_1d,dum_1d)
1277  vac%sec_X = nint(dum_1d)
1278  deallocate(dum_1d)
1279  call dealloc_var_1d(var_1d)
1280 
1281  ! norm
1282  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1283  &'norm',rich_lvl=rich_lvl_loc)
1284  chckerr('')
1285  call conv_1d2nd(var_1d,dum_2d)
1286  vac%norm = dum_2d
1287  call dealloc_var_1d(var_1d)
1288 
1289  ! x_vec
1290  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1291  &'x_vec',rich_lvl=rich_lvl_loc)
1292  chckerr('')
1293  call conv_1d2nd(var_1d,dum_2d)
1294  vac%x_vec = dum_2d
1295  call dealloc_var_1d(var_1d)
1296 
1297  ! copy variables specific to style
1298  select case (vac%style)
1299  case (1) ! field-line 3-D
1300  ! h_fac
1301  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1302  &'h_fac',rich_lvl=rich_lvl_loc)
1303  chckerr('')
1304  call conv_1d2nd(var_1d,dum_2d)
1305  vac%h_fac = dum_2d
1306  call dealloc_var_1d(var_1d)
1307  case (2) ! axisymmetric
1308  ! ang
1309  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1310  &'ang',rich_lvl=rich_lvl_loc)
1311  chckerr('')
1312  call conv_1d2nd(var_1d,dum_2d)
1313  vac%ang = dum_2d
1314  call dealloc_var_1d(var_1d)
1315 
1316  ! dnorm
1317  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1318  &'dnorm',rich_lvl=rich_lvl_loc)
1319  chckerr('')
1320  call conv_1d2nd(var_1d,dum_2d)
1321  vac%dnorm = dum_2d
1322  call dealloc_var_1d(var_1d)
1323  end select
1324 
1325  if (rank.eq.n_procs-1) then
1326  ! RE_res
1327  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1328  &'RE_res',rich_lvl=rich_lvl_loc)
1329  chckerr('')
1330  call conv_1d2nd(var_1d,dum_2d)
1331  vac%res = dum_2d
1332  call dealloc_var_1d(var_1d)
1333 
1334  ! IM_res
1335  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),&
1336  &'IM_res',rich_lvl=rich_lvl_loc)
1337  chckerr('')
1338  call conv_1d2nd(var_1d,dum_2d)
1339  vac%res = vac%res + iu*dum_2d
1340  call dealloc_var_1d(var_1d)
1341  end if
1342  end function reconstruct_pb3d_vac
1343 
1357  integer function reconstruct_pb3d_sol(mds,grid_sol,sol,data_name,rich_lvl,&
1358  &lim_sec_sol,lim_pos) result(ierr)
1359  use num_vars, only: pb3d_name
1360  use hdf5_ops, only: read_hdf5_arr
1361  use x_vars, only: n_mod_x
1362  use pb3d_utilities, only: conv_1d2nd
1363 
1364  character(*), parameter :: rout_name = 'reconstruct_PB3D_sol'
1365 
1366  ! input / output
1367  type(modes_type), intent(in) :: mds
1368  type(grid_type), intent(in) :: grid_sol
1369  type(sol_type), intent(inout) :: sol
1370  character(len=*), intent(in) :: data_name
1371  integer, intent(in), optional :: rich_lvl
1372  integer, intent(in), optional :: lim_sec_sol(2)
1373  integer, intent(in), optional :: lim_pos(1,2)
1374 
1375  ! local variables
1376  type(var_1d_type) :: var_1d ! 1D variable
1377  integer :: lim_sec_sol_loc(2) ! local version of lim_sec_sol
1378  integer :: lim_mem(3,2) ! memory limits for variables
1379  integer :: n_ev ! nr. of Eigenvalues
1380 
1381  ! initialize ierr
1382  ierr = 0
1383 
1384  ! prepare
1385 
1386  ! set n_EV
1387  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'RE_sol_vec',&
1388  &rich_lvl=rich_lvl)
1389  chckerr('')
1390  n_ev = var_1d%tot_i_max(3)-var_1d%tot_i_min(3)+1
1391  call dealloc_var_1d(var_1d)
1392 
1393  ! set up local lim_sec_sol
1394  lim_sec_sol_loc = [1,n_mod_x]
1395  if (present(lim_sec_sol)) lim_sec_sol_loc = lim_sec_sol
1396 
1397  ! create solution
1398  call sol%init(mds,grid_sol,n_ev,lim_sec_sol)
1399 
1400  ! set up local limits for HDF5 reconstruction
1401  lim_mem(1,:) = lim_sec_sol_loc
1402  lim_mem(2,:) = [grid_sol%i_min,grid_sol%i_max]
1403  lim_mem(3,:) = [-1,-1]
1404  if (present(lim_pos)) then
1405  lim_mem(2,:) = lim_mem(2,:) + lim_pos(1,1) - 1 ! take into account the grid limits (relative to position subset)
1406  end if
1407 
1408  ! restore variables
1409 
1410  ! RE_sol_val
1411  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'RE_sol_val',&
1412  &rich_lvl=rich_lvl)
1413  chckerr('')
1414  call conv_1d2nd(var_1d,dum_1d)
1415  sol%val = dum_1d
1416  deallocate(dum_1d)
1417  call dealloc_var_1d(var_1d)
1418 
1419  ! IM_sol_val
1420  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'IM_sol_val',&
1421  &rich_lvl=rich_lvl)
1422  chckerr('')
1423  call conv_1d2nd(var_1d,dum_1d)
1424  sol%val = sol%val + iu*dum_1d
1425  deallocate(dum_1d)
1426  call dealloc_var_1d(var_1d)
1427 
1428  ! RE_sol_vec
1429  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'RE_sol_vec',&
1430  &rich_lvl=rich_lvl,lim_loc=lim_mem)
1431  chckerr('')
1432  call conv_1d2nd(var_1d,dum_3d)
1433  sol%vec = dum_3d
1434  deallocate(dum_3d)
1435  call dealloc_var_1d(var_1d)
1436 
1437  ! IM_sol_vec
1438  ierr = read_hdf5_arr(var_1d,pb3d_name,trim(data_name),'IM_sol_vec',&
1439  &rich_lvl=rich_lvl,lim_loc=lim_mem)
1440  chckerr('')
1441  call conv_1d2nd(var_1d,dum_3d)
1442  sol%vec = sol%vec + iu*dum_3d
1443  deallocate(dum_3d)
1444  call dealloc_var_1d(var_1d)
1445  end function reconstruct_pb3d_sol
1446 
1452  integer function get_pb3d_grid_size(n,grid_name,rich_lvl,tot_rich) &
1453  &result(ierr)
1454  use num_vars, only: pb3d_name
1455  use pb3d_utilities, only: conv_1d2nd
1456  use hdf5_ops, only: read_hdf5_arr
1457 
1458  character(*), parameter :: rout_name = 'get_PB3D_grid_size'
1459 
1460  ! input / output
1461  integer, intent(inout) :: n(3)
1462  character(len=*), intent(in) :: grid_name
1463  integer, intent(in), optional :: rich_lvl
1464  logical, intent(in), optional :: tot_rich
1465 
1466  ! local variables
1467  type(var_1d_type) :: var_1d ! 1D variable
1468  real(dp), allocatable :: dum_1d(:) ! dummy variables
1469 
1470  ! initialize ierr
1471  ierr = 0
1472 
1473  ! set n from HDF5 for local Richardson level
1474  n = 0
1475  ierr = read_hdf5_arr(var_1d,pb3d_name,'grid_'//trim(grid_name),&
1476  &'n',rich_lvl=rich_lvl)
1477  chckerr('')
1478  ! loop over all equilibrium jobs
1479  call conv_1d2nd(var_1d,dum_1d)
1480  n = nint(dum_1d)
1481  call dealloc_var_1d(var_1d)
1482 
1483  ! possibly only half of the points were saved in local Richardson level
1484  if (present(tot_rich) .and. present(rich_lvl)) then
1485  if (tot_rich .and. rich_lvl.gt.1) n(1) = n(1)*2+1
1486  end if
1487  end function get_pb3d_grid_size
1488 end module pb3d_ops
pb3d_ops::reconstruct_pb3d_sol
integer function, public reconstruct_pb3d_sol(mds, grid_sol, sol, data_name, rich_lvl, lim_sec_sol, lim_pos)
Reconstructs the solution variables from PB3D HDF5 output.
Definition: PB3D_ops.f90:1359
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
num_vars::max_name_ln
integer, parameter, public max_name_ln
maximum length of filenames
Definition: num_vars.f90:51
hdf5_vars::dealloc_var_1d
Deallocates 1D variable.
Definition: HDF5_vars.f90:68
x_utilities::get_sec_x_range
subroutine, public get_sec_x_range(sec_X_range_loc, sec_X_range_tot, m, sym, lim_sec_X)
Gets one of the the local ranges of contiguous tensorial perturbation variables to be printed or read...
Definition: X_utilities.f90:128
helena_vars::rbphi_h
real(dp), dimension(:,:), allocatable, public rbphi_h
Definition: HELENA_vars.f90:29
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
num_vars::use_normalization
logical, public use_normalization
whether to use normalization or not
Definition: num_vars.f90:115
eq_vars::vac_perm
real(dp), public vac_perm
either usual mu_0 (default) or normalized
Definition: eq_vars.f90:48
num_vars::bc_style
integer, dimension(2), public bc_style
style for BC left and right
Definition: num_vars.f90:94
num_vars
Numerical variables used by most other modules.
Definition: num_vars.f90:4
vmec_vars::mn_v
integer, dimension(:,:), allocatable, public mn_v
m and n of modes
Definition: VMEC_vars.f90:32
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
helena_vars::rot_t_h
real(dp), dimension(:,:), allocatable, public rot_t_h
rotational transform
Definition: HELENA_vars.f90:28
num_vars::norm_disc_prec_eq
integer, public norm_disc_prec_eq
precision for normal discretization for equilibrium
Definition: num_vars.f90:120
helena_vars::q_saf_h
real(dp), dimension(:,:), allocatable, public q_saf_h
safety factor
Definition: HELENA_vars.f90:27
hdf5_ops
Operations on HDF5 and XDMF variables.
Definition: HDF5_ops.f90:27
helena_vars::flux_p_h
real(dp), dimension(:,:), allocatable, public flux_p_h
poloidal flux
Definition: HELENA_vars.f90:24
num_vars::max_njq_change
real(dp), public max_njq_change
maximum change of prim. mode number times saf. fac. / rot. transf. when using X_style 2 (fast)
Definition: num_vars.f90:119
vmec_vars::b_v_c
real(dp), dimension(:,:), allocatable, public b_v_c
Coeff. of magnitude of B in sine series (HM and FM)
Definition: VMEC_vars.f90:50
x_vars::min_r_sol
real(dp), public min_r_sol
min. normal range for pert.
Definition: X_vars.f90:135
vmec_vars::j_v_sup_int
real(dp), dimension(:,:), allocatable, public j_v_sup_int
Integrated poloidal and toroidal current (FM)
Definition: VMEC_vars.f90:52
hdf5_ops::read_hdf5_arr
integer function, public read_hdf5_arr(var, PB3D_name, head_name, var_name, rich_lvl, disp_info, lim_loc)
Reads a PB3D output file in HDF5 format.
Definition: HDF5_ops.f90:1530
eq_vars::t_0
real(dp), public t_0
derived normalization constant for nondimensionalization
Definition: eq_vars.f90:47
num_vars::iu
complex(dp), parameter, public iu
complex unit
Definition: num_vars.f90:85
num_vars::n_procs
integer, public n_procs
nr. of MPI processes
Definition: num_vars.f90:69
num_vars::norm_disc_prec_sol
integer, public norm_disc_prec_sol
precision for normal discretization for solution
Definition: num_vars.f90:122
hdf5_vars
Variables pertaining to HDF5 and XDMF.
Definition: HDF5_vars.f90:4
grid_vars::n_r_eq
integer, public n_r_eq
nr. of normal points in equilibrium grid
Definition: grid_vars.f90:20
vmec_vars::jac_v_s
real(dp), dimension(:,:,:), allocatable, public jac_v_s
Coeff. of in cosine series (HM and FM) and norm. deriv.
Definition: VMEC_vars.f90:46
num_vars::x_style
integer, public x_style
style for secondary mode numbers (1: prescribed, 2: fast)
Definition: num_vars.f90:95
num_vars::debug_version
logical, public debug_version
debug version used
Definition: num_vars.f90:62
str_utilities
Operations on strings.
Definition: str_utilities.f90:4
vmec_vars::b_v_sub_s
real(dp), dimension(:,:,:), allocatable, public b_v_sub_s
Coeff. of B_i in cosine series (r,theta,phi) (FM)
Definition: VMEC_vars.f90:48
vmec_vars::flux_p_v
real(dp), dimension(:,:), allocatable, public flux_p_v
poloidal flux
Definition: VMEC_vars.f90:35
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_vars::grid_type
Type for grids.
Definition: grid_vars.f90:59
eq_vars::max_flux_e
real(dp), public max_flux_e
max. flux in Equilibrium coordinates, set in calc_norm_range_PB3D_in
Definition: eq_vars.f90:49
str_utilities::r2strt
elemental character(len=max_str_ln) function, public r2strt(k)
Convert a real (double) to string.
Definition: str_utilities.f90:54
num_vars::prog_style
integer, public prog_style
program style (1: PB3D, 2: PB3D_POST)
Definition: num_vars.f90:53
pb3d_ops::reconstruct_pb3d_in
integer function, public reconstruct_pb3d_in(data_name)
Reconstructs the input variables from PB3D HDF5 output.
Definition: PB3D_ops.f90:41
eq_vars::r_0
real(dp), public r_0
independent normalization constant for nondimensionalization
Definition: eq_vars.f90:42
pb3d_ops::reconstruct_pb3d_vac
integer function, public reconstruct_pb3d_vac(vac, data_name, rich_lvl)
Reconstructs the vacuum variables from PB3D HDF5 output.
Definition: PB3D_ops.f90:1228
eq_vars::eq_1_type
flux equilibrium type
Definition: eq_vars.f90:63
pb3d_utilities::setup_rich_id
integer function, dimension(2), public setup_rich_id(rich_lvl_max, tot_rich)
Returns richardson id.
Definition: PB3D_utilities.f90:234
vac_vars
Variables pertaining to the vacuum quantities.
Definition: vac_vars.f90:4
sol_vars::sol_type
solution type
Definition: sol_vars.f90:30
helena_vars::h_h_33
real(dp), dimension(:,:), allocatable, public h_h_33
upper metric factor (1 / gem12)
Definition: HELENA_vars.f90:32
helena_vars::r_h
real(dp), dimension(:,:), allocatable, public r_h
major radius (xout)
Definition: HELENA_vars.f90:33
grid_vars::n_r_sol
integer, public n_r_sol
nr. of normal points in solution grid
Definition: grid_vars.f90:22
pb3d_ops::reconstruct_pb3d_eq_2
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.
Definition: PB3D_ops.f90:695
pb3d_utilities::setup_par_id
integer function, dimension(3), public setup_par_id(grid, rich_lvl_max, rich_lvl_loc, tot_rich, par_lim, par_id_mem)
Setup parallel id.
Definition: PB3D_utilities.f90:181
helena_vars::flux_t_h
real(dp), dimension(:,:), allocatable, public flux_t_h
toroidal flux
Definition: HELENA_vars.f90:25
num_vars::norm_disc_style_sol
integer, public norm_disc_style_sol
style for normal discretization for solution (1: central fin. diff., 2: left fin. diff....
Definition: num_vars.f90:123
helena_vars::h_h_12
real(dp), dimension(:,:), allocatable, public h_h_12
upper metric factor (gem12)
Definition: HELENA_vars.f90:31
vmec_vars::b_v_sub_c
real(dp), dimension(:,:,:), allocatable, public b_v_sub_c
Coeff. of B_i in sine series (r,theta,phi) (FM)
Definition: VMEC_vars.f90:47
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::ev_bc
real(dp), public ev_bc
value of artificial Eigenvalue for boundary condition
Definition: num_vars.f90:116
helena_vars::z_h
real(dp), dimension(:,:), allocatable, public z_h
height (yout)
Definition: HELENA_vars.f90:34
hdf5_vars::var_1d_type
1D equivalent of multidimensional variables, used for internal HDF5 storage.
Definition: HDF5_vars.f90:48
grid_vars::n_r_in
integer, public n_r_in
nr. of normal points in input grid
Definition: grid_vars.f90:19
helena_vars::pres_h
real(dp), dimension(:,:), allocatable, public pres_h
pressure profile
Definition: HELENA_vars.f90:26
x_vars::x_2_type
tensorial perturbation type
Definition: X_vars.f90:81
vmec_vars::flux_t_v
real(dp), dimension(:,:), allocatable, public flux_t_v
toroidal flux
Definition: VMEC_vars.f90:34
num_vars::eq_style
integer, public eq_style
either 1 (VMEC) or 2 (HELENA)
Definition: num_vars.f90:89
vmec_vars::pres_v
real(dp), dimension(:,:), allocatable, public pres_v
pressure
Definition: VMEC_vars.f90:36
num_vars::rho_style
integer, public rho_style
style for equilibrium density profile
Definition: num_vars.f90:90
x_vars
Variables pertaining to the perturbation quantities.
Definition: X_vars.f90:4
num_vars::use_pol_flux_e
logical, public use_pol_flux_e
whether poloidal flux is used in E coords.
Definition: num_vars.f90:113
eq_vars::max_flux_f
real(dp), public max_flux_f
max. flux in Flux coordinates, set in calc_norm_range_PB3D_in
Definition: eq_vars.f90:50
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
num_vars::norm_style
integer, public norm_style
style for normalization
Definition: num_vars.f90:92
helena_vars::h_h_11
real(dp), dimension(:,:), allocatable, public h_h_11
upper metric factor (gem11)
Definition: HELENA_vars.f90:30
eq_vars::rho_0
real(dp), public rho_0
independent normalization constant for nondimensionalization
Definition: eq_vars.f90:44
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
num_vars::max_deriv
integer, parameter, public max_deriv
highest derivatives for metric factors in Flux coords.
Definition: num_vars.f90:52
num_vars::alpha_style
integer, public alpha_style
style for alpha (1: one field line, many turns, 2: many field lines, one turn)
Definition: num_vars.f90:100
vmec_vars::rot_t_v
real(dp), dimension(:,:), allocatable, public rot_t_v
rotational transform
Definition: VMEC_vars.f90:37
vmec_vars::b_v_s
real(dp), dimension(:,:), allocatable, public b_v_s
Coeff. of magnitude of B in cosine series (HM and FM)
Definition: VMEC_vars.f90:51
num_vars::k_style
integer, public k_style
style for kinetic energy
Definition: num_vars.f90:93
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
vmec_vars::z_v_s
real(dp), dimension(:,:,:), allocatable, public z_v_s
Coeff. of in cosine series (FM) and norm. deriv.
Definition: VMEC_vars.f90:42
helena_vars
Variables that have to do with HELENA quantities.
Definition: HELENA_vars.f90:4
pb3d_ops
Operations on PB3D output.
Definition: PB3D_ops.f90:8
num_vars::pi
real(dp), parameter, public pi
Definition: num_vars.f90:83
grid_vars::max_alpha
real(dp), public max_alpha
max. of field-line label [ ] in field-aligned grid
Definition: grid_vars.f90:27
pb3d_ops::reconstruct_pb3d_eq_1
integer function, public reconstruct_pb3d_eq_1(grid_eq, eq, data_name, lim_pos)
Reconstructs the equilibrium variables from PB3D HDF5 output.
Definition: PB3D_ops.f90:597
vmec_vars::q_saf_v
real(dp), dimension(:,:), allocatable, public q_saf_v
safety factor
Definition: VMEC_vars.f90:38
num_vars::pb3d_name
character(len=max_str_ln), public pb3d_name
name of PB3D output file
Definition: num_vars.f90:139
vac_vars::vac_type
vacuum type
Definition: vac_vars.f90:46
vmec_vars::r_v_c
real(dp), dimension(:,:,:), allocatable, public r_v_c
Coeff. of in sine series (FM) and norm. deriv.
Definition: VMEC_vars.f90:39
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
helena_vars::ias
integer, public ias
0 if top-bottom symmetric, 1 if not
Definition: HELENA_vars.f90:36
messages::lvl_ud
subroutine, public lvl_ud(inc)
Increases/decreases lvl of output.
Definition: messages.f90:254
num_vars::min_pb3d_version
real(dp), parameter, public min_pb3d_version
minimum PB3D version for POST
Definition: num_vars.f90:60
pb3d_utilities::conv_1d2nd
Converts 1-D to n-D variables.
Definition: PB3D_utilities.f90:21
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
vmec_vars
Variables that concern the output of VMEC.
Definition: VMEC_vars.f90:4
vmec_vars::l_v_c
real(dp), dimension(:,:,:), allocatable, public l_v_c
Coeff. of in sine series (HM) and norm. deriv.
Definition: VMEC_vars.f90:43
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
eq_vars::psi_0
real(dp), public psi_0
derived normalization constant for nondimensionalization
Definition: eq_vars.f90:46
vmec_vars::r_v_s
real(dp), dimension(:,:,:), allocatable, public r_v_s
Coeff. of in cosine series (FM) and norm. deriv.
Definition: VMEC_vars.f90:40
vmec_vars::jac_v_c
real(dp), dimension(:,:,:), allocatable, public jac_v_c
Coeff. of in sine series (HM and FM) and norm. deriv.
Definition: VMEC_vars.f90:45
x_vars::x_1_type
vectorial perturbation type
Definition: X_vars.f90:51
sol_vars
Variables pertaining to the solution quantities.
Definition: sol_vars.f90:4
grid_vars::min_alpha
real(dp), public min_alpha
min. of field-line label [ ] in field-aligned grid
Definition: grid_vars.f90:26
helena_vars::chi_h
real(dp), dimension(:), allocatable, public chi_h
poloidal angle
Definition: HELENA_vars.f90:23
grid_vars::n_alpha
integer, public n_alpha
nr. of field-lines
Definition: grid_vars.f90:23
output_ops
Operations concerning giving output, on the screen as well as in output files.
Definition: output_ops.f90:5
eq_vars::pres_0
real(dp), public pres_0
independent normalization constant for nondimensionalization
Definition: eq_vars.f90:43
num_vars::rank
integer, public rank
MPI rank.
Definition: num_vars.f90:68
pb3d_ops::get_pb3d_grid_size
integer function, public get_pb3d_grid_size(n, grid_name, rich_lvl, tot_rich)
get grid size
Definition: PB3D_ops.f90:1454
vmec_vars::l_v_s
real(dp), dimension(:,:,:), allocatable, public l_v_s
Coeff. of in cosine series (HM) and norm. deriv.
Definition: VMEC_vars.f90:44
num_vars::matrix_slepc_style
integer, public matrix_slepc_style
style for matrix storage (1: sparse, 2: shell)
Definition: num_vars.f90:96
eq_vars::eq_2_type
metric equilibrium type
Definition: eq_vars.f90:114
helena_vars::nchi
integer, public nchi
nr. of poloidal points
Definition: HELENA_vars.f90:35
vmec_vars::z_v_c
real(dp), dimension(:,:,:), allocatable, public z_v_c
Coeff. of in sine series (FM) and norm. deriv.
Definition: VMEC_vars.f90:41
pb3d_utilities
Numerical utilities related to PB3D operations.
Definition: PB3D_utilities.f90:4
num_vars::ev_style
integer, public ev_style
determines the method used for solving an EV problem
Definition: num_vars.f90:88
eq_vars::b_0
real(dp), public b_0
derived normalization constant for nondimensionalization
Definition: eq_vars.f90:45
num_vars::norm_disc_prec_x
integer, public norm_disc_prec_x
precision for normal discretization for perturbation
Definition: num_vars.f90:121