PB3D [2.47]
Ideal linear high-n MHD stability in 3-D
Loading...
Searching...
No Matches
VMEC_ops.f90
Go to the documentation of this file.
1!------------------------------------------------------------------------------!
2!> Operations that concern the output of VMEC.
3!------------------------------------------------------------------------------!
4module vmec_ops
5#include <PB3D_macros.h>
7 use output_ops
8 use messages
9 use num_vars, only: &
11 use read_wout_mod, only: read_wout_file, read_wout_deallocate, & ! from LIBSTELL
12 &n_r_vmec => ns, & ! n_r
13 &xn, xm, & ! xn, xm
14 &lrfp, & ! whether or not the poloidal flux is used as radial variable
15 &phi, phipf, & ! toroidal flux (FM), norm. deriv. of toroidal flux (FM)
16 &iotaf, & ! rot. transf. (tor. flux) or saf. fac. (pol. flux) (FM)
17 &presf, & ! pressure (FM)
18 &bsubumns, bsubumnc, bsubvmns, bsubvmnc, bsubsmns, bsubsmnc, & ! B_theta (HM), B_zeta (HM), B_r (FM)
19 &jcuru, jcurv, & ! poloidal (FM) and toroidal (FM) current
20 &bmns, bmnc, & ! magnitude of B (HM)
21 &lmns, lmnc, rmns, rmnc, zmns, zmnc, & ! lambda (HM), R (FM), Z(FM)
22 &gmnc, gmns ! Jacobian in VMEC coordinates
23 use vmec_vars
24
25 implicit none
26 private
28
29contains
30 !> Reads the VMEC equilibrium data.
31 !!
32 !! \return ierr
33 integer function read_vmec(n_r_in,use_pol_flux_V) result(ierr)
34 use num_utilities, only: calc_int, spline
36
37 character(*), parameter :: rout_name = 'read_VMEC'
38
39 ! input / output
40 integer, intent(inout) :: n_r_in !< nr. of normal points in input grid
41 logical, intent(inout) :: use_pol_flux_v !< .true. if VMEC equilibrium is based on poloidal flux
42
43 ! local variables
44 integer :: id, kd ! counters
45 real(dp), allocatable :: r_v(:) ! normal coordinate
46 real(dp), allocatable :: l_c_h(:,:,:) ! temporary HM variable
47 real(dp), allocatable :: l_s_h(:,:,:) ! temporary HM variable
48 real(dp), allocatable :: jac_c_h(:,:,:) ! temporary HM variable
49 real(dp), allocatable :: jac_s_h(:,:,:) ! temporary HM variable
50 character(len=max_str_ln) :: err_msg ! error message
51 character(len=8) :: flux_name ! either poloidal or toroidal
52 real(dp), allocatable :: b_v_sub_c_m(:,:,:), b_v_sub_s_m(:,:,:) ! Coeff. of B_i in (co)sine series (r,theta,phi) (FM, HM, HM)
53#if ldebug
54 real(dp), allocatable :: b_v_c_h(:,:), b_v_s_h(:,:) ! Coeff. of magnitude of B (HM)
55#endif
56
57 ! initialize ierr
58 ierr = 0
59
60 call writo('Reading data from VMEC output "' &
61 &// trim(eq_name) // '"')
62 call lvl_ud(1)
63
64 ! read VMEC output using LIBSTELL
65 if (scan(eq_name, ".") .ne. scan(eq_name, ".", .true.)) &
66 &call writo('There seems to be a dot "." in "'//trim(eq_name)//&
67 &'". This can create problems.', alert=.true.)
68 call read_wout_file(eq_name,ierr) ! read the VMEC file
69 chckerr('Failed to read the VMEC file')
70
71 ! set some variables
72 n_r_in = n_r_vmec
73 allocate(flux_t_v(n_r_in,0:max_deriv+1))
74 allocate(flux_p_v(n_r_in,0:max_deriv+1))
75 allocate(rot_t_v(n_r_in,0:max_deriv+1))
76 allocate(q_saf_v(n_r_in,0:max_deriv+1))
77 allocate(pres_v(n_r_in,0:max_deriv+1))
78 use_pol_flux_v = lrfp
79 flux_t_v(:,0) = phi
80 flux_t_v(:,1) = phipf
81 flux_p_v(:,1) = iotaf*phipf
82 ierr = calc_int(flux_p_v(:,1),1.0_dp/(n_r_in-1.0_dp),flux_p_v(:,0))
83 chckerr('')
84 pres_v(:,0) = presf
85 rot_t_v(:,0) = iotaf
86 q_saf_v(:,0) = 1._dp/iotaf
87 if (use_pol_flux_v) then
88 flux_name = 'poloidal'
89 else
90 flux_name = 'toroidal'
91 end if
92 b_0_v = sum(1.5*bmnc(:,2)-0.5_dp*bmnc(:,3)) ! convert to full mesh
93
94 call writo('VMEC version is ' // trim(r2str(vmec_version)))
95 if (is_freeb_v) then
96 call writo('Free boundary VMEC')
97 else
98 call writo('Fixed boundary VMEC')
99 end if
100 if (is_asym_v) then
101 call writo('No stellarator symmetry')
102 else
103 call writo('Stellarator symmetry')
104 end if
105 call writo('VMEC has '//trim(i2str(mpol_v))//' poloidal and '&
106 &//trim(i2str(ntor_v))//' toroidal modes,')
107 if (ntor_v.gt.0) call writo('(basic toroidal field period NFP = '//&
108 &trim(i2str(nfp_v))//')')
109 call writo('defined on '//trim(i2str(n_r_in))//' '//flux_name//&
110 &' flux surfaces')
111
112 call writo('Running tests')
113 call lvl_ud(1)
114 if (mnmax_v.ne.((ntor_v+1)+(2*ntor_v+1)*(mpol_v-1))) then ! for mpol_V = 0, only ntor_V+1 values needed
115 err_msg = 'Inconsistency in ntor_V, mpol_V and mnmax_V'
116 ierr = 1
117 chckerr(err_msg)
118 end if
119 call lvl_ud(-1)
120
121 call writo('Updating variables')
122 allocate(mn_v(mnmax_v,2))
123 mn_v(:,1) = nint(xm(1:mnmax_v))
124 mn_v(:,2) = nint(xn(1:mnmax_v))
125
126 call lvl_ud(-1)
127 call writo('Data from VMEC output successfully read')
128
129 call writo('Saving VMEC output to PB3D')
130 call lvl_ud(1)
131
132 ! Allocate the Fourier coefficients
133 allocate(r_v_c(mnmax_v,n_r_in,0:max_deriv+1))
134 allocate(r_v_s(mnmax_v,n_r_in,0:max_deriv+1))
135 allocate(z_v_c(mnmax_v,n_r_in,0:max_deriv+1))
136 allocate(z_v_s(mnmax_v,n_r_in,0:max_deriv+1))
137 allocate(l_v_c(mnmax_v,n_r_in,0:max_deriv+1))
138 allocate(l_v_s(mnmax_v,n_r_in,0:max_deriv+1))
139 allocate(jac_v_c(mnmax_v,n_r_in,0:max_deriv))
140 allocate(jac_v_s(mnmax_v,n_r_in,0:max_deriv))
141 allocate(l_c_h(mnmax_v,n_r_in,0:0))
142 allocate(l_s_h(mnmax_v,n_r_in,0:0))
143 allocate(jac_c_h(mnmax_v,n_r_in,0:0))
144 allocate(jac_s_h(mnmax_v,n_r_in,0:0))
145
146 ! factors R_V_c,s; Z_V_c,s and L_C,s and HM varieties
147 r_v_c(:,:,0) = rmnc(1:mnmax_v,:)
148 z_v_s(:,:,0) = zmns(1:mnmax_v,:)
149 l_s_h(:,:,0) = lmns(1:mnmax_v,:)
150 jac_c_h(:,:,0) = gmnc(1:mnmax_v,:)
151 if (is_asym_v) then ! following only needed in asymmetric situations
152 r_v_s(:,:,0) = rmns(1:mnmax_v,:)
153 z_v_c(:,:,0) = zmnc(1:mnmax_v,:)
154 l_c_h(:,:,0) = lmnc(1:mnmax_v,:)
155 jac_s_h(:,:,0) = gmns(1:mnmax_v,:)
156 else
157 r_v_s = 0._dp
158 z_v_c = 0._dp
159 l_c_h = 0._dp
160 jac_s_h = 0._dp
161 end if
162
163 ! calculate data for normal derivatives with the toroidal (or poloidal)
164 ! flux, normalized wrt. to the maximum flux, equidistantly, so the step
165 ! size is 1/(n_r_in-1).
166 allocate(r_v(n_r_in))
167 r_v = [((kd-1._dp)/(n_r_in-1),kd=1,n_r_in)]
168
169 ! up to order 2, full mesh
170 do kd = 1,2
171 ierr = spline(r_v,q_saf_v(:,0),r_v,q_saf_v(:,kd),&
172 &ord=norm_disc_prec_eq,deriv=kd)
173 chckerr('')
174 ierr = spline(r_v,rot_t_v(:,0),r_v,rot_t_v(:,kd),&
175 &ord=norm_disc_prec_eq,deriv=kd)
176 chckerr('')
177 ierr = spline(r_v,pres_v(:,0),r_v,pres_v(:,kd),&
178 &ord=norm_disc_prec_eq,deriv=kd)
179 chckerr('')
180 ierr = spline(r_v,flux_t_v(:,0),r_v,flux_t_v(:,kd),&
181 &ord=norm_disc_prec_eq,deriv=kd)
182 chckerr('')
183 ierr = spline(r_v,flux_p_v(:,0),r_v,flux_p_v(:,kd),&
184 &ord=norm_disc_prec_eq,deriv=kd)
185 chckerr('')
186 end do
187
188 ! up to order 2, half mesh: also extrapolate
189 do kd = 0,2
190 do id = 1,mnmax_v
191 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
192 &jac_c_h(id,2:n_r_in,0),r_v,jac_v_c(id,:,kd),&
193 &ord=norm_disc_prec_eq,deriv=kd,extrap=.true.)
194 chckerr('')
195 if (is_asym_v) then
196 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
197 &jac_s_h(id,2:n_r_in,0),r_v,jac_v_s(id,:,kd),&
198 &ord=norm_disc_prec_eq,deriv=kd,extrap=.true.)
199 chckerr('')
200 end if
201 end do
202 end do
203
204 ! up to order 2, full mesh
205 do kd = 1,3
206 do id = 1,mnmax_v
207 ! full mesh
208 ierr = spline(r_v,r_v_c(id,:,0),r_v,r_v_c(id,:,kd),&
209 &ord=norm_disc_prec_eq,deriv=kd)
210 chckerr('')
211 ierr = spline(r_v,z_v_s(id,:,0),r_v,z_v_s(id,:,kd),&
212 &ord=norm_disc_prec_eq,deriv=kd)
213 chckerr('')
214 if (is_asym_v) then
215 ierr = spline(r_v,r_v_s(id,:,0),r_v,r_v_s(id,:,kd),&
216 &ord=norm_disc_prec_eq,deriv=kd)
217 chckerr('')
218 ierr = spline(r_v,z_v_c(id,:,0),r_v,z_v_c(id,:,kd),&
219 &ord=norm_disc_prec_eq,deriv=kd)
220 chckerr('')
221 end if
222 end do
223 end do
224
225 ! up to order 2, half mesh: also extrapolate
226 do kd = 0,3
227 do id = 1,mnmax_v
228 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
229 &l_s_h(id,2:n_r_in,0),r_v,l_v_s(id,:,kd),&
230 &ord=norm_disc_prec_eq,deriv=kd,extrap=.true.)
231 chckerr('')
232 if (is_asym_v) then
233 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
234 &l_c_h(id,2:n_r_in,0),r_v,l_v_c(id,:,kd),&
235 &ord=norm_disc_prec_eq,deriv=kd,extrap=.true.)
236 chckerr('')
237 end if
238 end do
239 end do
240
241 ! These are only necessary if the compiler searches for uninitialized
242 ! memory
243 flux_t_v(:,3:) = 0._dp
244 flux_p_v(:,3:) = 0._dp
245 q_saf_v(:,3:) = 0._dp
246 rot_t_v(:,3:) = 0._dp
247 pres_v(:,3:) = 0._dp
248
249 ! allocate helper variables
250 allocate(b_v_sub_c_m(mnmax_v,n_r_in,3)); b_v_sub_c_m = 0._dp
251 allocate(b_v_sub_s_m(mnmax_v,n_r_in,3)); b_v_sub_s_m = 0._dp
252#if ldebug
253 allocate(b_v_c_h(mnmax_v,n_r_in)); b_v_c_h = 0._dp
254 allocate(b_v_s_h(mnmax_v,n_r_in)); b_v_s_h = 0._dp
255#endif
256
257 ! store in helper variables
258 b_v_sub_s_m(:,:,1) = bsubsmns(1:mnmax_v,:)
259 b_v_sub_c_m(:,:,2) = bsubumnc(1:mnmax_v,:)
260 b_v_sub_c_m(:,:,3) = bsubvmnc(1:mnmax_v,:)
261 if (is_asym_v) then ! following only needed in asymmetric situations
262 b_v_sub_c_m(:,:,1) = bsubsmnc(1:mnmax_v,:)
263 b_v_sub_s_m(:,:,2) = bsubumns(1:mnmax_v,:)
264 b_v_sub_s_m(:,:,3) = bsubvmns(1:mnmax_v,:)
265 else
266 b_v_sub_c_m(:,:,1) = 0._dp
267 b_v_sub_s_m(:,:,2) = 0._dp
268 b_v_sub_s_m(:,:,3) = 0._dp
269 end if
270#if ldebug
271 b_v_c_h(:,:) = bmnc(1:mnmax_v,:)
272 if (is_asym_v) then ! following only needed in asymmetric situations
273 b_v_s_h(:,:) = bmns(1:mnmax_v,:)
274 else
275 b_v_s_h(:,:) = 0._dp
276 end if
277#endif
278
279 ! allocate FM variables
280 allocate(b_v_sub_c(mnmax_v,n_r_in,3))
281 allocate(b_v_sub_s(mnmax_v,n_r_in,3))
282#if ldebug
283 allocate(b_v_c(mnmax_v,n_r_in))
284 allocate(b_v_s(mnmax_v,n_r_in))
285 allocate(j_v_sup_int(n_r_in,2))
286#endif
287
288 ! half mesh: extrapolate
289 do id = 1,mnmax_v
290 do kd = 2,3
291 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
292 &b_v_sub_c_m(id,2:n_r_in,kd),r_v,b_v_sub_c(id,:,kd),&
293 &ord=norm_disc_prec_eq,deriv=0,extrap=.true.)
294 chckerr('')
295 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
296 &b_v_sub_s_m(id,2:n_r_in,kd),r_v,b_v_sub_s(id,:,kd),&
297 &ord=norm_disc_prec_eq,deriv=0,extrap=.true.)
298 chckerr('')
299 end do
300 end do
301#if ldebug
302 do id = 1,mnmax_v
303 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
304 &b_v_c_h(id,2:n_r_in),r_v,b_v_c(id,:),&
305 &ord=norm_disc_prec_eq,deriv=0,extrap=.true.)
306 chckerr('')
307 ierr = spline(-0.5_dp/n_r_in+r_v(2:n_r_in),&
308 &b_v_s_h(id,2:n_r_in),r_v,b_v_s(id,:),&
309 &ord=norm_disc_prec_eq,deriv=0,extrap=.true.)
310 chckerr('')
311 end do
312 j_v_sup_int(:,1) = jcuru
313 j_v_sup_int(:,2) = jcurv
314#endif
315
316 ! deallocate repacked variables
317 call read_wout_deallocate()
318
319 call lvl_ud(-1)
320 call writo('Conversion complete')
321 end function read_vmec
322
323 !> Normalizes VMEC input.
324 !!
325 !! \note The normal VMEC coordinate runs from 0 to 1, whatever the
326 !! normalization.
327 subroutine normalize_vmec
328 use eq_vars, only: pres_0, psi_0, r_0, b_0
329
330 ! scale the VMEC quantities
334 r_v_c = r_v_c/r_0
335 r_v_s = r_v_s/r_0
336 z_v_c = z_v_c/r_0
337 z_v_s = z_v_s/r_0
338 l_v_c = l_v_c
339 l_v_s = l_v_s
340 jac_v_c = jac_v_c/(r_0**3)
341 jac_v_s = jac_v_s/(r_0**3)
344#if ldebug
345 b_v_c = b_v_c/b_0
346 b_v_s = b_v_s/b_0
348#endif
349 end subroutine normalize_vmec
350end module vmec_ops
Integrates a function using the trapezoidal rule.
Wrapper to the pspline library, making it easier to use for 1-D applications where speed is not the m...
Variables that have to do with equilibrium quantities and the grid used in the calculations:
Definition eq_vars.f90:27
real(dp), public r_0
independent normalization constant for nondimensionalization
Definition eq_vars.f90:42
real(dp), public psi_0
derived normalization constant for nondimensionalization
Definition eq_vars.f90:46
real(dp), public pres_0
independent normalization constant for nondimensionalization
Definition eq_vars.f90:43
real(dp), public b_0
derived normalization constant for nondimensionalization
Definition eq_vars.f90:45
Numerical utilities related to giving output.
Definition messages.f90:4
subroutine, public lvl_ud(inc)
Increases/decreases lvl of output.
Definition messages.f90:254
subroutine, public writo(input_str, persistent, error, warning, alert)
Write output to file identified by output_i.
Definition messages.f90:275
Numerical utilities.
Numerical variables used by most other modules.
Definition num_vars.f90:4
integer, parameter, public dp
double precision
Definition num_vars.f90:46
real(dp), parameter, public pi
Definition num_vars.f90:83
character(len=max_str_ln), public eq_name
name of equilibrium file from VMEC or HELENA
Definition num_vars.f90:138
integer, parameter, public max_str_ln
maximum length of strings
Definition num_vars.f90:50
integer, parameter, public max_deriv
highest derivatives for metric factors in Flux coords.
Definition num_vars.f90:52
integer, public norm_disc_prec_eq
precision for normal discretization for equilibrium
Definition num_vars.f90:120
Operations concerning giving output, on the screen as well as in output files.
Definition output_ops.f90:5
Operations on strings.
elemental character(len=max_str_ln) function, public i2str(k)
Convert an integer to string.
elemental character(len=max_str_ln) function, public r2str(k)
Convert a real (double) to string.
Operations that concern the output of VMEC.
Definition VMEC_ops.f90:4
integer function, public read_vmec(n_r_in, use_pol_flux_v)
Reads the VMEC equilibrium data.
Definition VMEC_ops.f90:34
subroutine, public normalize_vmec
Normalizes VMEC input.
Definition VMEC_ops.f90:328
Variables that concern the output of VMEC.
Definition VMEC_vars.f90:4
real(dp), dimension(:,:), allocatable, public q_saf_v
safety factor
Definition VMEC_vars.f90:38
real(dp), dimension(:,:,:), allocatable, public jac_v_c
Coeff. of in sine series (HM and FM) and norm. deriv.
Definition VMEC_vars.f90:45
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
integer, dimension(:,:), allocatable, public mn_v
m and n of modes
Definition VMEC_vars.f90:32
real(dp), dimension(:,:,:), allocatable, public l_v_s
Coeff. of in cosine series (HM) and norm. deriv.
Definition VMEC_vars.f90:44
real(dp), dimension(:,:,:), allocatable, public z_v_c
Coeff. of in sine series (FM) and norm. deriv.
Definition VMEC_vars.f90:41
real(dp), dimension(:,:), allocatable, public rot_t_v
rotational transform
Definition VMEC_vars.f90:37
real(dp), public b_0_v
the magnitude of B at the magnetic axis,
Definition VMEC_vars.f90:33
real(dp), dimension(:,:,:), allocatable, public jac_v_s
Coeff. of in cosine series (HM and FM) and norm. deriv.
Definition VMEC_vars.f90:46
real(dp), dimension(:,:,:), allocatable, public r_v_c
Coeff. of in sine series (FM) and norm. deriv.
Definition VMEC_vars.f90:39
real(dp), dimension(:,:), allocatable, public pres_v
pressure
Definition VMEC_vars.f90:36
real(dp), dimension(:,:), allocatable, public j_v_sup_int
Integrated poloidal and toroidal current (FM).
Definition VMEC_vars.f90:52
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
real(dp), dimension(:,:), allocatable, public b_v_s
Coeff. of magnitude of B in cosine series (HM and FM).
Definition VMEC_vars.f90:51
real(dp), dimension(:,:), allocatable, public flux_t_v
toroidal flux
Definition VMEC_vars.f90:34
real(dp), dimension(:,:,:), allocatable, public z_v_s
Coeff. of in cosine series (FM) and norm. deriv.
Definition VMEC_vars.f90:42
real(dp), dimension(:,:,:), allocatable, public r_v_s
Coeff. of in cosine series (FM) and norm. deriv.
Definition VMEC_vars.f90:40
real(dp), dimension(:,:,:), allocatable, public l_v_c
Coeff. of in sine series (HM) and norm. deriv.
Definition VMEC_vars.f90:43
real(dp), dimension(:,:), allocatable, public b_v_c
Coeff. of magnitude of B in sine series (HM and FM).
Definition VMEC_vars.f90:50
real(dp), dimension(:,:), allocatable, public flux_p_v
poloidal flux
Definition VMEC_vars.f90:35