GCC Code Coverage Report


Directory: src/
Coverage: low: ≥ 0% medium: ≥ 75.0% high: ≥ 90.0%
Coverage Exec / Excl / Total
Lines: 95.1% 77 / 0 / 81
Functions: 90.9% 10 / 0 / 11
Branches: 78.1% 50 / 0 / 64

coord_transform/s_coord_transform_toroidal_bspline.f90
Line Branch Exec Source
1 !> @file s_coord_transform_toroidal_bspline.f90
2 !> @brief Submodule implementing ToroidalBSplineTransform methods
3 submodule(m_coord_transform_toroidal_bspline) s_coord_transform_toroidal_bspline
4 implicit none
5
6 contains
7
8 7 module function init_ToroidalBSplineTransform(major_radius, minor_radius, bspline) result(this)
9 use m_bspline_linalg, only: l2_projection
10
11 real(wp), intent(in) :: major_radius, minor_radius
12 type(BSplineSpace), intent(in) :: bspline
13 type(ToroidalBSplineTransform) :: this
14
15 7 this%major_radius = major_radius
16 7 this%minor_radius = minor_radius
17 7 this%is_orthogonal = .false.
18 7 this%has_polar_xy_singularity = .true.
19
20 ! TODO why not allow for clamped B-splines?
21 7 if (.not. bspline%is_periodic) then
22 write (*, '(A)') "ERROR in ToroidalBSplineTransform: B-spline basis must be periodic. " &
23 //"Clamped B-splines are not supported for toroidal B-spline transformations."
24 error stop 1
25 end if
26
27 7 call l2_projection(this%cos_bspline, bspline, cosine)
28 7 call l2_projection(this%sin_bspline, bspline, sine)
29
30 ! NOTE: L2 projection should preserve the integral value, but we enforce zero mean anyways (in agreement with the polar
31 ! NOTE regularity constraint)
32 this%cos_bspline%data = this%cos_bspline%data - sum(this%cos_bspline%data(0:this%cos_bspline%bspline%nr_bsplines - 1)) / &
33
4/4
✓ Branch 10 → 11 taken 88 times.
✓ Branch 10 → 12 taken 7 times.
✓ Branch 13 → 14 taken 103 times.
✓ Branch 13 → 15 taken 7 times.
198 this%cos_bspline%bspline%nr_bsplines
34 this%sin_bspline%data = this%sin_bspline%data - sum(this%sin_bspline%data(0:this%sin_bspline%bspline%nr_bsplines - 1)) / &
35
5/6
✗ Branch 2 → 3 not taken.
✓ Branch 2 → 7 taken 7 times.
✓ Branch 16 → 17 taken 88 times.
✓ Branch 16 → 18 taken 7 times.
✓ Branch 19 → 20 taken 103 times.
✓ Branch 19 → 21 taken 7 times.
205 this%sin_bspline%bspline%nr_bsplines
36 contains
37 292 pure real(wp) function cosine(yp) result(ans)
38 real(wp), intent(in) :: yp
39 292 ans = cos(2 * pi * yp)
40 292 end function cosine
41
42 292 pure real(wp) function sine(yp) result(ans)
43 real(wp), intent(in) :: yp
44 292 ans = sin(2 * pi * yp)
45 292 end function sine
46 end function init_ToroidalBSplineTransform
47
48 1800 module procedure toroidal_bspline_transform
49 !$acc routine seq
50 2400 select case (i)
51 case (1)
52 600 ans = (this%major_radius + this%minor_radius * xp * evaluate(this%sin_bspline, yp)) * cos(2 * pi * zp)
53 case (2)
54 600 ans = (this%major_radius + this%minor_radius * xp * evaluate(this%sin_bspline, yp)) * sin(2 * pi * zp)
55 case (3)
56
3/4
✓ Branch 2 → 3 taken 600 times.
✓ Branch 2 → 4 taken 600 times.
✓ Branch 2 → 5 taken 600 times.
✗ Branch 2 → 6 not taken.
1800 ans = this%minor_radius * xp * evaluate(this%cos_bspline, yp)
57 end select
58 1800 end procedure toroidal_bspline_transform
59
60 module procedure toroidal_bspline_inverse_transform
61 !$acc routine seq
62 ! character(len=256) :: error_msg
63
64 ! write (error_msg, '(A)') "ERROR in ToroidalBSplineTransform::inverse_transform: inverse transform is not implemented. " &
65 ! //"Consider using the analytical ToroidalTransform instead if inverse transformation is required."
66 ! error stop error_msg
67 ans = 0._wp
68 end procedure toroidal_bspline_inverse_transform
69
70 70326964 module procedure toroidal_bspline_jacobian
71 !$acc routine seq
72 ans = xp * 2 * pi * this%minor_radius**2 * toroidal_twopi(this, yp) * (this%major_radius + this%minor_radius * xp *&
73 70326964 & evaluate(this%sin_bspline, yp))
74 70326964 end procedure toroidal_bspline_jacobian
75
76 96595316 module procedure toroidal_twopi
77 !$acc routine seq
78 ans = evaluate(this%cos_bspline, yp) * evaluate(this%sin_bspline, yp, derivative=.true.) &
79 96595316 - evaluate(this%sin_bspline, yp) * evaluate(this%cos_bspline, yp, derivative=.true.)
80 96595316 end procedure toroidal_twopi
81
82 5100 module procedure toroidal_bspline_jacobian_matrix
83 !$acc routine seq
84 6800 select case (i)
85 case (1)
86 4000 select case (j)
87 case (1)
88 600 ans = this%minor_radius * evaluate(this%sin_bspline, yp) * cos(2 * pi * zp)
89 case (2)
90 600 ans = this%minor_radius * xp * evaluate(this%sin_bspline, yp, derivative=.true.) * cos(2 * pi * zp)
91 case (3)
92
3/4
✗ Branch 4 → 3 not taken.
✓ Branch 4 → 5 taken 600 times.
✓ Branch 4 → 6 taken 600 times.
✓ Branch 4 → 7 taken 500 times.
1700 ans = -2 * pi * (this%major_radius + this%minor_radius * xp * evaluate(this%sin_bspline, yp)) * sin(2 * pi * zp)
93 end select
94 case (2)
95 4000 select case (j)
96 case (1)
97 600 ans = this%minor_radius * evaluate(this%sin_bspline, yp) * sin(2 * pi * zp)
98 case (2)
99 600 ans = this%minor_radius * xp * evaluate(this%sin_bspline, yp, derivative=.true.) * sin(2 * pi * zp)
100 case (3)
101
3/4
✗ Branch 8 → 3 not taken.
✓ Branch 8 → 9 taken 600 times.
✓ Branch 8 → 10 taken 600 times.
✓ Branch 8 → 11 taken 500 times.
1700 ans = 2 * pi * (this%major_radius + this%minor_radius * xp * evaluate(this%sin_bspline, yp)) * cos(2 * pi * zp)
102 end select
103 case (3)
104
3/4
✗ Branch 2 → 3 not taken.
✓ Branch 2 → 4 taken 1700 times.
✓ Branch 2 → 8 taken 1700 times.
✓ Branch 2 → 12 taken 1700 times.
5100 select case (j)
105 case (1)
106 600 ans = this%minor_radius * evaluate(this%cos_bspline, yp)
107 case (2)
108 600 ans = this%minor_radius * xp * evaluate(this%cos_bspline, yp, derivative=.true.)
109 case (3)
110
3/4
✗ Branch 12 → 3 not taken.
✓ Branch 12 → 13 taken 600 times.
✓ Branch 12 → 14 taken 600 times.
✓ Branch 12 → 15 taken 500 times.
1700 ans = 0._wp
111 end select
112 end select
113 5100 end procedure toroidal_bspline_jacobian_matrix
114
115 900 module procedure toroidal_bspline_jacobian_matrix_inv
116 !$acc routine seq
117 integer :: k
118 real(wp) :: tmp
119
120 ! We use the identity: J^-1 = G^-1 * J^T
121 ans = 0._wp
122
2/2
✓ Branch 3 → 4 taken 2700 times.
✓ Branch 3 → 7 taken 900 times.
3600 do k = 1, 3
123 2700 tmp = toroidal_bspline_G_matrix_inv(this, i, k, xp, yp, zp)
124
2/2
✓ Branch 4 → 5 taken 1500 times.
✓ Branch 4 → 6 taken 1200 times.
3600 if (tmp /= 0._wp) then
125 1500 ans = ans + tmp * toroidal_bspline_jacobian_matrix(this, j, k, xp, yp, zp)
126 end if
127 end do
128 900 end procedure toroidal_bspline_jacobian_matrix_inv
129
130 50375556 module procedure toroidal_bspline_G_matrix
131 !$acc routine seq
132 69811824 select case (i)
133 case (1)
134 42449884 select case (j)
135 case (1)
136 6221764 ans = this%minor_radius**2 * (evaluate(this%cos_bspline, yp)**2 + evaluate(this%sin_bspline, yp)**2)
137 case (2)
138 ans = this%minor_radius**2 * xp * &
139 & (evaluate(this%cos_bspline, yp) * evaluate(this%cos_bspline, yp, derivative=.true.) + &
140 6607252 & evaluate(this%sin_bspline, yp) * evaluate(this%sin_bspline, yp, derivative=.true.))
141 case (3)
142
3/4
✓ Branch 3 → 4 taken 6221764 times.
✓ Branch 3 → 5 taken 6607252 times.
✓ Branch 3 → 10 taken 6607252 times.
✗ Branch 3 → 12 not taken.
19436268 ans = 0._wp
143 end select
144 case (2)
145 34902124 select case (j)
146 case (1)
147 ans = this%minor_radius**2 * xp * &
148 & (evaluate(this%cos_bspline, yp) * evaluate(this%cos_bspline, yp, derivative=.true.) + &
149 3962836 & evaluate(this%sin_bspline, yp) * evaluate(this%sin_bspline, yp, derivative=.true.))
150 case (2)
151 ans = (this%minor_radius * xp)**2 * (evaluate(this%cos_bspline, yp, derivative=.true.)**2 + &
152 6221764 & evaluate(this%sin_bspline, yp, derivative=.true.)**2)
153 case (3)
154
3/4
✓ Branch 6 → 7 taken 3962836 times.
✓ Branch 6 → 8 taken 6221764 times.
✓ Branch 6 → 10 taken 6607252 times.
✗ Branch 6 → 12 not taken.
16791852 ans = 0._wp
155 end select
156 case (3)
157
3/4
✓ Branch 2 → 3 taken 19436268 times.
✓ Branch 2 → 6 taken 16791852 times.
✓ Branch 2 → 9 taken 14147436 times.
✗ Branch 2 → 12 not taken.
50375556 select case (j)
158 case (1)
159 ans = 0._wp
160 case (2)
161 6221764 ans = 0._wp
162 case (3)
163
2/3
✓ Branch 9 → 10 taken 7925672 times.
✓ Branch 9 → 11 taken 6221764 times.
✗ Branch 9 → 12 not taken.
14147436 ans = (2 * pi * (this%major_radius + this%minor_radius * xp * evaluate(this%sin_bspline, yp)))**2
164 end select
165 end select
166 50375556 end procedure toroidal_bspline_G_matrix
167
168 57825936 module procedure toroidal_bspline_G_matrix_inv
169 !$acc routine seq
170 79745664 select case (i)
171 case (1)
172 48188080 select case (j)
173 case (1)
174 ans = (evaluate(this%cos_bspline, yp, derivative=.true.)**2 + &
175 & evaluate(this%sin_bspline, yp, derivative=.true.)**2) / &
176 6993040 & (toroidal_twopi(this, yp)**2 * this%minor_radius**2)
177 case (2)
178 ans = -(evaluate(this%cos_bspline, yp) * evaluate(this%cos_bspline, yp, derivative=.true.) + evaluate(this%sin_bspline, yp)&
179 7463344 & * evaluate(this%sin_bspline, yp, derivative=.true.)) / (toroidal_twopi(this, yp)**2 * this%minor_radius**2 * xp)
180 case (3)
181
3/4
✓ Branch 3 → 4 taken 6993040 times.
✓ Branch 3 → 5 taken 7463344 times.
✓ Branch 3 → 10 taken 7463344 times.
✗ Branch 3 → 12 not taken.
21919728 ans = 0._wp
182 end select
183 case (2)
184 40725136 select case (j)
185 case (1)
186 ans = -(evaluate(this%cos_bspline, yp) * evaluate(this%cos_bspline, yp, derivative=.true.) + &
187 & evaluate(this%sin_bspline, yp) *&
188 4818928 & evaluate(this%sin_bspline, yp, derivative=.true.)) / (toroidal_twopi(this, yp)**2 * this%minor_radius**2 * xp)
189 case (2)
190 ans = (evaluate(this%cos_bspline, yp)**2 + evaluate(this%sin_bspline, yp)**2) / &
191 6993040 & (toroidal_twopi(this, yp)**2 * this%minor_radius**2 * xp**2)
192 case (3)
193
3/4
✓ Branch 6 → 7 taken 4818928 times.
✓ Branch 6 → 8 taken 6993040 times.
✓ Branch 6 → 10 taken 7463344 times.
✗ Branch 6 → 12 not taken.
19275312 ans = 0._wp
194 end select
195 case (3)
196
3/4
✓ Branch 2 → 3 taken 21919728 times.
✓ Branch 2 → 6 taken 19275312 times.
✓ Branch 2 → 9 taken 16630896 times.
✗ Branch 2 → 12 not taken.
57825936 select case (j)
197 case (1)
198 ans = 0._wp
199 case (2)
200 6993040 ans = 0._wp
201 case (3)
202
2/3
✓ Branch 9 → 10 taken 9637856 times.
✓ Branch 9 → 11 taken 6993040 times.
✗ Branch 9 → 12 not taken.
16630896 ans = 1 / (2 * pi * (this%major_radius + this%minor_radius * xp * evaluate(this%sin_bspline, yp)))**2
203 end select
204 end select
205 57825936 end procedure toroidal_bspline_G_matrix_inv
206
207 end submodule s_coord_transform_toroidal_bspline
208