-
Notifications
You must be signed in to change notification settings - Fork 1
/
vlm.py
1771 lines (1412 loc) · 64.8 KB
/
vlm.py
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
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""
Define the aerodynamic analysis component using a vortex lattice method.
We input a nodal mesh and properties of the airflow to calculate the
circulations of the horseshoe vortices. We then compute the forces, lift,
and drag acting on the lifting surfaces. Currently we can compute the induced
and viscous drag.
"""
from __future__ import division, print_function
import numpy as np
from openmdao.api import Component, Group, AnalysisError
from scipy.linalg import lu_factor, lu_solve
def view_mat(mat):
""" Helper function used to visually examine matrices. """
import matplotlib.pyplot as plt
if len(mat.shape) > 2:
mat = np.sum(mat, axis=2)
im = plt.imshow(mat.real, interpolation='none')
plt.colorbar(im, orientation='horizontal')
plt.show()
try:
import OAS_API
fortran_flag = True
data_type = float
except:
fortran_flag = False
data_type = complex
def norm(vec):
""" Finds the 2-norm of a vector. """
return np.sqrt(np.sum(vec**2))
def _calc_vorticity(A, B, P):
""" Calculates the influence coefficient for a vortex filament.
Parameters
----------
A[3] : numpy array
Coordinates for the start point of the filament.
B[3] : numpy array
Coordinates for the end point of the filament.
P[3] : numpy array
Coordinates for the collocation point where the influence coefficient
is computed.
Returns
-------
out[3] : numpy array
Influence coefficient contribution for the described filament.
"""
r1 = P - A
r2 = P - B
r1_mag = norm(r1)
r2_mag = norm(r2)
return (r1_mag + r2_mag) * np.cross(r1, r2) / \
(r1_mag * r2_mag * (r1_mag * r2_mag + r1.dot(r2)))
def _assemble_AIC_mtx(mtx, params, surfaces, skip=False):
"""
Compute the aerodynamic influence coefficient matrix
for either solving the linear system or solving for the drag.
We use a nested for loop structure to loop through the lifting surfaces to
obtain the corresponding mesh, then for each mesh we again loop through
the lifting surfaces to obtain the collocation points used to compute
the horseshoe vortex influence coefficients.
This creates mtx with blocks corresponding to each lifting surface's
effects on other lifting surfaces. The block diagonal portions
correspond to each lifting surface's influencen on itself. For a single
lifting surface, this is the entire mtx.
Parameters
----------
mtx[(nx-1)*(ny-1), (nx-1)*(ny-1), 3] : numpy array
Aerodynamic influence coefficient (AIC) matrix, or the
derivative of v w.r.t. circulations.
params : dictionary
OpenMDAO params dictionary for a given aero problem
surfaces : dictionary
Dictionary containing all surfaces in an aero problem.
skip : boolean
If false, the bound vortex contributions on the collocation point
corresponding to the same panel are not included. Used for the drag
computation.
Returns
-------
mtx[tot_panels, tot_panels, 3] : numpy array
Aerodynamic influence coefficient (AIC) matrix, or the
derivative of v w.r.t. circulations.
"""
alpha = params['alpha']
mtx[:, :, :] = 0.0
cosa = np.cos(alpha * np.pi / 180.)
sina = np.sin(alpha * np.pi / 180.)
u = np.array([cosa, 0, sina])
# Set the u-vector to the x-vector to match AVL
# u[0] = 1.
# u[2] = 0.
i_ = 0
i_bpts_ = 0
i_panels_ = 0
# Loop over the lifting surfaces to compute their influence on the flow
# velocity at the collocation points
for surface_ in surfaces:
# Variable names with a trailing underscore correspond to the lifting
# surface being examined, not the collocation point
name_ = surface_['name']
nx_ = surface_['num_x']
ny_ = surface_['num_y']
n_ = nx_ * ny_
n_bpts_ = (nx_ - 1) * ny_
n_panels_ = (nx_ - 1) * (ny_ - 1)
# Obtain the lifting surface mesh in the form expected by the solver,
# with shape [nx_, ny_, 3]
mesh = params[name_+'def_mesh']
bpts = params[name_+'b_pts']
# Set counters to know where to index the sub-matrix within the full mtx
i = 0
i_bpts = 0
i_panels = 0
for surface in surfaces:
# These variables correspond to the collocation points
name = surface['name']
nx = surface['num_x']
ny = surface['num_y']
n = nx * ny
n_bpts = (nx - 1) * ny
n_panels = (nx - 1) * (ny - 1)
symmetry = surface['symmetry']
# Obtain the collocation points used to compute the AIC mtx.
# If setting up the AIC mtx, we use the collocation points (c_pts),
# but if setting up the matrix to solve for drag, we use the
# midpoints of the bound vortices.
if skip:
# Find the midpoints of the bound points, used in drag computations
pts = (params[name+'b_pts'][:, 1:, :] + \
params[name+'b_pts'][:, :-1, :]) / 2
else:
pts = params[name+'c_pts']
# Initialize sub-matrix to populate within full mtx
small_mat = np.zeros((n_panels, n_panels_, 3), dtype=data_type)
# Dense fortran assembly for the AIC matrix
if fortran_flag:
small_mat[:, :, :] = OAS_API.oas_api.assembleaeromtx(alpha, pts, bpts,
mesh, skip, symmetry)
# Python matrix assembly
else:
# Spanwise loop through horseshoe elements
for el_j in range(ny_ - 1):
el_loc_j = el_j * (nx_ - 1)
C_te = mesh[-1, el_j + 1, :]
D_te = mesh[-1, el_j + 0, :]
# Mirror the horseshoe vortex points
if symmetry:
C_te_sym = C_te.copy()
D_te_sym = D_te.copy()
C_te_sym[1] = -C_te_sym[1]
D_te_sym[1] = -D_te_sym[1]
# Spanwise loop through control points
for cp_j in range(ny - 1):
cp_loc_j = cp_j * (nx - 1)
# Chordwise loop through control points
for cp_i in range(nx - 1):
cp_loc = cp_i + cp_loc_j
P = pts[cp_i, cp_j]
r1 = P - D_te
r2 = P - C_te
r1_mag = norm(r1)
r2_mag = norm(r2)
t1 = np.cross(u, r2) / \
(r2_mag * (r2_mag - u.dot(r2)))
t3 = np.cross(u, r1) / \
(r1_mag * (r1_mag - u.dot(r1)))
# AIC contribution from trailing vortex filaments
# coming off the trailing edge
trailing = t1 - t3
# Calculate the effects across the symmetry plane
if symmetry:
r1 = P - D_te_sym
r2 = P - C_te_sym
r1_mag = norm(r1)
r2_mag = norm(r2)
t1 = np.cross(u, r2) / \
(r2_mag * (r2_mag - u.dot(r2)))
t3 = np.cross(u, r1) / \
(r1_mag * (r1_mag - u.dot(r1)))
trailing += t3 - t1
edges = 0
# Chordwise loop through horseshoe elements in
# reversed order, starting with the panel closest
# to the leading edge. This is done to sum the
# AIC contributions from the side vortex filaments
# as we loop through the elements
for el_i in reversed(range(nx_ - 1)):
el_loc = el_i + el_loc_j
A = bpts[el_i, el_j + 0, :]
B = bpts[el_i, el_j + 1, :]
# Check if this is the last panel; if so, use
# the trailing edge mesh points for C & D, else
# use the directly aft panel's bound points
# for C & D
if el_i == nx_ - 2:
C = mesh[-1, el_j + 1, :]
D = mesh[-1, el_j + 0, :]
else:
C = bpts[el_i + 1, el_j + 1, :]
D = bpts[el_i + 1, el_j + 0, :]
# Calculate and store the contributions from
# the vortex filaments on the sides of the
# panels, adding as we progress through the
# panels
edges += _calc_vorticity(B, C, P)
edges += _calc_vorticity(D, A, P)
# Mirror the horseshoe vortex points and
# calculate the effects across
# the symmetry plane
if symmetry:
A_sym = A.copy()
B_sym = B.copy()
C_sym = C.copy()
D_sym = D.copy()
A_sym[1] = -A_sym[1]
B_sym[1] = -B_sym[1]
C_sym[1] = -C_sym[1]
D_sym[1] = -D_sym[1]
edges += _calc_vorticity(C_sym, B_sym, P)
edges += _calc_vorticity(A_sym, D_sym, P)
# If skip, do not include the contributions
# from the panel's bound vortex filament, as
# this causes a singularity when we're taking
# the influence of a panel on its own
# collocation point. This true for the drag
# computation and false for circulation
# computation, due to the different collocation
# points.
if skip and el_loc == cp_loc:
if symmetry:
bound = _calc_vorticity(B_sym, A_sym, P)
else:
bound = np.zeros((3))
small_mat[cp_loc, el_loc, :] = \
trailing + edges + bound
else:
bound = _calc_vorticity(A, B, P)
# Account for symmetry by including the
# mirrored bound vortex
if symmetry:
bound += _calc_vorticity(B_sym, A_sym, P)
small_mat[cp_loc, el_loc, :] = \
trailing + edges + bound
# Populate the full-size matrix with these surface-surface AICs
mtx[i_panels:i_panels+n_panels,
i_panels_:i_panels_+n_panels_, :] = small_mat
i += n
i_bpts += n_bpts
i_panels += n_panels
i_ += n_
i_bpts_ += n_bpts_
i_panels_ += n_panels_
mtx /= 4 * np.pi
def _assemble_AIC_mtx_d(mtxd, params, dparams, dunknowns, dresids, surfaces, skip=False):
"""
Differentiated code to get the forward mode seeds for the AIC matrix assembly.
"""
alpha = params['alpha']
alphad = dparams['alpha']
i_ = 0
i_bpts_ = 0
i_panels_ = 0
# Loop over the lifting surfaces to compute their influence on the flow
# velocity at the collocation points
for surface_ in surfaces:
# Variable names with a trailing underscore correspond to the lifting
# surface being examined, not the collocation point
name_ = surface_['name']
nx_ = surface_['num_x']
ny_ = surface_['num_y']
n_ = nx_ * ny_
n_bpts_ = (nx_ - 1) * ny_
n_panels_ = (nx_ - 1) * (ny_ - 1)
# Obtain the lifting surface mesh in the form expected by the solver,
# with shape [nx_, ny_, 3]
mesh = params[name_+'def_mesh']
bpts = params[name_+'b_pts']
meshd = dparams[name_+'def_mesh']
bptsd = dparams[name_+'b_pts']
# Set counters to know where to index the sub-matrix within the full mtx
i = 0
i_bpts = 0
i_panels = 0
for surface in surfaces:
# These variables correspond to the collocation points
name = surface['name']
nx = surface['num_x']
ny = surface['num_y']
n = nx * ny
n_bpts = (nx - 1) * ny
n_panels = (nx - 1) * (ny - 1)
symmetry = surface['symmetry']
# Obtain the collocation points used to compute the AIC mtx.
# If setting up the AIC mtx, we use the collocation points (c_pts),
# but if setting up the matrix to solve for drag, we use the
# midpoints of the bound vortices.
if skip:
# Find the midpoints of the bound points, used in drag computations
pts = (params[name+'b_pts'][:, 1:, :] + \
params[name+'b_pts'][:, :-1, :]) / 2
ptsd = (dparams[name+'b_pts'][:, 1:, :] + \
dparams[name+'b_pts'][:, :-1, :]) / 2
else:
pts = params[name+'c_pts']
ptsd = dparams[name+'c_pts']
_, small_mat = OAS_API.oas_api.assembleaeromtx_d(alpha, alphad, pts, ptsd,
bpts, bptsd, mesh, meshd,
skip, symmetry)
# Populate the full-size matrix with these surface-surface AICs
mtxd[i_panels:i_panels+n_panels,
i_panels_:i_panels_+n_panels_, :] = small_mat
i += n
i_bpts += n_bpts
i_panels += n_panels
i_ += n_
i_bpts_ += n_bpts_
i_panels_ += n_panels_
mtxd /= 4 * np.pi
def _assemble_AIC_mtx_b(mtxb, params, dparams, dunknowns, dresids, surfaces, skip=False):
"""
Differentiated code to get the reverse mode seeds for the AIC matrix assembly.
"""
alpha = params['alpha']
mtxb /= 4 * np.pi
i_ = 0
i_bpts_ = 0
i_panels_ = 0
# Loop over the lifting surfaces to compute their influence on the flow
# velocity at the collocation points
for surface_ in surfaces:
# Variable names with a trailing underscore correspond to the lifting
# surface being examined, not the collocation point
name_ = surface_['name']
nx_ = surface_['num_x']
ny_ = surface_['num_y']
n_ = nx_ * ny_
n_bpts_ = (nx_ - 1) * ny_
n_panels_ = (nx_ - 1) * (ny_ - 1)
# Obtain the lifting surface mesh in the form expected by the solver,
# with shape [nx_, ny_, 3]
mesh = params[name_+'def_mesh']
bpts = params[name_+'b_pts']
# Set counters to know where to index the sub-matrix within the full mtx
i = 0
i_bpts = 0
i_panels = 0
for surface in surfaces:
# These variables correspond to the collocation points
name = surface['name']
nx = surface['num_x']
ny = surface['num_y']
n = nx * ny
n_bpts = (nx - 1) * ny
n_panels = (nx - 1) * (ny - 1)
symmetry = surface['symmetry']
# Obtain the collocation points used to compute the AIC mtx.
# If setting up the AIC mtx, we use the collocation points (c_pts),
# but if setting up the matrix to solve for drag, we use the
# midpoints of the bound vortices.
if skip:
# Find the midpoints of the bound points, used in drag computations
pts = (params[name+'b_pts'][:, 1:, :] + \
params[name+'b_pts'][:, :-1, :]) / 2
else:
pts = params[name+'c_pts']
small_mtxb = mtxb[i_panels:i_panels+n_panels, i_panels_:i_panels_+n_panels_, :]
alphab, ptsb, bptsb, meshb, mtx = OAS_API.oas_api.assembleaeromtx_b(alpha, pts, bpts,
mesh, skip, symmetry, small_mtxb)
dparams[name_+'def_mesh'] += meshb.real
dparams[name_+'b_pts'] += bptsb.real
if skip:
dparams[name+'b_pts'][:, 1:, :] += ptsb.real / 2
dparams[name+'b_pts'][:, :-1, :] += ptsb.real / 2
else:
dparams[name+'c_pts'] += ptsb.real
dparams['alpha'] += alphab
i += n
i_bpts += n_bpts
i_panels += n_panels
i_ += n_
i_bpts_ += n_bpts_
i_panels_ += n_panels_
class VLMGeometry(Component):
""" Compute various geometric properties for VLM analysis.
Parameters
----------
def_mesh[nx, ny, 3] : numpy array
Array defining the nodal coordinates of the lifting surface.
Returns
-------
b_pts[nx-1, ny, 3] : numpy array
Bound points for the horseshoe vortices, found along the 1/4 chord.
c_pts[nx-1, ny-1, 3] : numpy array
Collocation points on the 3/4 chord line where the flow tangency
condition is satisfed. Used to set up the linear system.
widths[ny-1] : numpy array
The spanwise widths of each individual panel.
lengths[ny] : numpy array
The chordwise length of the entire airfoil following the camber line.
chords[ny] : numpy array
The chordwise distance between the leading and trailing edges.
normals[nx-1, ny-1, 3] : numpy array
The normal vector for each panel, computed as the cross of the two
diagonals from the mesh points.
S_ref : float
The reference area of the lifting surface.
"""
def __init__(self, surface):
super(VLMGeometry, self).__init__()
self.surface = surface
self.ny = surface['num_y']
self.nx = surface['num_x']
# self.fem_origin = surface['fem_origin']
# self.fem_origin = (surface['data_x_upper'][0] + surface['data_x_upper'][-1]) / 2.
self.fem_origin = (surface['data_x_upper'][0] *(surface['data_y_upper'][0]-surface['data_y_lower'][0]) + \
surface['data_x_upper'][-1]*(surface['data_y_upper'][-1]-surface['data_y_lower'][-1])) / \
( (surface['data_y_upper'][0]-surface['data_y_lower'][0]) + (surface['data_y_upper'][-1]-surface['data_y_lower'][-1]))
self.add_param('def_mesh', val=np.zeros((self.nx, self.ny, 3),
dtype=data_type))
self.add_output('b_pts', val=np.zeros((self.nx-1, self.ny, 3),
dtype=data_type))
self.add_output('c_pts', val=np.zeros((self.nx-1, self.ny-1, 3)))
self.add_output('widths', val=np.zeros((self.ny-1)))
self.add_output('cos_sweep', val=np.zeros((self.ny-1)))
self.add_output('lengths', val=np.zeros((self.ny)))
self.add_output('chords', val=np.zeros((self.ny)))
self.add_output('normals', val=np.zeros((self.nx-1, self.ny-1, 3)))
self.add_output('S_ref', val=0.)
def solve_nonlinear(self, params, unknowns, resids):
mesh = params['def_mesh']
# Compute the bound points at 1/4 chord
b_pts = mesh[:-1, :, :] * .75 + mesh[1:, :, :] * .25
# Compute the collocation points at the midpoints of each
# panel's 3/4 chord line
c_pts = 0.5 * 0.25 * mesh[:-1, :-1, :] + \
0.5 * 0.75 * mesh[1:, :-1, :] + \
0.5 * 0.25 * mesh[:-1, 1:, :] + \
0.5 * 0.75 * mesh[1:, 1:, :]
# Compute the widths of each panel at the quarter-chord line
quarter_chord = 0.25 * mesh[-1] + 0.75 * mesh[0]
widths = np.linalg.norm(quarter_chord[1:, :] - quarter_chord[:-1, :], axis=1)
# Compute the numerator of the cosine of the sweep angle of each panel
# (we need this for the viscous drag dependence on sweep, and we only compute
# the numerator because the denominator of the cosine fraction is the width,
# which we have already computed. They are combined in the viscous drag
# calculation.)
cos_sweep = np.linalg.norm(quarter_chord[1:, [1,2]] - quarter_chord[:-1, [1,2]], axis=1)
# Compute the length of each chordwise set of mesh points through the camber line.
dx = mesh[1:, :, 0] - mesh[:-1, :, 0]
dy = mesh[1:, :, 1] - mesh[:-1, :, 1]
dz = mesh[1:, :, 2] - mesh[:-1, :, 2]
lengths = np.sum(np.sqrt(dx**2 + dy**2 + dz**2), axis=0)
# Compute the normal of each panel by taking the cross-product of
# its diagonals. Note that this could be a nonplanar surface
normals = np.cross(
mesh[:-1, 1:, :] - mesh[1:, :-1, :],
mesh[:-1, :-1, :] - mesh[1:, 1:, :],
axis=2)
# Normalize the normal vectors
norms = np.sqrt(np.sum(normals**2, axis=2))
for j in range(3):
normals[:, :, j] /= norms
# If the user provides an S_ref, we use that value.
# This value is multiplied by 2 if the problem is symmetric.
if self.surface['S_ref'] is not None:
S_ref = self.surface['S_ref']
else:
# Compute the wetted surface area
if self.surface['S_ref_type'] == 'wetted':
S_ref = 0.5 * np.sum(norms)
# Compute the projected surface area
elif self.surface['S_ref_type'] == 'projected':
proj_mesh = mesh.copy()
proj_mesh[: , :, 2] = 0.
proj_normals = np.cross(
proj_mesh[:-1, 1:, :] - proj_mesh[1:, :-1, :],
proj_mesh[:-1, :-1, :] - proj_mesh[1:, 1:, :],
axis=2)
proj_norms = np.sqrt(np.sum(proj_normals**2, axis=2))
for j in xrange(3):
proj_normals[:, :, j] /= proj_norms
S_ref = 0.5 * np.sum(proj_norms)
# Multiply the surface area by 2 if symmetric to get consistent area measures.
# However, only do this if we compute the area, not if the user inputs
# the reference area.
if self.surface['symmetry']:
S_ref *= 2
chords = np.linalg.norm(mesh[0, :, :] - mesh[-1, :, :], axis=1)
# Store each array in the unknowns dict
unknowns['b_pts'] = b_pts
unknowns['c_pts'] = c_pts
unknowns['widths'] = widths
unknowns['cos_sweep'] = cos_sweep
unknowns['lengths'] = lengths
unknowns['normals'] = normals
unknowns['S_ref'] = S_ref
unknowns['chords'] = chords
def linearize(self, params, unknowns, resids):
""" Jacobian for VLM geometry."""
jac = self.alloc_jacobian()
name = self.surface['name']
mesh = params['def_mesh']
nx = self.surface['num_x']
ny = self.surface['num_y']
if fortran_flag:
normalsb = np.zeros(unknowns['normals'].shape)
for i in range(nx-1):
for j in range(ny-1):
for ind in range(3):
normalsb[:, :, :] = 0.
normalsb[i, j, ind] = 1.
meshb, _, _ = OAS_API.oas_api.compute_normals_b(mesh, normalsb, 0.)
jac['normals', 'def_mesh'][i*(ny-1)*3 + j*3 + ind, :] = meshb.flatten()
normalsb[:, :, :] = 0.
if self.surface['S_ref_type'] == 'wetted':
seed_mesh = mesh
elif self.surface['S_ref_type'] == 'projected':
seed_mesh = mesh.copy()
seed_mesh[:, :, 2] = 0.
meshb, _, _ = OAS_API.oas_api.compute_normals_b(seed_mesh, normalsb, 1.)
jac['S_ref', 'def_mesh'] = np.atleast_2d(meshb.flatten())
if self.surface['symmetry']:
jac['S_ref', 'def_mesh'] *= 2
if self.surface['S_ref'] is not None:
jac['S_ref', 'def_mesh'][:] = 0.
else:
cs_jac = self.complex_step_jacobian(params, unknowns, resids,
fd_params=['def_mesh'],
fd_unknowns=['normals', 'S_ref'],
fd_states=[])
jac.update(cs_jac)
cs_jac = self.fd_jacobian(params, unknowns, resids,
fd_params=['def_mesh'],
fd_unknowns=['chords'],
fd_states=[])
jac.update(cs_jac)
for iz, v in zip((0, ny*3), (.75, .25)):
np.fill_diagonal(jac['b_pts', 'def_mesh'][:, iz:], v)
for iz, v in zip((0, 3, ny*3, (ny+1)*3),
(.125, .125, .375, .375)):
for ix in range(nx-1):
np.fill_diagonal(jac['c_pts', 'def_mesh']
[(ix*(ny-1))*3:((ix+1)*(ny-1))*3, iz+ix*ny*3:], v)
jac['widths', 'def_mesh'] = np.zeros_like(jac['widths', 'def_mesh'])
jac['cos_sweep', 'def_mesh'] = np.zeros_like(jac['cos_sweep', 'def_mesh'])
qc = 0.25 * mesh[-1] + 0.75 * mesh[0]
gap = [0, (nx-1)*ny*3]
factor = [0.75, 0.25]
for i in range(ny-1):
w = unknowns['widths'][i]
cos_sweep = unknowns['cos_sweep'][i]
dx = (qc[i+1, 0] - qc[i, 0])
dy = (qc[i+1, 1] - qc[i, 1])
dz = (qc[i+1, 2] - qc[i, 2])
for j in range(2):
jac['widths', 'def_mesh'][i, i*3+gap[j]] -= dx * factor[j] / w
jac['widths', 'def_mesh'][i, (i+1)*3+gap[j]] += dx * factor[j] / w
jac['widths', 'def_mesh'][i, i*3+1+gap[j]] -= dy * factor[j] / w
jac['widths', 'def_mesh'][i, (i+1)*3+1+gap[j]] += dy * factor[j] / w
jac['widths', 'def_mesh'][i, i*3+2+gap[j]] -= dz * factor[j] / w
jac['widths', 'def_mesh'][i, (i+1)*3+2+gap[j]] += dz * factor[j] / w
jac['cos_sweep', 'def_mesh'][i, i*3+1+gap[j]] -= dy / cos_sweep * factor[j]
jac['cos_sweep', 'def_mesh'][i, (i+1)*3+1+gap[j]] += dy / cos_sweep * factor[j]
jac['cos_sweep', 'def_mesh'][i, i*3+2+gap[j]] -= dz / cos_sweep * factor[j]
jac['cos_sweep', 'def_mesh'][i, (i+1)*3+2+gap[j]] += dz / cos_sweep * factor[j]
jac['lengths', 'def_mesh'] = np.zeros_like(jac['lengths', 'def_mesh'])
for i in range(ny):
dx = mesh[1:, i, 0] - mesh[:-1, i, 0]
dy = mesh[1:, i, 1] - mesh[:-1, i, 1]
dz = mesh[1:, i, 2] - mesh[:-1, i, 2]
for j in range(nx-1):
l = np.sqrt(dx[j]**2 + dy[j]**2 + dz[j]**2)
jac['lengths', 'def_mesh'][i, (j*ny+i)*3] -= dx[j] / l
jac['lengths', 'def_mesh'][i, ((j+1)*ny+i)*3] += dx[j] / l
jac['lengths', 'def_mesh'][i, (j*ny+i)*3 + 1] -= dy[j] / l
jac['lengths', 'def_mesh'][i, ((j+1)*ny+i)*3 + 1] += dy[j] / l
jac['lengths', 'def_mesh'][i, (j*ny+i)*3 + 2] -= dz[j] / l
jac['lengths', 'def_mesh'][i, ((j+1)*ny+i)*3 + 2] += dz[j] / l
return jac
class AssembleAIC(Component):
"""
Compute the circulations based on the AIC matrix and the panel velocities.
Note that the flow tangency condition is enforced at the 3/4 chord point.
There are multiple versions of the first four parameters with one
for each surface defined.
Each of these four parameters has the name of the surface prepended on the
actual parameter name.
Parameters
----------
def_mesh[nx, ny, 3] : numpy array
Array defining the nodal coordinates of the lifting surface.
b_pts[nx-1, ny, 3] : numpy array
Bound points for the horseshoe vortices, found along the 1/4 chord.
c_pts[nx-1, ny-1, 3] : numpy array
Collocation points on the 3/4 chord line where the flow tangency
condition is satisfed. Used to set up the linear system.
normals[nx-1, ny-1, 3] : numpy array
The normal vector for each panel, computed as the cross of the two
diagonals from the mesh points.
v : float
Freestream air velocity in m/s.
alpha : float
Angle of attack in degrees.
Returns
-------
AIC[(nx-1)*(ny-1), (nx-1)*(ny-1)] : numpy array
The aerodynamic influence coefficient matrix. Solving the linear system
of AIC * circulations = n * v gives us the circulations for each of the
horseshoe vortices.
rhs[(nx-1)*(ny-1)] : numpy array
The right-hand-side of the linear system that yields the circulations.
"""
def __init__(self, surfaces):
super(AssembleAIC, self).__init__()
self.surfaces = surfaces
tot_panels = 0
for surface in surfaces:
self.surface = surface
ny = surface['num_y']
nx = surface['num_x']
name = surface['name']
self.add_param(name+'def_mesh', val=np.zeros((nx, ny, 3),
dtype=data_type))
self.add_param(name+'b_pts', val=np.zeros((nx-1, ny, 3),
dtype=data_type))
self.add_param(name+'c_pts', val=np.zeros((nx-1, ny-1, 3),
dtype=data_type))
self.add_param(name+'normals', val=np.zeros((nx-1, ny-1, 3)))
tot_panels += (nx - 1) * (ny - 1)
self.tot_panels = tot_panels
self.add_param('v', val=1.)
self.add_param('alpha', val=0.)
self.add_output('AIC', val=np.zeros((tot_panels, tot_panels), dtype=data_type))
self.add_output('rhs', val=np.zeros((tot_panels), dtype=data_type))
self.AIC_mtx = np.zeros((tot_panels, tot_panels, 3),
dtype=data_type)
self.mtx = np.zeros((tot_panels, tot_panels),
dtype=data_type)
if not fortran_flag:
self.deriv_options['type'] = 'cs'
self.deriv_options['form'] = 'central'
def solve_nonlinear(self, params, unknowns, resids):
# Actually assemble the AIC matrix
_assemble_AIC_mtx(self.AIC_mtx, params, self.surfaces)
# Construct an flattened array with the normals of each surface in order
# so we can do the normals with velocities to set up the right-hand-side
# of the system.
flattened_normals = np.zeros((self.tot_panels, 3), dtype=data_type)
i = 0
for surface in self.surfaces:
name = surface['name']
num_panels = (surface['num_x'] - 1) * (surface['num_y'] - 1)
flattened_normals[i:i+num_panels, :] = params[name+'normals'].reshape(-1, 3, order='F')
i += num_panels
# Construct a matrix that is the AIC_mtx dotted by the normals at each
# collocation point. This is used to compute the circulations
self.mtx[:, :] = 0.
for ind in range(3):
self.mtx[:, :] += (self.AIC_mtx[:, :, ind].T *
flattened_normals[:, ind]).T
# Obtain the freestream velocity direction and magnitude by taking
# alpha into account
alpha = params['alpha'] * np.pi / 180.
cosa = np.cos(alpha)
sina = np.sin(alpha)
v_inf = params['v'] * np.array([cosa, 0., sina], dtype=data_type)
# Populate the right-hand side of the linear system with the
# expected velocities at each collocation point
unknowns['rhs'] = -flattened_normals.\
reshape(-1, flattened_normals.shape[-1], order='F').dot(v_inf)
unknowns['AIC'] = self.mtx
def apply_linear(self, params, unknowns, dparams, dunknowns, dresids, mode):
if mode == 'fwd':
AIC_mtxd = np.zeros(self.AIC_mtx.shape)
# Actually assemble the AIC matrix
_assemble_AIC_mtx_d(AIC_mtxd, params, dparams, dunknowns, dresids, self.surfaces)
# Construct an flattened array with the normals of each surface in order
# so we can do the normals with velocities to set up the right-hand-side
# of the system.
flattened_normals = np.zeros((self.tot_panels, 3))
flattened_normalsd = np.zeros((self.tot_panels, 3))
i = 0
for surface in self.surfaces:
name = surface['name']
num_panels = (surface['num_x'] - 1) * (surface['num_y'] - 1)
flattened_normals[i:i+num_panels, :] = params[name+'normals'].reshape(-1, 3, order='F')
flattened_normalsd[i:i+num_panels, :] = dparams[name+'normals'].reshape(-1, 3, order='F')
i += num_panels
# Construct a matrix that is the AIC_mtx dotted by the normals at each
# collocation point. This is used to compute the circulations
self.mtx[:, :] = 0.
for ind in range(3):
self.mtx[:, :] += (AIC_mtxd[:, :, ind].T *
flattened_normals[:, ind]).T
self.mtx[:, :] += (self.AIC_mtx[:, :, ind].T *
flattened_normalsd[:, ind]).T
# Obtain the freestream velocity direction and magnitude by taking
# alpha into account
alpha = params['alpha'] * np.pi / 180.
alphad = dparams['alpha'] * np.pi / 180.
cosa = np.cos(alpha)
sina = np.sin(alpha)
cosad = -sina * alphad
sinad = cosa * alphad
freestream_direction = np.array([cosa, 0., sina])
v_inf = params['v'] * freestream_direction
v_infd = dparams['v'] * freestream_direction
v_infd += params['v'] * np.array([cosad, 0., sinad])
# Populate the right-hand side of the linear system with the
# expected velocities at each collocation point
dresids['rhs'] = -flattened_normalsd.\
reshape(-1, 3, order='F').dot(v_inf)
dresids['rhs'] += -flattened_normals.\
reshape(-1, 3, order='F').dot(v_infd)
dresids['AIC'] = self.mtx
if mode == 'rev':
# Construct an flattened array with the normals of each surface in order
# so we can do the normals with velocities to set up the right-hand-side
# of the system.
flattened_normals = np.zeros((self.tot_panels, 3))
i = 0
for surface in self.surfaces:
name = surface['name']
num_panels = (surface['num_x'] - 1) * (surface['num_y'] - 1)
flattened_normals[i:i+num_panels, :] = params[name+'normals'].reshape(-1, 3, order='F')
i += num_panels
AIC_mtxb = np.zeros((self.tot_panels, self.tot_panels, 3))
flattened_normalsb = np.zeros(flattened_normals.shape)
for ind in range(3):
AIC_mtxb[:, :, ind] = (dresids['AIC'].T * flattened_normals[:, ind]).T
flattened_normalsb[:, ind] += np.sum(self.AIC_mtx[:, :, ind].real * dresids['AIC'], axis=1).T
# Actually assemble the AIC matrix
_assemble_AIC_mtx_b(AIC_mtxb, params, dparams, dunknowns, dresids, self.surfaces)
# Obtain the freestream velocity direction and magnitude by taking
# alpha into account
alpha = params['alpha'] * np.pi / 180.
cosa = np.cos(alpha)
sina = np.sin(alpha)
arr = np.array([cosa, 0., sina])
v_inf = params['v'] * arr
fn = flattened_normals
fnb = np.zeros(fn.shape)
rhsb = dresids['rhs']
v_infb = 0.
for ind in reversed(range(self.tot_panels)):
fnb[ind, :] -= v_inf * rhsb[ind]
v_infb -= fn[ind, :] * rhsb[ind]
dparams['v'] += sum(arr * v_infb)
arrb = params['v'] * v_infb
alphab = np.cos(alpha) * arrb[2]
alphab -= np.sin(alpha) * arrb[0]
alphab *= np.pi / 180.
dparams['alpha'] += alphab
i = 0
for surface in self.surfaces:
name = surface['name']
nx = surface['num_x']
ny = surface['num_y']
num_panels = (nx - 1) * (ny - 1)
dparams[name+'normals'] += flattened_normalsb[i:i+num_panels, :].reshape(nx-1, ny-1, 3, order='F')
dparams[name+'normals'] += fnb[i:i+num_panels, :].reshape(nx-1, ny-1, 3, order='F')
i += num_panels
class AeroCirculations(Component):
"""
Compute the circulation strengths of the horseshoe vortices by solving the
linear system AIC * circulations = n * v.
This component is copied from OpenMDAO's LinearSystem component with the
names of the parameters and outputs changed to match our problem formulation.
Parameters
----------
AIC[(nx-1)*(ny-1), (nx-1)*(ny-1)] : numpy array
The aerodynamic influence coefficient matrix. Solving the linear system
of AIC * circulations = n * v gives us the circulations for each of the
horseshoe vortices.
rhs[(nx-1)*(ny-1)] : numpy array
The right-hand-side of the linear system that yields the circulations.
Returns
-------
circulations[(nx-1)*(ny-1)] : numpy array
Augmented displacement array. Obtained by solving the system
AIC * circulations = n * v.
"""
def __init__(self, size):
super(AeroCirculations, self).__init__()
self.add_param('AIC', val=np.zeros((size, size), dtype=data_type))
self.add_param('rhs', val=np.zeros((size), dtype=data_type))
self.add_state('circulations', val=np.zeros((size), dtype=data_type))
self.size = size
# cache
self.lup = None
self.rhs_cache = None
def solve_nonlinear(self, params, unknowns, resids):
if np.any(np.isnan(params['AIC'])) or np.any(np.isinf(params['AIC'])):
raise AnalysisError
# lu factorization for use with solve_linear
self.lup = lu_factor(params['AIC'])
unknowns['circulations'] = lu_solve(self.lup, params['rhs'])
resids['circulations'] = params['AIC'].dot(unknowns['circulations']) - params['rhs']
def apply_nonlinear(self, params, unknowns, resids):
"""Evaluating residual for given state."""
resids['circulations'] = params['AIC'].dot(unknowns['circulations']) - params['rhs']
def apply_linear(self, params, unknowns, dparams, dunknowns, dresids, mode):
""" Apply the derivative of state variable with respect to
everything."""
if mode == 'fwd':
if 'circulations' in dunknowns:
dresids['circulations'] += params['AIC'].dot(dunknowns['circulations'])
if 'AIC' in dparams:
dresids['circulations'] += dparams['AIC'].dot(unknowns['circulations'])
if 'rhs' in dparams:
dresids['circulations'] -= dparams['rhs']
elif mode == 'rev':
if 'circulations' in dunknowns:
dunknowns['circulations'] += params['AIC'].T.dot(dresids['circulations'])
if 'AIC' in dparams: