PB3D  [2.45]
Ideal linear high-n MHD stability in 3-D
eq_vars.f90
Go to the documentation of this file.
1 !------------------------------------------------------------------------------!
26 !------------------------------------------------------------------------------!
27 module eq_vars
28 #include <PB3D_macros.h>
29  use str_utilities
30  use messages
32  use grid_vars, only: grid_type
33 
34  implicit none
35  private
37 #if ldebug
39 #endif
40 
41  ! global variables
42  real(dp) :: r_0
43  real(dp) :: pres_0
44  real(dp) :: rho_0
45  real(dp) :: b_0
46  real(dp) :: psi_0
47  real(dp) :: t_0
49  real(dp) :: max_flux_e
50  real(dp) :: max_flux_f
51 #if ldebug
52 
53  integer :: n_alloc_eq_1s
55  integer :: n_alloc_eq_2s
56 #endif
57 
63  type, public :: eq_1_type
64  real(dp), allocatable :: pres_e(:,:)
65  real(dp), allocatable :: q_saf_e(:,:)
66  real(dp), allocatable :: rot_t_e(:,:)
67  real(dp), allocatable :: flux_p_e(:,:)
68  real(dp), allocatable :: flux_t_e(:,:)
69  real(dp), allocatable :: pres_fd(:,:)
70  real(dp), allocatable :: q_saf_fd(:,:)
71  real(dp), allocatable :: rot_t_fd(:,:)
72  real(dp), allocatable :: flux_p_fd(:,:)
73  real(dp), allocatable :: flux_t_fd(:,:)
74  real(dp), allocatable :: rho(:)
75 #if ldebug
76  real(dp) :: estim_mem_usage(2)
77 #endif
78  contains
80  procedure :: init => init_eq_1
82  procedure :: copy => copy_eq_1
84  procedure :: dealloc => dealloc_eq_1
85  end type
86 
114  type, public :: eq_2_type
115  ! coordinate variables R, Z and L (VMEC)
116  real(dp), allocatable :: r_e(:,:,:,:,:,:)
117  real(dp), allocatable :: z_e(:,:,:,:,:,:)
118  real(dp), allocatable :: l_e(:,:,:,:,:,:)
119  ! upper (h) and lower (g) metric factors
120  real(dp), allocatable :: g_c(:,:,:,:,:,:,:)
121  real(dp), allocatable :: g_e(:,:,:,:,:,:,:)
122  real(dp), allocatable :: h_e(:,:,:,:,:,:,:)
123  real(dp), allocatable :: g_f(:,:,:,:,:,:,:)
124  real(dp), allocatable :: h_f(:,:,:,:,:,:,:)
125  real(dp), allocatable :: g_fd(:,:,:,:,:,:,:)
126  real(dp), allocatable :: h_fd(:,:,:,:,:,:,:)
127  ! transformation matrices
128  real(dp), allocatable :: t_vc(:,:,:,:,:,:,:)
129  real(dp), allocatable :: t_ef(:,:,:,:,:,:,:)
130  real(dp), allocatable :: t_fe(:,:,:,:,:,:,:)
131  ! determinants of transformation matrices
132  real(dp), allocatable :: det_t_ef(:,:,:,:,:,:)
133  real(dp), allocatable :: det_t_fe(:,:,:,:,:,:)
134  ! Jacobians
135  real(dp), allocatable :: jac_e(:,:,:,:,:,:)
136  real(dp), allocatable :: jac_f(:,:,:,:,:,:)
137  real(dp), allocatable :: jac_fd(:,:,:,:,:,:)
138  ! derived variables
139  real(dp), allocatable :: s(:,:,:)
140  real(dp), allocatable :: kappa_n(:,:,:)
141  real(dp), allocatable :: kappa_g(:,:,:)
142  real(dp), allocatable :: sigma(:,:,:)
143 #if ldebug
144  real(dp) :: estim_mem_usage(2)
145 #endif
146  contains
148  procedure :: init => init_eq_2
150  procedure :: copy => copy_eq_2
152  procedure :: dealloc => dealloc_eq_2
153  end type
154 
155 contains
173  subroutine init_eq_1(eq,grid,setup_E,setup_F)
174  use num_vars, only: max_deriv
175 #if ldebug
176  use num_vars, only: print_mem_usage, rank
177 #endif
178 
179  ! input / output
180  class(eq_1_type), intent(inout) :: eq
181  type(grid_type), intent(in) :: grid
182  logical, intent(in), optional :: setup_E
183  logical, intent(in), optional :: setup_F
184 
185  ! local variables
186  integer :: loc_n_r, n ! local and total nr. of normal points
187  logical :: setup_E_loc, setup_F_loc ! local versions of setup_E and setup_F
188 #if ldebug
189  ! initialize memory usage
190  if (print_mem_usage) eq%estim_mem_usage = 0._dp
191 #endif
192 
193  ! setup local setup_E and setup_F
194  setup_e_loc = .true.
195  if (present(setup_e)) setup_e_loc = setup_e
196  setup_f_loc = .true.
197  if (present(setup_f)) setup_f_loc = setup_f
198 
199  ! set local variables
200  loc_n_r = grid%loc_n_r
201  n = grid%n(3)
202 
203  if (setup_e_loc) then
204  ! pres_E
205  allocate(eq%pres_E(loc_n_r,0:max_deriv+1))
206 
207  ! q_saf_E
208  allocate(eq%q_saf_E(loc_n_r,0:max_deriv+1))
209 
210  ! rot_t_E
211  allocate(eq%rot_t_E(loc_n_r,0:max_deriv+1))
212 
213  ! flux_p_E
214  allocate(eq%flux_p_E(loc_n_r,0:max_deriv+1))
215 
216  ! flux_t_E
217  allocate(eq%flux_t_E(loc_n_r,0:max_deriv+1))
218 
219 #if ldebug
220  ! set estimated memory usage
221  if (print_mem_usage) eq%estim_mem_usage(1) = &
222  &eq%estim_mem_usage(1) + loc_n_r*(max_deriv+2)*5
223 #endif
224  end if
225 
226  if (setup_f_loc) then
227  ! pres_FD
228  allocate(eq%pres_FD(loc_n_r,0:max_deriv+1))
229 
230  ! flux_p_FD
231  allocate(eq%flux_p_FD(loc_n_r,0:max_deriv+1))
232 
233  ! flux_t_FD
234  allocate(eq%flux_t_FD(loc_n_r,0:max_deriv+1))
235 
236  ! q_saf_FD
237  allocate(eq%q_saf_FD(loc_n_r,0:max_deriv+1))
238 
239  ! rot_t_FD
240  allocate(eq%rot_t_FD(loc_n_r,0:max_deriv+1))
241 
242  ! rho
243  allocate(eq%rho(grid%loc_n_r))
244 
245 #if ldebug
246  ! set estimated memory usage
247  if (print_mem_usage) eq%estim_mem_usage(2) = &
248  &eq%estim_mem_usage(2) + loc_n_r*((max_deriv+2)*5+1)
249 #endif
250  end if
251 
252 #if ldebug
253  ! increment n_alloc_eq_1s
255 
256  ! print memory usage
257  if (print_mem_usage) call writo('[rank '//trim(i2str(rank))//&
258  &' - Expected memory usage of eq_1: '//&
259  &trim(r2strt(eq%estim_mem_usage(1)*0.008))//' kB for E and '//&
260  &trim(r2strt(eq%estim_mem_usage(2)*0.008))//' kB for E]',&
261  &alert=.true.)
262 #endif
263  end subroutine init_eq_1
264 
272  subroutine init_eq_2(eq,grid,setup_E,setup_F) ! metric version
274 #if ldebug
275  use num_vars, only: print_mem_usage, rank
276 #endif
277 
278  ! input / output
279  class(eq_2_type), intent(inout) :: eq
280  type(grid_type), intent(in) :: grid
281  logical, intent(in), optional :: setup_E
282  logical, intent(in), optional :: setup_F
283 
284  ! local variables
285  logical :: setup_E_loc, setup_F_loc ! local versions of setup_E and setup_F
286  integer :: dims(3) ! dimensions
287 
288  ! set up local variables
289  dims = [grid%n(1),grid%n(2),grid%loc_n_r]
290 
291 #if ldebug
292  ! initialize memory usage
293  if (print_mem_usage) eq%estim_mem_usage = 0._dp
294 #endif
295 
296  ! setup local setup_E and setup_F
297  setup_e_loc = .true.
298  if (present(setup_e)) setup_e_loc = setup_e
299  setup_f_loc = .true.
300  if (present(setup_f)) setup_f_loc = setup_f
301 
302  if (setup_e_loc) then
303  ! initialize variables that are used for all equilibrium styles
304  ! g_E
305  allocate(eq%g_E(dims(1),dims(2),dims(3),6,&
306  &0:max_deriv,0:max_deriv,0:max_deriv))
307 
308  ! h_E
309  allocate(eq%h_E(dims(1),dims(2),dims(3),6,&
310  &0:max_deriv,0:max_deriv,0:max_deriv))
311 
312  ! g_F
313  allocate(eq%g_F(dims(1),dims(2),dims(3),6,&
314  &0:max_deriv,0:max_deriv,0:max_deriv))
315 
316  ! h_F
317  allocate(eq%h_F(dims(1),dims(2),dims(3),6,&
318  &0:max_deriv,0:max_deriv,0:max_deriv))
319 
320  ! jac_E
321  allocate(eq%jac_E(dims(1),dims(2),dims(3),&
322  &0:max_deriv,0:max_deriv,0:max_deriv))
323 
324  ! jac_F
325  allocate(eq%jac_F(dims(1),dims(2),dims(3),&
326  &0:max_deriv,0:max_deriv,0:max_deriv))
327 
328  ! T_EF
329  allocate(eq%T_EF(dims(1),dims(2),dims(3),9,&
330  &0:max_deriv,0:max_deriv,0:max_deriv))
331 
332  ! T_FE
333  allocate(eq%T_FE(dims(1),dims(2),dims(3),9,&
334  &0:max_deriv,0:max_deriv,0:max_deriv))
335 
336  ! det_T_EF
337  allocate(eq%det_T_EF(dims(1),dims(2),dims(3),&
338  &0:max_deriv,0:max_deriv,0:max_deriv))
339 
340  ! det_T_FE
341  allocate(eq%det_T_FE(dims(1),dims(2),dims(3),&
342  &0:max_deriv,0:max_deriv,0:max_deriv))
343 
344 #if ldebug
345  ! set estimated memory usage
346  if (print_mem_usage) eq%estim_mem_usage(1) = &
347  &eq%estim_mem_usage(1) + product(dims)*&
348  &((max_deriv+1)**3*(4*6+2*9+4))
349 #endif
350 
351  ! initialize variables that are specificic to which equilibrium
352  ! style is being used:
353  ! 1: VMEC
354  ! 2: HELENA
355  select case (eq_style)
356  case (1) ! VMEC
357  ! g_C
358  allocate(eq%g_C(dims(1),dims(2),dims(3),6,&
359  &0:max_deriv,0:max_deriv,0:max_deriv))
360 
361  ! T_VC
362  allocate(eq%T_VC(dims(1),dims(2),dims(3),9,&
363  &0:max_deriv,0:max_deriv,0:max_deriv))
364 
365  ! R
366  allocate(eq%R_E(dims(1),dims(2),dims(3),&
367  &0:max_deriv+1,0:max_deriv+1,0:max_deriv+1))
368 
369  ! Z
370  allocate(eq%Z_E(dims(1),dims(2),dims(3),&
371  &0:max_deriv+1,0:max_deriv+1,0:max_deriv+1))
372 
373  ! lambda
374  allocate(eq%L_E(dims(1),dims(2),dims(3),&
375  &0:max_deriv+1,0:max_deriv+1,0:max_deriv+1))
376 
377 #if ldebug
378  ! set estimated memory usage
379  if (print_mem_usage) eq%estim_mem_usage(1) = &
380  &eq%estim_mem_usage(1) + product(dims)*(&
381  &(max_deriv+1)**3*(1*6+1*9+2) + &
382  &(max_deriv+2)**3*(3))
383 #endif
384  case (2) ! HELENA
385  ! do nothing
386  end select
387  end if
388 
389  if (setup_f_loc) then
390  ! initialize variables that are used for all equilibrium styles
391  ! g_FD
392  allocate(eq%g_FD(dims(1),dims(2),dims(3),6,&
393  &0:max_deriv,0:max_deriv,0:max_deriv))
394 
395  ! h_FD
396  allocate(eq%h_FD(dims(1),dims(2),dims(3),6,&
397  &0:max_deriv,0:max_deriv,0:max_deriv))
398 
399  ! jac_FD
400  allocate(eq%jac_FD(dims(1),dims(2),dims(3),&
401  &0:max_deriv,0:max_deriv,0:max_deriv))
402 
403  ! magnetic shear
404  allocate(eq%S(dims(1),dims(2),dims(3)))
405 
406  ! normal curvature
407  allocate(eq%kappa_n(dims(1),dims(2),dims(3)))
408 
409  ! geodesic curvature
410  allocate(eq%kappa_g(dims(1),dims(2),dims(3)))
411 
412  ! parallel current
413  allocate(eq%sigma(dims(1),dims(2),dims(3)))
414 
415 #if ldebug
416  ! set estimated memory usage
417  if (print_mem_usage) eq%estim_mem_usage(2) = &
418  &eq%estim_mem_usage(2) + product(dims)*(&
419  &(max_deriv+1)**3*(2*6+1) + &
420  &4)
421 #endif
422  end if
423 
424 #if ldebug
425  ! increment n_alloc_eq_2s
427 
428  ! print memory usage
429  if (print_mem_usage) call writo('[rank '//trim(i2str(rank))//&
430  &' - Expected memory usage of eq_2: '//&
431  &trim(r2strt(eq%estim_mem_usage(1)*0.008))//' kB for E and '//&
432  &trim(r2strt(eq%estim_mem_usage(2)*0.008))//' kB for E]',&
433  &alert=.true.)
434 #endif
435  end subroutine init_eq_2
436 
438  subroutine copy_eq_1(eq_i,grid_i,eq_o)
439  use grid_vars, only: grid_type
440 
441  ! input / output
442  class(eq_1_type), intent(in) :: eq_i
443  type(grid_type), intent(in) :: grid_i
444  type(eq_1_type), intent(inout) :: eq_o
445 
446  ! local variables
447  logical :: setup_E
448  logical :: setup_F
449 
450  setup_e = allocated(eq_i%pres_E)
451  setup_f = allocated(eq_i%pres_FD)
452 
453  call eq_o%init(grid_i,setup_e=setup_e,setup_f=setup_f)
454 
455  if (setup_e) then
456  if (allocated(eq_i%pres_E)) eq_o%pres_E = eq_i%pres_E
457  if (allocated(eq_i%q_saf_E)) eq_o%q_saf_E = eq_i%q_saf_E
458  if (allocated(eq_i%rot_t_E)) eq_o%rot_t_E = eq_i%rot_t_E
459  if (allocated(eq_i%flux_p_E)) eq_o%flux_p_E = eq_i%flux_p_E
460  if (allocated(eq_i%flux_t_E)) eq_o%flux_t_E = eq_i%flux_t_E
461  end if
462 
463  if (setup_f) then
464  if (allocated(eq_i%pres_FD)) eq_o%pres_FD = eq_i%pres_FD
465  if (allocated(eq_i%flux_p_FD)) eq_o%flux_p_FD = eq_i%flux_p_FD
466  if (allocated(eq_i%flux_t_FD)) eq_o%flux_t_FD = eq_i%flux_t_FD
467  if (allocated(eq_i%q_saf_FD)) eq_o%q_saf_FD = eq_i%q_saf_FD
468  if (allocated(eq_i%rot_t_FD)) eq_o%rot_t_FD = eq_i%rot_t_FD
469  if (allocated(eq_i%rho)) eq_o%rho = eq_i%rho
470  end if
471  end subroutine copy_eq_1
472 
476  subroutine copy_eq_2(eq_i,grid_i,eq_o)
477  use num_vars, only: eq_style
478  use grid_vars, only: grid_type
479 
480  ! input / output
481  class(eq_2_type), intent(in) :: eq_i
482  type(grid_type), intent(in) :: grid_i
483  type(eq_2_type), intent(inout) :: eq_o
484 
485  ! local variables
486  logical :: setup_E
487  logical :: setup_F
488 
489  setup_e = allocated(eq_i%jac_E)
490  setup_f = allocated(eq_i%jac_FD)
491 
492  call eq_o%init(grid_i,setup_e=setup_e,setup_f=setup_f)
493 
494  if (setup_e) then
495  if (allocated(eq_i%g_E)) eq_o%g_E = eq_i%g_E
496  if (allocated(eq_i%h_E)) eq_o%h_E = eq_i%h_E
497  if (allocated(eq_i%g_F)) eq_o%g_F = eq_i%g_F
498  if (allocated(eq_i%h_F)) eq_o%h_F = eq_i%h_F
499  if (allocated(eq_i%jac_E)) eq_o%jac_E = eq_i%jac_E
500  if (allocated(eq_i%jac_F)) eq_o%jac_F = eq_i%jac_F
501  if (allocated(eq_i%T_EF)) eq_o%T_EF = eq_i%T_EF
502  if (allocated(eq_i%T_FE)) eq_o%T_FE = eq_i%T_FE
503  if (allocated(eq_i%det_T_EF)) eq_o%det_T_EF = eq_i%det_T_EF
504  if (allocated(eq_i%det_T_FE)) eq_o%det_T_FE = eq_i%det_T_FE
505 
506  select case (eq_style)
507  case (1) ! VMEC
508  if (allocated(eq_i%g_C)) eq_o%g_C = eq_i%g_C
509  if (allocated(eq_i%T_VC)) eq_o%T_VC = eq_i%T_VC
510  if (allocated(eq_i%R_E)) eq_o%R_E = eq_i%R_E
511  if (allocated(eq_i%Z_E)) eq_o%Z_E = eq_i%Z_E
512  if (allocated(eq_i%l_E)) eq_o%l_E = eq_i%l_E
513  case (2) ! HELENA
514  ! do nothing
515  end select
516  end if
517 
518  if (setup_f) then
519  if (allocated(eq_i%g_FD)) eq_o%g_FD = eq_i%g_FD
520  if (allocated(eq_i%h_FD)) eq_o%h_FD = eq_i%h_FD
521  if (allocated(eq_i%jac_FD)) eq_o%jac_FD = eq_i%jac_FD
522  if (allocated(eq_i%S)) eq_o%S = eq_i%S
523  if (allocated(eq_i%kappa_n)) eq_o%kappa_n = eq_i%kappa_n
524  if (allocated(eq_i%kappa_g)) eq_o%kappa_g = eq_i%kappa_g
525  if (allocated(eq_i%sigma)) eq_o%sigma = eq_i%sigma
526  end if
527  end subroutine copy_eq_2
528 
530  subroutine dealloc_eq_1(eq)
531 #if ldebug
532  use num_vars, only: rank, print_mem_usage
533 #endif
534 
535  ! input / output
536  class(eq_1_type), intent(inout) :: eq
537 
538 #if ldebug
539  ! local variables
540  integer :: mem_diff ! difference in memory
541  real(dp) :: estim_mem_usage ! estimated memory usage
542 
543  ! memory usage before deallocation
544  if (print_mem_usage) then
545  mem_diff = get_mem_usage()
546  estim_mem_usage = sum(eq%estim_mem_usage)
547  end if
548 #endif
549 
550  ! deallocate allocatable variables
551  call dealloc_eq_1_final(eq)
552 
553 #if ldebug
554  ! decrement n_alloc_eq_1s
556 
557  ! memory usage difference after deallocation
558  if (print_mem_usage) then
559  mem_diff = mem_diff - get_mem_usage()
560  call writo('[Rank '//trim(i2str(rank))//' - liberated '//&
561  &trim(i2str(mem_diff))//'kB deallocating eq_1 ('//&
562  &trim(i2str(nint(100*mem_diff/&
563  &(estim_mem_usage*weight_dp))))//&
564  &'% of estimated)]',alert=.true.)
565  end if
566 #endif
567  contains
568  ! Note: intent(out) automatically deallocates the variable
570  subroutine dealloc_eq_1_final(eq)
571  ! input / output
572  type(eq_1_type), intent(out) :: eq ! equilibrium to be deallocated
573  end subroutine dealloc_eq_1_final
574  end subroutine dealloc_eq_1
575 
577  subroutine dealloc_eq_2(eq)
578 #if ldebug
579  use num_vars, only: rank, print_mem_usage
580 #endif
581  ! input / output
582  class(eq_2_type), intent(inout) :: eq
583 
584 #if ldebug
585  ! local variables
586  integer :: mem_diff ! difference in memory
587  real(dp) :: estim_mem_usage ! estimated memory usage
588 
589  ! memory usage before deallocation
590  if (print_mem_usage) then
591  mem_diff = get_mem_usage()
592  estim_mem_usage = sum(eq%estim_mem_usage)
593  end if
594 #endif
595 
596  ! deallocate allocatable variables
597  call dealloc_eq_2_final(eq)
598 
599 #if ldebug
600  ! decrement n_alloc_eq_2s
602 
603  ! memory usage difference after deallocation
604  if (print_mem_usage) then
605  mem_diff = mem_diff - get_mem_usage()
606  call writo('[Rank '//trim(i2str(rank))//' - liberated '//&
607  &trim(i2str(mem_diff))//'kB deallocating eq_2 ('//&
608  &trim(i2str(nint(100*mem_diff/&
609  &(estim_mem_usage*weight_dp))))//&
610  &'% of estimated)]',alert=.true.)
611  end if
612 #endif
613  contains
614  ! Note: intent(out) automatically deallocates the variable
616  subroutine dealloc_eq_2_final(eq)
617  ! input / output
618  type(eq_2_type), intent(out) :: eq ! equilibrium to be deallocated
619  end subroutine dealloc_eq_2_final
620  end subroutine dealloc_eq_2
621 end module eq_vars
num_vars::mu_0_original
real(dp), parameter, public mu_0_original
permeability of free space
Definition: num_vars.f90:84
eq_vars::init_eq_2
subroutine init_eq_2(eq, grid, setup_E, setup_F)
Initializes new metric equilibrium.
Definition: eq_vars.f90:273
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
eq_vars::vac_perm
real(dp), public vac_perm
either usual mu_0 (default) or normalized
Definition: eq_vars.f90:48
num_vars
Numerical variables used by most other modules.
Definition: num_vars.f90:4
num_vars::max_str_ln
integer, parameter, public max_str_ln
maximum length of strings
Definition: num_vars.f90:50
eq_vars::copy_eq_2
subroutine copy_eq_2(eq_i, grid_i, eq_o)
Deep copy of metric equilibrium variables.
Definition: eq_vars.f90:477
str_utilities::i2str
elemental character(len=max_str_ln) function, public i2str(k)
Convert an integer to string.
Definition: str_utilities.f90:18
messages::get_mem_usage
integer function, public get_mem_usage()
Returns the memory usage in kilobytes.
Definition: messages.f90:554
eq_vars::t_0
real(dp), public t_0
derived normalization constant for nondimensionalization
Definition: eq_vars.f90:47
eq_vars::n_alloc_eq_1s
integer, public n_alloc_eq_1s
nr. of allocated eq_1 variables
Definition: eq_vars.f90:53
str_utilities
Operations on strings.
Definition: str_utilities.f90:4
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
eq_vars::r_0
real(dp), public r_0
independent normalization constant for nondimensionalization
Definition: eq_vars.f90:42
eq_vars::eq_1_type
flux equilibrium type
Definition: eq_vars.f90:63
num_vars::weight_dp
real(dp), parameter, public weight_dp
size of double precision in kB
Definition: num_vars.f90:49
num_vars::print_mem_usage
logical, public print_mem_usage
print memory usage is printed
Definition: num_vars.f90:149
num_vars::eq_style
integer, public eq_style
either 1 (VMEC) or 2 (HELENA)
Definition: num_vars.f90:89
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
eq_vars::rho_0
real(dp), public rho_0
independent normalization constant for nondimensionalization
Definition: eq_vars.f90:44
num_vars::max_deriv
integer, parameter, public max_deriv
highest derivatives for metric factors in Flux coords.
Definition: num_vars.f90:52
messages::writo
subroutine, public writo(input_str, persistent, error, warning, alert)
Write output to file identified by output_i.
Definition: messages.f90:275
messages
Numerical utilities related to giving output.
Definition: messages.f90:4
num_vars::pi
real(dp), parameter, public pi
Definition: num_vars.f90:83
grid_vars
Variables pertaining to the different grids used.
Definition: grid_vars.f90:4
eq_vars::dealloc_eq_2
subroutine dealloc_eq_2(eq)
Deallocates metric equilibrium quantities.
Definition: eq_vars.f90:578
eq_vars::copy_eq_1
subroutine copy_eq_1(eq_i, grid_i, eq_o)
Deep copy of flux equilibrium variables.
Definition: eq_vars.f90:439
eq_vars::psi_0
real(dp), public psi_0
derived normalization constant for nondimensionalization
Definition: eq_vars.f90:46
eq_vars::dealloc_eq_1
subroutine dealloc_eq_1(eq)
Deallocates flux equilibrium quantities.
Definition: eq_vars.f90:531
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
eq_vars::init_eq_1
subroutine init_eq_1(eq, grid, setup_E, setup_F)
Initializes new flux equilibrium.
Definition: eq_vars.f90:174
eq_vars::n_alloc_eq_2s
integer, public n_alloc_eq_2s
nr. of allocated eq_2 variables
Definition: eq_vars.f90:55
eq_vars::eq_2_type
metric equilibrium type
Definition: eq_vars.f90:114
eq_vars::b_0
real(dp), public b_0
derived normalization constant for nondimensionalization
Definition: eq_vars.f90:45