forked from landreman/regcoil_pm
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathregcoil_evaluate_coil_surface.f90
312 lines (279 loc) · 16.4 KB
/
regcoil_evaluate_coil_surface.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
subroutine regcoil_evaluate_coil_surface()
! This subroutine takes the arrays rmnc_coil, zmns_coil, etc, and evaluates the position vector r_coil
! and its first and second derivatives with respect to theta and zeta.
use regcoil_Gaussian_quadrature
use regcoil_variables
use stel_kinds
implicit none
integer :: imn, m, n, iflag, j, js
real(dp) :: rmnc, rmns, zmnc, zmns
real(dp) :: angle, sinangle, cosangle, dsinangledtheta, dcosangledtheta, dsinangledzeta, dcosangledzeta
real(dp) :: d2sinangledtheta2, d2sinangledthetadzeta, d2sinangledzeta2, d2cosangledtheta2, d2cosangledthetadzeta, d2cosangledzeta2
real(dp) :: angle2, sinangle2, cosangle2, dsinangle2dzeta, dcosangle2dzeta, d2sinangle2dzeta2, d2cosangle2dzeta2
real(dp), dimension(:,:), allocatable :: major_R_squared
integer :: itheta, izeta
integer :: tic, toc, countrate
real(dp), dimension(:,:), allocatable :: sin_m_theta, cos_m_theta, sin_n_zeta, cos_n_zeta
real(dp), dimension(:,:), allocatable :: fundamental_form_E, fundamental_form_F, fundamental_form_G
real(dp), dimension(:,:), allocatable :: fundamental_form_L, fundamental_form_M, fundamental_form_P
real(dp), dimension(:,:), allocatable :: temp_matrix
real(dp) :: a, temp, d1, d2
call system_clock(tic,countrate)
allocate(sin_m_theta(mnmax_coil,ntheta_coil))
allocate(cos_m_theta(mnmax_coil,ntheta_coil))
allocate(sin_n_zeta( mnmax_coil,nzetal_coil))
allocate(cos_n_zeta( mnmax_coil,nzetal_coil))
do itheta = 1,ntheta_coil
sin_m_theta(:,itheta) = sin(xm_coil * theta_coil(itheta))
cos_m_theta(:,itheta) = cos(xm_coil * theta_coil(itheta))
end do
do izeta = 1,nzetal_coil
sin_n_zeta(:, izeta) = sin(xn_coil * zetal_coil(izeta))
cos_n_zeta(:, izeta) = cos(xn_coil * zetal_coil(izeta))
end do
r_coil = 0
drdtheta_coil = 0
drdzeta_coil = 0
d2rdtheta2_coil = 0
d2rdthetadzeta_coil = 0
d2rdzeta2_coil = 0
do izeta = 1,nzetal_coil
angle2 = zetal_coil(izeta)
sinangle2 = sin(angle2)
cosangle2 = cos(angle2)
dsinangle2dzeta = cosangle2
dcosangle2dzeta = -sinangle2
d2sinangle2dzeta2 = -sinangle2
d2cosangle2dzeta2 = -cosangle2
do imn = 1, mnmax_coil
m = xm_coil(imn)
n = xn_coil(imn)
rmnc = rmnc_coil(imn)
rmns = rmns_coil(imn)
zmnc = zmnc_coil(imn)
zmns = zmns_coil(imn)
do itheta = 1,ntheta_coil
!angle = m*theta_coil(itheta) - n*zetal_coil(izeta)
!sinangle = sin(angle)
!cosangle = cos(angle)
! Trig angle sum formulae for angle = m*theta_coil(itheta) - n*zetal_coil(izeta):
sinangle = sin_m_theta(imn,itheta) * cos_n_zeta(imn,izeta) - cos_m_theta(imn,itheta) * sin_n_zeta(imn,izeta)
cosangle = cos_m_theta(imn,itheta) * cos_n_zeta(imn,izeta) + sin_m_theta(imn,itheta) * sin_n_zeta(imn,izeta)
!if (abs(sinangle - sin(angle)) > 1d-10) stop "Error sin"
!if (abs(cosangle - cos(angle)) > 1d-10) stop "Error cos"
dsinangledtheta = cosangle*m
dcosangledtheta = -sinangle*m
dsinangledzeta = -cosangle*n
dcosangledzeta = sinangle*n
r_coil(1,itheta,izeta) = r_coil(1,itheta,izeta) + rmnc * cosangle * cosangle2 + rmns * sinangle * cosangle2
r_coil(2,itheta,izeta) = r_coil(2,itheta,izeta) + rmnc * cosangle * sinangle2 + rmns * sinangle * sinangle2
r_coil(3,itheta,izeta) = r_coil(3,itheta,izeta) + zmns * sinangle + zmnc * cosangle
drdtheta_coil(1,itheta,izeta) = drdtheta_coil(1,itheta,izeta) + rmnc * dcosangledtheta * cosangle2 + rmns * dsinangledtheta * cosangle2
drdtheta_coil(2,itheta,izeta) = drdtheta_coil(2,itheta,izeta) + rmnc * dcosangledtheta * sinangle2 + rmns * dsinangledtheta * sinangle2
drdtheta_coil(3,itheta,izeta) = drdtheta_coil(3,itheta,izeta) + zmns * dsinangledtheta + zmnc * dcosangledtheta
drdzeta_coil(1,itheta,izeta) = drdzeta_coil(1,itheta,izeta) + rmnc * (dcosangledzeta * cosangle2 + cosangle * dcosangle2dzeta) &
+ rmns * (dsinangledzeta * cosangle2 + sinangle * dcosangle2dzeta)
drdzeta_coil(2,itheta,izeta) = drdzeta_coil(2,itheta,izeta) + rmnc * (dcosangledzeta * sinangle2 + cosangle * dsinangle2dzeta) &
+ rmns * (dsinangledzeta * sinangle2 + sinangle * dsinangle2dzeta)
drdzeta_coil(3,itheta,izeta) = drdzeta_coil(3,itheta,izeta) + zmns * dsinangledzeta + zmnc * dcosangledzeta
! 2nd derivatives are only constructed on 1 field period.
if (izeta > nzeta_coil) cycle
d2sinangledtheta2 = -m*m*sinangle
d2sinangledthetadzeta = m*n*sinangle
d2sinangledzeta2 = -n*n*sinangle
d2cosangledtheta2 = -m*m*cosangle
d2cosangledthetadzeta = m*n*cosangle
d2cosangledzeta2 = -n*n*cosangle
d2rdtheta2_coil(1,itheta,izeta) = d2rdtheta2_coil(1,itheta,izeta) + rmnc * d2cosangledtheta2 * cosangle2 + rmns * d2sinangledtheta2 * cosangle2
d2rdtheta2_coil(2,itheta,izeta) = d2rdtheta2_coil(2,itheta,izeta) + rmnc * d2cosangledtheta2 * sinangle2 + rmns * d2sinangledtheta2 * sinangle2
d2rdtheta2_coil(3,itheta,izeta) = d2rdtheta2_coil(3,itheta,izeta) + zmns * d2sinangledtheta2 + zmnc * d2cosangledtheta2
d2rdthetadzeta_coil(1,itheta,izeta) = d2rdthetadzeta_coil(1,itheta,izeta) + rmnc * (d2cosangledthetadzeta * cosangle2 + dcosangledtheta * dcosangle2dzeta) &
+ rmns * (d2sinangledthetadzeta * cosangle2 + dsinangledtheta * dcosangle2dzeta)
d2rdthetadzeta_coil(2,itheta,izeta) = d2rdthetadzeta_coil(2,itheta,izeta) + rmnc * (d2cosangledthetadzeta * sinangle2 + dcosangledtheta * dsinangle2dzeta) &
+ rmns * (d2sinangledthetadzeta * sinangle2 + dsinangledtheta * dsinangle2dzeta)
d2rdthetadzeta_coil(3,itheta,izeta) = d2rdthetadzeta_coil(3,itheta,izeta) + zmns * d2sinangledthetadzeta + zmnc * d2cosangledthetadzeta
d2rdzeta2_coil(1,itheta,izeta) = d2rdzeta2_coil(1,itheta,izeta) + rmnc * (d2cosangledzeta2 * cosangle2 + dcosangledzeta * dcosangle2dzeta &
+ dcosangledzeta * dcosangle2dzeta + cosangle * d2cosangle2dzeta2) &
+ rmns * (d2sinangledzeta2 * cosangle2 + dsinangledzeta * dcosangle2dzeta &
+ dsinangledzeta * dcosangle2dzeta + sinangle * d2cosangle2dzeta2)
d2rdzeta2_coil(2,itheta,izeta) = d2rdzeta2_coil(2,itheta,izeta) + rmnc * (d2cosangledzeta2 * sinangle2 + dcosangledzeta * dsinangle2dzeta &
+ dcosangledzeta * dsinangle2dzeta + cosangle * d2sinangle2dzeta2) &
+ rmns * (d2sinangledzeta2 * sinangle2 + dsinangledzeta * dsinangle2dzeta &
+ dsinangledzeta * dsinangle2dzeta + sinangle * d2sinangle2dzeta2)
d2rdzeta2_coil(3,itheta,izeta) = d2rdzeta2_coil(3,itheta,izeta) + zmns * d2sinangledzeta2 + zmnc * d2cosangledzeta2
end do
end do
end do
deallocate(sin_m_theta, cos_m_theta, sin_n_zeta, cos_n_zeta)
call system_clock(toc)
if (verbose) print *," Evaluating coil surface & derivatives:",real(toc-tic)/countrate," sec."
!!$ print *,"mnmax_coil:",mnmax_coil
!!$ print *,"xm_coil:",xm_coil
!!$ print *,"xn_coil:",xn_coil
!!$ print *,"rmnc_coil:",rmnc_coil
!!$ print *,"zmns_coil:",zmns_coil
!!$ print *,"r_coil(1,:,:):"
!!$ do j = 1,ntheta_coil
!!$ print *,r_coil(1,j,:)
!!$ end do
!!$ print *,"r_coil(2,:,:):"
!!$ do j = 1,ntheta_coil
!!$ print *,r_coil(2,j,:)
!!$ end do
!!$ print *,"r_coil(3,:,:):"
!!$ do j = 1,ntheta_coil
!!$ print *,r_coil(3,j,:)
!!$ end do
!!$ print *,"d2xdtheta2:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdtheta2_coil(1,j,:)
!!$ end do
!!$ print *,"d2xdthetadzeta:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdthetadzeta_coil(1,j,:)
!!$ end do
!!$ print *,"d2xdzeta2:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdzeta2_coil(1,j,:)
!!$ end do
!!$
!!$ print *,"d2ydtheta2:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdtheta2_coil(2,j,:)
!!$ end do
!!$ print *,"d2ydthetadzeta:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdthetadzeta_coil(2,j,:)
!!$ end do
!!$ print *,"d2ydzeta2:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdzeta2_coil(2,j,:)
!!$ end do
!!$
!!$ print *,"d2zdtheta2:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdtheta2_coil(3,j,:)
!!$ end do
!!$ print *,"d2zdthetadzeta:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdthetadzeta_coil(3,j,:)
!!$ end do
!!$ print *,"d2zdzeta2:"
!!$ do j = 1,ntheta_coil
!!$ print *,d2rdzeta2_coil(3,j,:)
!!$ end do
! Evaluate cross product:
normal_coil(1,:,:) = drdzeta_coil(2,:,:) * drdtheta_coil(3,:,:) - drdtheta_coil(2,:,:) * drdzeta_coil(3,:,:)
normal_coil(2,:,:) = drdzeta_coil(3,:,:) * drdtheta_coil(1,:,:) - drdtheta_coil(3,:,:) * drdzeta_coil(1,:,:)
normal_coil(3,:,:) = drdzeta_coil(1,:,:) * drdtheta_coil(2,:,:) - drdtheta_coil(1,:,:) * drdzeta_coil(2,:,:)
if (allocated(norm_normal_coil)) deallocate(norm_normal_coil)
allocate(norm_normal_coil(ntheta_coil, nzeta_coil),stat=iflag)
if (iflag .ne. 0) stop 'Allocation error! regcoil_init_coil_surface 11'
norm_normal_coil = sqrt(normal_coil(1,:,1:nzeta_coil)**2 + normal_coil(2,:,1:nzeta_coil)**2 &
+ normal_coil(3,:,1:nzeta_coil)**2)
area_coil = nfp * dtheta_coil * dzeta_coil * sum(norm_normal_coil)
! Evaluate coefficients of the first and second fundamental forms:
allocate(fundamental_form_E(ntheta_coil, nzeta_coil))
allocate(fundamental_form_F(ntheta_coil, nzeta_coil))
allocate(fundamental_form_G(ntheta_coil, nzeta_coil))
allocate(fundamental_form_L(ntheta_coil, nzeta_coil))
allocate(fundamental_form_M(ntheta_coil, nzeta_coil))
allocate(fundamental_form_P(ntheta_coil, nzeta_coil))
fundamental_form_E = drdtheta_coil(1,:,1:nzeta_coil) * drdtheta_coil(1,:,1:nzeta_coil) + drdtheta_coil(2,:,1:nzeta_coil) * drdtheta_coil(2,:,1:nzeta_coil) + drdtheta_coil(3,:,1:nzeta_coil) * drdtheta_coil(3,:,1:nzeta_coil)
fundamental_form_F = drdtheta_coil(1,:,1:nzeta_coil) * drdzeta_coil(1,:,1:nzeta_coil) + drdtheta_coil(2,:,1:nzeta_coil) * drdzeta_coil(2,:,1:nzeta_coil) + drdtheta_coil(3,:,1:nzeta_coil) * drdzeta_coil(3,:,1:nzeta_coil)
fundamental_form_G = drdzeta_coil(1,:,1:nzeta_coil) * drdzeta_coil(1,:,1:nzeta_coil) + drdzeta_coil(2,:,1:nzeta_coil) * drdzeta_coil(2,:,1:nzeta_coil) + drdzeta_coil(3,:,1:nzeta_coil) * drdzeta_coil(3,:,1:nzeta_coil)
fundamental_form_L = (normal_coil(1,:,1:nzeta_coil) * d2rdtheta2_coil(1,:,:) + normal_coil(2,:,1:nzeta_coil) * d2rdtheta2_coil(2,:,:) + normal_coil(3,:,1:nzeta_coil) * d2rdtheta2_coil(3,:,:)) / norm_normal_coil
fundamental_form_M = (normal_coil(1,:,1:nzeta_coil) * d2rdthetadzeta_coil(1,:,:) + normal_coil(2,:,1:nzeta_coil) * d2rdthetadzeta_coil(2,:,:) + normal_coil(3,:,1:nzeta_coil) * d2rdthetadzeta_coil(3,:,:)) / norm_normal_coil
fundamental_form_P = (normal_coil(1,:,1:nzeta_coil) * d2rdzeta2_coil(1,:,:) + normal_coil(2,:,1:nzeta_coil) * d2rdzeta2_coil(2,:,:) + normal_coil(3,:,1:nzeta_coil) * d2rdzeta2_coil(3,:,:)) / norm_normal_coil
mean_curvature_coil = (fundamental_form_L * fundamental_form_G + fundamental_form_P * fundamental_form_E - 2 * fundamental_form_M * fundamental_form_F) &
/ (2 * (fundamental_form_E * fundamental_form_G - fundamental_form_F * fundamental_form_F))
! Initialize s arrays
allocate(s_integration(ns_integration))
allocate(s_weights(ns_integration))
allocate(s_magnetization(ns_magnetization))
allocate(temp_matrix(ns_integration,ns_integration))
select case (trim(s_integration_option))
case (s_integration_option_Chebyshev)
call regcoil_Chebyshev_grid(ns_integration, 0.0_dp, 1.0_dp, s_integration, s_weights, temp_matrix)
case (s_integration_option_Gaussian)
call get_legendre_grids(0.0_dp, 1.0_dp, s_integration, s_weights)
case default
stop "Unrecognized s_integration_option"
end select
deallocate(temp_matrix)
allocate(temp_matrix(ns_magnetization,ns_magnetization))
allocate(s_magnetization_weights(ns_magnetization))
call regcoil_Chebyshev_grid(ns_magnetization, 0.0_dp, 1.0_dp, s_magnetization, s_magnetization_weights, temp_matrix)
deallocate(temp_matrix)
allocate(interpolate_magnetization_to_integration(ns_integration,ns_magnetization))
call regcoil_Chebyshev_interpolation_matrix(ns_magnetization, ns_integration, s_magnetization, s_integration, interpolate_magnetization_to_integration)
print *,"s_integration:",s_integration
print *,"s_weights:",s_weights
print *,"s_magnetization",s_magnetization
print *,"s_magnetization_weights",s_magnetization_weights
print *,"interpolate_magnetization_to_integration:"
do j = 1,ns_integration
print *,interpolate_magnetization_to_integration(j,:)
end do
! For now, just set d=constant. We'll improve this eventually.
allocate(d(ntheta_coil, nzeta_coil))
d = d_initial
d0 = d_initial
allocate(Jacobian_ssquared_term(ntheta_coil, nzeta_coil))
! We wll need this quantity later to generate the Jacobian of the (s, theta, zeta) coordinates in the magnetization region:
Jacobian_ssquared_term = (fundamental_form_M * fundamental_form_M - fundamental_form_L * fundamental_form_P) / norm_normal_coil
allocate(max_d_before_singularity(ntheta_coil, nzeta_coil))
do izeta = 1, nzeta_coil
do itheta = 1, ntheta_coil
! Use quadratic formula to find when Jacobian=0 for s=1:
a = (fundamental_form_M(itheta,izeta) * fundamental_form_M(itheta,izeta) - fundamental_form_L(itheta,izeta) * fundamental_form_P(itheta,izeta)) &
/ (norm_normal_coil(itheta,izeta) * norm_normal_coil(itheta,izeta))
if (abs(a) < 1e-12) then
! Handle the a=0 case separately
d1 = 1 / (2 * sign_normal * mean_curvature_coil(itheta,izeta))
if (d1 > 0) then
max_d_before_singularity(itheta,izeta) = d1
else
max_d_before_singularity(itheta,izeta) = 1.0d+200
end if
!print *,"a=0:",a," mean_curv:",mean_curvature_coil(itheta,izeta)," temp:",temp," d1:",d1," d2:",d2," max_d:",max_d_before_singularity(itheta,izeta)
cycle
end if
temp = sqrt(4 * mean_curvature_coil(itheta,izeta) * mean_curvature_coil(itheta,izeta) + 4 * a)
d1 = (-2 * sign_normal * mean_curvature_coil(itheta,izeta) + temp) / (2 * a)
d2 = (-2 * sign_normal * mean_curvature_coil(itheta,izeta) - temp) / (2 * a)
if (d1 < 0) then
if (d2 >= 0) then
max_d_before_singularity(itheta,izeta) = d2
else
max_d_before_singularity(itheta,izeta) = 1.0d+200
end if
else
if ((d2 < d1) .and. (d2 >= 0)) then
max_d_before_singularity(itheta,izeta) = d2
else
max_d_before_singularity(itheta,izeta) = d1
end if
end if
!print "(3(a,es10.2))","d1:",d1," d2:",d2," max_d:",max_d_before_singularity(itheta,izeta)
!print *,"a:",a," mean_curv:",mean_curvature_coil(itheta,izeta)," temp:",temp," d1:",d1," d2:",d2," max_d:",max_d_before_singularity(itheta,izeta)
end do
end do
print "(a,es12.4)"," Maximum allowed uniform d before singularity:",minval(max_d_before_singularity)
if (d_initial >= minval(max_d_before_singularity)) stop "Error! d_initial is too large. Jacobian is singular."
deallocate(fundamental_form_E, fundamental_form_F, fundamental_form_G, fundamental_form_L, fundamental_form_M, fundamental_form_P)
! Compute coil surface volume using \int (1/2) R^2 dZ dzeta.
! These quantities will be evaluated on the half theta grid, which is the natural grid for dZ,
! but we will need to interpolate R^2 from the full to half grid.
allocate(major_R_squared(ntheta_coil,nzetal_coil))
major_R_squared = r_coil(1,:,:)*r_coil(1,:,:) + r_coil(2,:,:)*r_coil(2,:,:)
! First handle the interior of the theta grid:
volume_coil = sum((major_R_squared(1:ntheta_coil-1,:) + major_R_squared(2:ntheta_coil,:)) * (0.5d+0) & ! R^2, interpolated from full to half grid
* (r_coil(3,2:ntheta_coil,:)-r_coil(3,1:ntheta_coil-1,:))) ! dZ
! Add the contribution from the ends of the theta grid:
volume_coil = volume_coil + sum((major_R_squared(1,:) + major_R_squared(ntheta_coil,:)) * (0.5d+0) & ! R^2, interpolated from full to half grid
* (r_coil(3,1,:)-r_coil(3,ntheta_coil,:))) ! dZ
volume_coil = abs(volume_coil * dzeta_coil / 2) ! r includes all nfp periods already, so no factor of nfp needed.
deallocate(major_R_squared)
if (verbose) print "(a,es10.3,a,es10.3,a)"," Coil surface area:",area_coil," m^2. Volume:",volume_coil," m^3."
end subroutine regcoil_evaluate_coil_surface