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
| #include "filaments.h"
double test_vorticity_filament0(coord pcar, int n_seg, double a, double* t0, coord* c, coord* tvec, coord* nvec, coord* bvec, int period){
/* Compute the local coordinates required for the vorticity field.
Each position P(x,y,z) is projected into the local Frenet-Serret frame to
obtain a set of local coordinates, such that:
i) (P - X . T = 0
ii) (P - X) . N = x_n
ii) (P - X) . B = x_b
This requires finding the value of X(t0) along each space curve that verifies
i) through a minization process.
Then, we use the local coordinates (x_n, x_b) to define a radial coordinate
rho required to compute the vorticity of a Lamb-Oseen vortex as
omega = Gamma/(pi a²) exp(-rho²/a²) . T
where Gamma is the circulation and a the core size.
*/
/* First, we approximate the minimal distance between the point P and each
segment of the curve. */
double dmin = 1e30, tmin=0;
double rho_loc;
for (int i = 0; i < n_seg; i++){
if (vecdist2(pcar, c[i]) < dmin){
dmin = vecdist2(pcar, c[i]);
tmin = t0[i];
}
}
if (period != 0)
tmin = fmod(tmin + period*2*pi, period*2*pi);
// If P is close to the vortex, we refine the initial guess
coord ploc, ccar, frenet[3];
double tq = frenet_projection_min(n_seg, a, t0, c, tvec, nvec, bvec, pcar, tmin);
ccar = gsl_interp1d( n_seg, t0, c, tq);
// Then, compute the local coordinates for the vortex
frenet[0] = gsl_interp1d( n_seg, t0, tvec, tq);
frenet[1] = gsl_interp1d( n_seg, t0, nvec, tq);
frenet[2] = gsl_interp1d( n_seg, t0, bvec, tq);
ploc.x = vecdot(vecdiff(pcar, ccar), frenet[0]);
ploc.y = vecdot(vecdiff(pcar, ccar), frenet[1]);
ploc.z = vecdot(vecdiff(pcar, ccar), frenet[2]);
rho_loc = sqrt(vecdot(ploc, ploc));
return rho_loc;
}
coord test_vorticity_filament1(coord pcar, int n_seg, double a, double* t0, coord* c, coord* tvec, coord* nvec, coord* bvec, int period){
double dmin = 1e30, tmin=0;
for (int i = 0; i < n_seg; i++){
if (vecdist2(pcar, c[i]) < dmin){
dmin = vecdist2(pcar, c[i]);
tmin = t0[i];
}
}
if (period != 0)
tmin = fmod(tmin + period*2*pi, period*2*pi);
// If P is close to the vortex, we refine the initial guess
coord ccar, frenet[1];
double tq = frenet_projection_min(n_seg, a, t0, c, tvec, nvec, bvec, pcar, tmin);
ccar = gsl_interp1d( n_seg, t0, c, tq);
// Then, compute the local coordinates for the vortex
frenet[0] = gsl_interp1d( n_seg, t0, tvec, tq);
return (coord) {frenet[0].x, frenet[0].y, frenet[0].z};
}
coord get_vorticity_filament(coord pcar, int n_seg, double a, double* t0, coord* c, coord* tvec, coord* nvec, coord* bvec, int period){
/* Compute the local coordinates required for the vorticity field.
Each position P(x,y,z) is projected into the local Frenet-Serret frame to
obtain a set of local coordinates, such that:
i) (P - X . T = 0
ii) (P - X) . N = x_n
ii) (P - X) . B = x_b
This requires finding the value of X(t0) along each space curve that verifies
i) through a minization process.
Then, we use the local coordinates (x_n, x_b) to define a radial coordinate
rho required to compute the vorticity of a Lamb-Oseen vortex as
omega = Gamma/(pi a²) exp(-rho²/a²) . T
where Gamma is the circulation and a the core size.
*/
/* First, we approximate the minimal distance between the point P and each
segment of the curve. */
coord omega;
double dmin = 1e30, tmin=0;
double rho_loc, omega_mag;
for (int i = 0; i < n_seg; i++){
if (vecdist2(pcar, c[i]) < dmin){
dmin = vecdist2(pcar, c[i]);
tmin = t0[i];
}
}
if (period != 0)
tmin = fmod(tmin + period*2*pi, period*2*pi);
if (dmin < L0){
// If P is close to the vortex, we refine the initial guess
coord ploc, ccar, frenet[3];
double tq = frenet_projection_min(n_seg, a, t0, c, tvec, nvec, bvec, pcar, tmin);
ccar = gsl_interp1d( n_seg, t0, c, tq);
// Then, compute the local coordinates for the vortex
frenet[0] = gsl_interp1d( n_seg, t0, tvec, tq);
frenet[1] = gsl_interp1d( n_seg, t0, nvec, tq);
frenet[2] = gsl_interp1d( n_seg, t0, bvec, tq);
ploc.x = vecdot(vecdiff(pcar, ccar), frenet[0]);
ploc.y = vecdot(vecdiff(pcar, ccar), frenet[1]);
ploc.z = vecdot(vecdiff(pcar, ccar), frenet[2]);
rho_loc = sqrt(vecdot(ploc, ploc));
// Last, we compute the vorticity for the vortex 1
omega_mag = exp(-sq(rho_loc)/sq(a))/(pi*sq(a));
omega = (coord) {omega_mag * frenet[0].x, omega_mag * frenet[0].y, omega_mag * frenet[0].z};
}
else {
// Otherwise, if the point is too far, we set the vorticity to zero.
foreach_dimension()
omega.x = 0;
}
return omega;
}
coord get_vorticity_filament2(double Gamma, double Uc, coord pcar, int n_seg, double a, double* t0, coord* c, coord* tvec, coord* nvec, coord* bvec, int period){
/* Compute the local coordinates required for the vorticity field.
Each position P(x,y,z) is projected into the local Frenet-Serret frame to
obtain a set of local coordinates, such that:
i) (P - X . T = 0
ii) (P - X) . N = x_n
ii) (P - X) . B = x_b
This requires finding the value of X(t0) along each space curve that verifies
i) through a minization process.
Then, we use the local coordinates (x_n, x_b) to define a radial coordinate
rho required to compute the vorticity of a Lamb-Oseen vortex as
omega = Gamma/(pi a²) exp(-rho²/a²) . T + 2 Uc/a (rho/a) exp(-rho²/a²)
where Gamma is the circulation and a the core size and Uc is the
axial centerline velocity.
*/
/* First, we approximate the minimal distance between the point P and each
segment of the curve. */
coord omega;
double dmin = 1e30, tmin=0;
double rho_loc, omega_mag;
for (int i = 0; i < n_seg; i++){
if (vecdist2(pcar, c[i]) < dmin){
dmin = vecdist2(pcar, c[i]);
tmin = t0[i];
}
}
if (period != 0)
tmin = fmod(tmin + period*2*pi, period*2*pi);
if (dmin < L0){
// If P is close to the vortex, we refine the initial guess
coord ploc, ccar, frenet[3];
double tq = frenet_projection_min(n_seg, a, t0, c, tvec, nvec, bvec, pcar, tmin);
ccar = gsl_interp1d( n_seg, t0, c, tq);
// Then, compute the local coordinates for the vortex
frenet[0] = gsl_interp1d( n_seg, t0, tvec, tq);
frenet[1] = gsl_interp1d( n_seg, t0, nvec, tq);
frenet[2] = gsl_interp1d( n_seg, t0, bvec, tq);
ploc.x = vecdot(vecdiff(pcar, ccar), frenet[0]);
ploc.y = vecdot(vecdiff(pcar, ccar), frenet[1]);
ploc.z = vecdot(vecdiff(pcar, ccar), frenet[2]);
rho_loc = sqrt(vecdot(ploc, ploc));
// Last, we compute the vorticity for the vortex 1
omega_mag = exp(-sq(rho_loc)/sq(a))/(sq(a));
omega.x = omega_mag * (Gamma/pi * frenet[0].x + (2.0 * Uc) * (-ploc.z * frenet[1].x + ploc.y * frenet[2].x));
omega.y = omega_mag * (Gamma/pi * frenet[0].y + (2.0 * Uc) * (-ploc.z * frenet[1].y + ploc.y * frenet[2].y));
omega.z = omega_mag * (Gamma/pi * frenet[0].z + (2.0 * Uc) * (-ploc.z * frenet[1].z + ploc.y * frenet[2].z));
}
else {
// Otherwise, if the point is too far, we set the vorticity to zero.
foreach_dimension()
omega.x = 0;
}
return omega;
}
|