PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Go to the documentation of this file.
5 #include <PB3D_macros.h>
16 character(len=max_str_ln),
allocatable :: command_arg(:)
19 character(len=max_str_ln),
allocatable ::
opt_args(:)
20 integer,
allocatable :: inc_args(:)
37 opt_args(11) =
'--invert_top_bottom_H'
38 opt_args(12) =
'-st_pc_factor_shift_type'
40 opt_args(14) =
'-st_pc_factor_mat_solver_package'
49 opt_args(23) =
'-st_pc_factor_mat_solver_package'
52 inc_args(7:25) = [0,0,0,0,0,1,1,1,1,0,1,1,1,0,1,1,1,0,1]
68 opt_args(5) =
'--do_execute_command_line'
70 inc_args(1:6) = [0,0,0,0,0,0]
81 character(*),
parameter :: rout_name =
'parse_args'
84 character(len=max_str_ln),
allocatable :: open_error(:)
85 character(len=max_str_ln),
allocatable :: open_help(:)
88 character(len=max_str_ln) :: dum_str
98 allocate(open_error(3))
99 allocate(open_help(6))
102 &
" USER_INPUT EQUILIBRIUM_INPUT [OPTIONS]"
103 open_error(3) =
"Try './"//
prog_name//
" --help' or './"//&
105 open_help(1) = open_error(2)
107 open_help(3) =
" USER_INPUT input provided by the &
109 open_help(4) =
" EQUILIBRIUM_INPUT output from VMEC or &
110 &HELENA, in .txt format"
111 open_help(5) =
" [OPTIONS] (optional) &
112 &command-line options"
116 allocate(open_error(3))
117 allocate(open_help(6))
120 &
" USER_INPUT PB3D_OUTPUT [OPTIONS]"
121 open_error(3) =
"Try './"//
prog_name//
" --help' or './"//&
123 open_help(1) = open_error(2)
125 open_help(3) =
" USER_INPUT input provided by the &
127 open_help(4) =
" PB3D_OUTPUT output from PB3D code, &
129 open_help(5) =
" [OPTIONS] (optional) &
130 &command-line options"
135 call writo(
"Parsing command line arguments")
138 numargs = command_argument_count()
140 if (numargs.gt.0)
call getarg(1,dum_str)
143 if (dum_str.eq.
'-h' .or. dum_str.eq.
'--help')
then
144 do id = 1,
size(open_help)
145 call writo(open_help(id))
149 else if (numargs.lt.min_args)
then
150 do id = 1,
size(open_error)
151 call writo(open_error(id))
156 allocate(command_arg(numargs))
158 command_arg(1) = dum_str
163 call getarg(iseq, command_arg(iseq))
166 call writo(
"Command line arguments parsed")
188 character(*),
parameter :: rout_name =
'open_input'
192 character(len=max_str_ln) :: err_msg
195 character(len=max_str_ln) :: file_ext
201 call writo(
"Opening files")
205 select case (prog_style)
208 input_name = command_arg(1)
209 open(unit=input_i,file=input_name,status=
'old',iostat=istat)
211 call writo(
'Input file "' // trim(input_name) &
214 call writo(
'No input file found. Default used')
220 eq_name = command_arg(2)
221 open(unit=eq_i,file=eq_name,status=
'old',iostat=ierr)
222 err_msg =
'No equilibrium file found.'
226 call writo(
'equilibrium file "' // trim(eq_name) &
235 do id = len(eq_name)-1,1,-1
236 if (eq_name(id:id).eq.
'.')
then
237 file_ext = eq_name(id+1:len(eq_name))
241 if (trim(file_ext).eq.
'nc')
then
246 read(eq_i,*,iostat=istat) first_int
253 if (eq_style.lt.1 .or. eq_style.gt.2)
then
255 call writo(
'Unable to recognize type of input file')
256 call writo(
'If you are trying to open HELENA files, &
257 &you need to have a version that outputs the &
258 &variable IAS. See tutorial!')
263 pb3d_name = prog_name//
'_'//trim(output_name)//
'.h5'
266 input_name = command_arg(1)
267 open(unit=input_i,file=input_name,status=
'old',iostat=istat)
269 call writo(
'Input file "' // trim(input_name) &
272 call writo(
'No input file found. Default used')
278 pb3d_name = command_arg(2)
279 open(unit=pb3d_i,file=pb3d_name,status=
'old',iostat=ierr)
280 err_msg =
'No PB3D file found.'
284 call writo(
'PB3D output file "' // trim(pb3d_name) &
290 if (numargs.gt.2)
then
291 call writo(
'applying options')
297 call writo(
'Files opened')
303 #if ( lwith_intel && !lwith_gnu)
308 logical :: opt_taken(size(
opt_args))
317 do while (id.le.numargs)
319 do while (jd.le.numopts .and. .not.opt_found)
320 if (trim(command_arg(id)).eq.trim(
opt_args(jd)))
then
322 if (opt_taken(jd))
then
323 call writo(
'Option "'//trim(command_arg(id))//&
324 &
'" already set, ignoring...',warning=.true.)
330 call writo(
'option test chosen')
333 call writo(
'option test not available. &
334 &Recompile with cpp flag &
335 &''ldebug''...',warning=.true.)
338 call writo(
'option no_plots chosen: &
342 call writo(
'option no_output chosen: &
346 call writo(
'option do_execute_command_line &
347 &chosen: execute_command_line enabled')
348 do_execute_command_line = .true.
351 call writo(
'option mem_info chosen: &
352 &memory usage is printed')
354 call writo(
'The PID of this process is: '//&
355 &trim(
i2str(getpid())))
356 print_mem_usage = .true.
359 call writo(
'option mem_info is only &
360 &available when compiled with &
361 &"ldebug"',warning=.true.)
365 select case (prog_style)
367 call apply_opt_pb3d(jd,id)
369 call apply_opt_post(jd,id)
372 opt_taken(jd) = .true.
374 id = id + inc_args(jd)
379 if (.not.opt_found)
then
380 call writo(
'Option "'//trim(command_arg(id))//
'" invalid',&
386 end subroutine read_opts
390 subroutine apply_opt_pb3d(opt_nr,arg_nr)
394 integer,
intent(in) :: opt_nr
395 integer,
intent(in) :: arg_nr
399 character(len=max_str_ln) :: opt_str
403 call writo(
'option no_guess chosen: Eigenfunction not &
404 &guessed from previous Richardson level')
407 call writo(
'option jump_to_sol chosen: Skip all possible &
408 &equilibrium and perturbation drivers for first &
412 if (eq_style.eq.2)
then
414 call writo(
'option export_HEL chosen: Exporting to &
418 call writo(
'Can only use export_HEL with 1 &
419 &process',warning=.true.)
423 call writo(
'Can only use export_HEL with HELENA',&
428 if (eq_style.eq.1)
then
429 call writo(
'option plot_VMEC_modes chosen: Plotting &
431 plot_vmec_modes = .true.
433 call writo(
'Can only output VMEC modes with VMEC',&
435 plot_vmec_modes = .false.
439 if (eq_style.eq.2)
then
440 call writo(
'option invert_top_bottom_H chosen: &
441 &inverting top and bottom of equilibrium')
442 invert_top_bottom_h = .true.
444 call writo(
'Can only invert top and bottom for HELENA',&
446 invert_top_bottom_h = .false.
449 call writo(
'Can only invert top and bottom for HELENA in &
450 &debug mode',warning=.true.)
451 invert_top_bottom_h = .false.
455 do id = 1,inc_args(opt_nr)
456 opt_str = trim(opt_str)//
' '//&
457 &trim(command_arg(arg_nr+id))
460 &trim(opt_str)//
'" passed to SLEPC')
462 call writo(
'Invalid option number',warning=.true.)
464 end subroutine apply_opt_pb3d
468 subroutine apply_opt_post(opt_nr,arg_nr)
472 integer,
intent(in) :: opt_nr
473 integer,
intent(in) :: arg_nr
480 call writo(
'option swap_angles chosen: theta and zeta are &
484 call writo(
'option compare_tor_pos chosen: Comparing &
485 &B, J and kappa at different toroidal positions')
487 compare_tor_pos = .true.
488 do id = 1,inc_args(opt_nr)
489 if (arg_nr+id.gt.
size(command_arg))
then
492 read(command_arg(arg_nr+id),*,iostat=ierr)
rz_0(id)
495 call writo(
'Did you provide an origin for the &
496 &geometrical poloidal angle in the form')
497 call writo(
'"--compare_tor_pos R_0, Z_0"?')
499 chckerr(
'failed to read RZ_0')
501 call writo(
'Origin for geometrical poloidal angle: ('//&
505 call writo(
'Invalid option number',warning=.true.)
507 end subroutine apply_opt_post
537 #if ( lwith_intel && !lwith_gnu)
541 character(*),
parameter :: rout_name =
'open_output'
546 logical :: file_exists
547 logical :: group_exists
548 character(len=max_str_ln) :: full_output_name
549 character(len=max_str_ln) :: err_msg
555 call writo(
'Attempting to open output files')
561 full_output_name = prog_name//
'_'//trim(output_name)//
'.txt'
564 if (rich_restart_lvl.gt.1 .and. prog_style.eq.1)
then
566 open(output_i,file=trim(full_output_name),status=
'old',&
567 &position=
'append',iostat=ierr)
568 chckerr(
'Failed to open output file')
571 call writo(
'log output file "'//trim(full_output_name)//&
572 &
'" reopened at number '//trim(
i2str(output_i)))
575 open(output_i,file=trim(full_output_name),status=
'replace',&
577 chckerr(
'Failed to open output file')
580 call writo(
'log output file "'//trim(full_output_name)//&
581 &
'" opened at number '//trim(
i2str(output_i)))
586 write(output_i,
'(A)',iostat=istat) trim(
temp_output(id))
592 full_output_name = prog_name//
'_'//trim(shell_commands_name)//
'.sh'
595 open(unit=shell_commands_i,file=trim(full_output_name),&
596 &status=
'replace',iostat=ierr)
597 chckerr(
'Failed to create shell command file')
600 write(shell_commands_i,
'(A)',iostat=istat)
'#!/bin/bash'
601 write(shell_commands_i,
'(A)',iostat=istat)
'# This file contains all &
602 &the shell commands from the '//trim(prog_name)//
' run'
603 close(shell_commands_i)
605 #if ( lwith_intel && !lwith_gnu)
606 istat = system(
'chmod +x '//trim(full_output_name))
608 call execute_command_line(
'chmod +x '//trim(full_output_name),&
613 call writo(
'shell commands script file "'//trim(full_output_name)//&
619 select case (prog_style)
622 inquire(file=pb3d_name,exist=file_exists,iostat=ierr)
623 chckerr(
'Failed to inquire about file')
624 if (file_exists)
then
626 &trim(
i2str(rich_restart_lvl)),group_exists)
629 group_exists = .false.
632 if (jump_to_sol)
then
633 if (.not.group_exists)
then
635 call writo(
'No integrated tensorial perturbation &
636 &quantitites found to jump over',alert=.true.)
637 err_msg =
'These quantitites need to be set up'
642 if (group_exists)
then
643 call writo(
'Integrated tensorial perturbation &
644 &quantitites were found for Richardson level '&
645 &//trim(
i2str(rich_restart_lvl))//
'.',&
648 call writo(
'Have you already calculated them and &
649 &are you re-running these simulations')
650 call writo(
'because you want to change the solution?')
651 call writo(
'If the solution method is the only &
652 &thing you are changing, you should')
653 call writo(
'consider using the command-line option &
659 if (rich_restart_lvl.eq.1)
then
675 open(
mem_usage_i,file=trim(full_output_name),status=
'replace',&
677 chckerr(
'Failed to open memory file')
680 write(
mem_usage_i,
'(A)',iostat=istat)
'# Rank [] Count [] &
681 & Time [s] Mem. [kB] Max. tot. Memory [kB] &
682 &Max. X. Memory [kB]'
686 call writo(
'memory usage data file "'//trim(full_output_name)//&
698 call writo(
'Output files opened')
706 call writo(
'Closing output files')
708 call writo(
'A log of the shell command is kept in the log file "'//&
711 &
call writo(
'Execute this script to do all the shell commands')
713 call writo(
'Output files closed')
logical, public swap_angles
swap angles theta and zeta in plots (only for POST)
logical, public invert_top_bottom_h
invert top and bottom for HELENA equilibria
integer, parameter, public dp
double precision
integer, parameter, public input_i
file number of input file
Numerical variables used by most other modules.
character(len=max_str_ln), public input_name
will hold the full name of the input file
integer, parameter, public max_str_ln
maximum length of strings
character(len=max_str_ln), dimension(:), allocatable, public opt_args
optional arguments that can be passed using –[name]
Variables concerning Richardson extrapolation.
elemental character(len=max_str_ln) function, public i2str(k)
Convert an integer to string.
Operations on HDF5 and XDMF variables.
integer, parameter, public shell_commands_i
file number of shell commands file
integer, public n_procs
nr. of MPI processes
Utilities pertaining to HDF5 and XDMF.
integer function, public open_input()
Open the input files.
integer, public lvl
determines the indenting. higher lvl = more indenting
integer, public prog_style
program style (1: PB3D, 2: PB3D_POST)
character(len=4), public prog_name
name of program, used for info
integer function, public create_output_hdf5(HDF5_name)
Creates an HDF5 output file.
logical, public ltest
whether or not to call the testing routines
integer function, public parse_args()
Parses the command line arguments.
character(len=14), parameter, public shell_commands_name
name of shell commands file
integer function, public probe_hdf5_group(HDF5_name, group_name, group_exists)
Probe HDF5 file for group existence.
logical, public compare_tor_pos
compare quantities at toroidal positions (only for POST)
subroutine, public init_files()
Initialize the variables for the module.
integer, public rich_restart_lvl
starting Richardson level (0: none [default])
logical, public export_hel
export HELENA
real(dp), dimension(2), public rz_0
origin of geometrical poloidal coordinate
logical, public print_mem_usage
print memory usage is printed
integer, parameter, public mem_usage_i
file number of memory usage file
integer, public eq_style
either 1 (VMEC) or 2 (HELENA)
logical, public no_guess
disable guessing Eigenfunction from previous level of Richardson
logical, public no_output
no output shown
character(len=max_str_ln), dimension(:), allocatable, public temp_output
temporary output, before output file is opened
elemental character(len=max_str_ln) function, public r2str(k)
Convert a real (double) to string.
subroutine, public writo(input_str, persistent, error, warning, alert)
Write output to file identified by output_i.
Numerical utilities related to giving output.
logical, public plot_vmec_modes
plot VMEC modes
Operations related to files !
character(len=9), parameter, public mem_usage_name
name of memory usage file
logical, public do_execute_command_line
call "execute_command_line" inside program
character(len=max_str_ln), public pb3d_name
name of PB3D output file
logical, public jump_to_sol
jump to solution
integer, parameter, public eq_i
file number of equilibrium file from VMEC or HELENA
subroutine, public lvl_ud(inc)
Increases/decreases lvl of output.
logical, public no_plots
no plots made
subroutine, public close_output()
Closes the output file.
integer, parameter, public pb3d_i
file number of PB3D output file
character(len=3), parameter, public output_name
name of output file
integer, public rank
MPI rank.
integer, parameter, public output_i
file number of output file
logical, public temp_output_active
true if temporary output is to be written in temp_output
character(len=max_str_ln), public eq_name
name of equilibrium file from VMEC or HELENA
integer function, public open_output()
Open the output files.