forked from ImperialCollegeLondon/sharpy
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathgenerate_hale.py
708 lines (627 loc) · 27.6 KB
/
generate_hale.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
#! /usr/bin/env python3
import h5py as h5
import numpy as np
import os
import sharpy.utils.algebra as algebra
case_name = 'hale'
route = os.path.dirname(os.path.realpath(__file__)) + '/'
# EXECUTION
flow = ['BeamLoader',
'AerogridLoader',
'StaticTrim',
'DynamicCoupled',
'BeamLoads'
]
# if free_flight is False, the motion of the centre of the wing is prescribed.
free_flight = True
if not free_flight:
case_name += '_prescribed'
amplitude = 0 * np.pi / 180
period = 3
case_name += '_amp_' + str(amplitude).replace('.', '') + '_period_' + str(period)
# FLIGHT CONDITIONS
# the simulation is set such that the aircraft flies at a u_inf velocity while
# the air is calm.
u_inf = 10
rho = 1.225
# trim sigma = 1.5
alpha = 4.31 * np.pi / 180
beta = 0
roll = 0
gravity = 'on'
cs_deflection = -2.08 * np.pi / 180
rudder_static_deflection = 0.0
rudder_step = 0.0 * np.pi / 180
thrust = 6.16
sigma = 1.5
lambda_dihedral = 20 * np.pi / 180
# gust settings
gust_intensity = 0.20
gust_length = 1 * u_inf
gust_offset = 0.2 * u_inf
# numerics
n_step = 5
structural_relaxation_factor = 0.6
relaxation_factor = 0.35
tolerance = 1e-6
fsi_tolerance = 1e-4
num_cores = 2
# MODEL GEOMETRY
# beam
span_main = 16.0
lambda_main = 0.25
ea_main = 0.3
ea = 1e7
ga = 1e5
gj = 1e4
eiy = 2e4
eiz = 4e6
m_bar_main = 0.75
j_bar_main = 0.075
length_fuselage = 10
offset_fuselage = 0
sigma_fuselage = 10
m_bar_fuselage = 0.2
j_bar_fuselage = 0.08
span_tail = 2.5
ea_tail = 0.5
fin_height = 2.5
ea_fin = 0.5
sigma_tail = 100
m_bar_tail = 0.3
j_bar_tail = 0.08
# lumped masses
n_lumped_mass = 1
lumped_mass_nodes = np.zeros((n_lumped_mass,), dtype=int)
lumped_mass = np.zeros((n_lumped_mass,))
lumped_mass[0] = 50
lumped_mass_inertia = np.zeros((n_lumped_mass, 3, 3))
lumped_mass_position = np.zeros((n_lumped_mass, 3))
# aero
chord_main = 1.0
chord_tail = 0.5
chord_fin = 0.5
# DISCRETISATION
# spatial discretisation
# chordiwse panels
m = 4
# spanwise elements
n_elem_multiplier = 2
n_elem_main = int(4 * n_elem_multiplier)
n_elem_tail = int(2 * n_elem_multiplier)
n_elem_fin = int(2 * n_elem_multiplier)
n_elem_fuselage = int(2 * n_elem_multiplier)
n_surfaces = 5
# temporal discretisation
physical_time = 30
tstep_factor = 1.
dt = 1.0 / m / u_inf * tstep_factor
n_tstep = 20
# END OF INPUT-----------------------------------------------------------------
# beam processing
n_node_elem = 3
span_main1 = (1.0 - lambda_main) * span_main
span_main2 = lambda_main * span_main
n_elem_main1 = round(n_elem_main * (1 - lambda_main))
n_elem_main2 = n_elem_main - n_elem_main1
# total number of elements
n_elem = 0
n_elem += n_elem_main1 + n_elem_main1
n_elem += n_elem_main2 + n_elem_main2
n_elem += n_elem_fuselage
n_elem += n_elem_fin
n_elem += n_elem_tail + n_elem_tail
# number of nodes per part
n_node_main1 = n_elem_main1 * (n_node_elem - 1) + 1
n_node_main2 = n_elem_main2 * (n_node_elem - 1) + 1
n_node_main = n_node_main1 + n_node_main2 - 1
n_node_fuselage = n_elem_fuselage * (n_node_elem - 1) + 1
n_node_fin = n_elem_fin * (n_node_elem - 1) + 1
n_node_tail = n_elem_tail * (n_node_elem - 1) + 1
# total number of nodes
n_node = 0
n_node += n_node_main1 + n_node_main1 - 1
n_node += n_node_main2 - 1 + n_node_main2 - 1
n_node += n_node_fuselage - 1
n_node += n_node_fin - 1
n_node += n_node_tail - 1
n_node += n_node_tail - 1
# stiffness and mass matrices
n_stiffness = 3
base_stiffness_main = sigma * np.diag([ea, ga, ga, gj, eiy, eiz])
base_stiffness_fuselage = base_stiffness_main.copy() * sigma_fuselage
base_stiffness_fuselage[4, 4] = base_stiffness_fuselage[5, 5]
base_stiffness_tail = base_stiffness_main.copy() * sigma_tail
base_stiffness_tail[4, 4] = base_stiffness_tail[5, 5]
n_mass = 3
base_mass_main = np.diag([m_bar_main, m_bar_main, m_bar_main, j_bar_main, 0.5 * j_bar_main, 0.5 * j_bar_main])
base_mass_fuselage = np.diag([m_bar_fuselage,
m_bar_fuselage,
m_bar_fuselage,
j_bar_fuselage,
j_bar_fuselage * 0.5,
j_bar_fuselage * 0.5])
base_mass_tail = np.diag([m_bar_tail,
m_bar_tail,
m_bar_tail,
j_bar_tail,
j_bar_tail * 0.5,
j_bar_tail * 0.5])
# PLACEHOLDERS
# beam
x = np.zeros((n_node,))
y = np.zeros((n_node,))
z = np.zeros((n_node,))
beam_number = np.zeros((n_elem,), dtype=int)
frame_of_reference_delta = np.zeros((n_elem, n_node_elem, 3))
structural_twist = np.zeros((n_elem, 3))
conn = np.zeros((n_elem, n_node_elem), dtype=int)
stiffness = np.zeros((n_stiffness, 6, 6))
elem_stiffness = np.zeros((n_elem,), dtype=int)
mass = np.zeros((n_mass, 6, 6))
elem_mass = np.zeros((n_elem,), dtype=int)
boundary_conditions = np.zeros((n_node,), dtype=int)
app_forces = np.zeros((n_node, 6))
# aero
airfoil_distribution = np.zeros((n_elem, n_node_elem), dtype=int)
surface_distribution = np.zeros((n_elem,), dtype=int) - 1
surface_m = np.zeros((n_surfaces,), dtype=int)
m_distribution = 'uniform'
aero_node = np.zeros((n_node,), dtype=bool)
twist = np.zeros((n_elem, n_node_elem))
sweep = np.zeros((n_elem, n_node_elem))
chord = np.zeros((n_elem, n_node_elem,))
elastic_axis = np.zeros((n_elem, n_node_elem,))
# FUNCTIONS-------------------------------------------------------------
def clean_test_files():
fem_file_name = route + '/' + case_name + '.fem.h5'
if os.path.isfile(fem_file_name):
os.remove(fem_file_name)
dyn_file_name = route + '/' + case_name + '.dyn.h5'
if os.path.isfile(dyn_file_name):
os.remove(dyn_file_name)
aero_file_name = route + '/' + case_name + '.aero.h5'
if os.path.isfile(aero_file_name):
os.remove(aero_file_name)
solver_file_name = route + '/' + case_name + '.sharpy'
if os.path.isfile(solver_file_name):
os.remove(solver_file_name)
flightcon_file_name = route + '/' + case_name + '.flightcon.txt'
if os.path.isfile(flightcon_file_name):
os.remove(flightcon_file_name)
def generate_fem():
stiffness[0, ...] = base_stiffness_main
stiffness[1, ...] = base_stiffness_fuselage
stiffness[2, ...] = base_stiffness_tail
mass[0, ...] = base_mass_main
mass[1, ...] = base_mass_fuselage
mass[2, ...] = base_mass_tail
we = 0
wn = 0
# inner right wing
beam_number[we:we + n_elem_main1] = 0
y[wn:wn + n_node_main1] = np.linspace(0.0, span_main1, n_node_main1)
for ielem in range(n_elem_main1):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [-1.0, 0.0, 0.0]
elem_stiffness[we:we + n_elem_main1] = 0
elem_mass[we:we + n_elem_main1] = 0
boundary_conditions[0] = 1
# remember this is in B FoR
app_forces[0] = [0, thrust, 0, 0, 0, 0]
we += n_elem_main1
wn += n_node_main1
# outer right wing
beam_number[we:we + n_elem_main1] = 0
y[wn:wn + n_node_main2 - 1] = y[wn - 1] + np.linspace(0.0, np.cos(lambda_dihedral) * span_main2, n_node_main2)[1:]
z[wn:wn + n_node_main2 - 1] = z[wn - 1] + np.linspace(0.0, np.sin(lambda_dihedral) * span_main2, n_node_main2)[1:]
for ielem in range(n_elem_main2):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [-1.0, 0.0, 0.0]
elem_stiffness[we:we + n_elem_main2] = 0
elem_mass[we:we + n_elem_main2] = 0
boundary_conditions[wn + n_node_main2 - 2] = -1
we += n_elem_main2
wn += n_node_main2 - 1
# inner left wing
beam_number[we:we + n_elem_main1 - 1] = 1
y[wn:wn + n_node_main1 - 1] = np.linspace(0.0, -span_main1, n_node_main1)[1:]
for ielem in range(n_elem_main1):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [1.0, 0.0, 0.0]
conn[we, 0] = 0
elem_stiffness[we:we + n_elem_main1] = 0
elem_mass[we:we + n_elem_main1] = 0
we += n_elem_main1
wn += n_node_main1 - 1
# outer left wing
beam_number[we:we + n_elem_main2] = 1
y[wn:wn + n_node_main2 - 1] = y[wn - 1] + np.linspace(0.0, -np.cos(lambda_dihedral) * span_main2, n_node_main2)[1:]
z[wn:wn + n_node_main2 - 1] = z[wn - 1] + np.linspace(0.0, np.sin(lambda_dihedral) * span_main2, n_node_main2)[1:]
for ielem in range(n_elem_main2):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [1.0, 0.0, 0.0]
elem_stiffness[we:we + n_elem_main2] = 0
elem_mass[we:we + n_elem_main2] = 0
boundary_conditions[wn + n_node_main2 - 2] = -1
we += n_elem_main2
wn += n_node_main2 - 1
# fuselage
beam_number[we:we + n_elem_fuselage] = 2
x[wn:wn + n_node_fuselage - 1] = np.linspace(0.0, length_fuselage, n_node_fuselage)[1:]
z[wn:wn + n_node_fuselage - 1] = np.linspace(0.0, offset_fuselage, n_node_fuselage)[1:]
for ielem in range(n_elem_fuselage):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [0.0, 1.0, 0.0]
conn[we, 0] = 0
elem_stiffness[we:we + n_elem_fuselage] = 1
elem_mass[we:we + n_elem_fuselage] = 1
we += n_elem_fuselage
wn += n_node_fuselage - 1
global end_of_fuselage_node
end_of_fuselage_node = wn - 1
# fin
beam_number[we:we + n_elem_fin] = 3
x[wn:wn + n_node_fin - 1] = x[end_of_fuselage_node]
z[wn:wn + n_node_fin - 1] = z[end_of_fuselage_node] + np.linspace(0.0, fin_height, n_node_fin)[1:]
for ielem in range(n_elem_fin):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [-1.0, 0.0, 0.0]
conn[we, 0] = end_of_fuselage_node
elem_stiffness[we:we + n_elem_fin] = 2
elem_mass[we:we + n_elem_fin] = 2
we += n_elem_fin
wn += n_node_fin - 1
end_of_fin_node = wn - 1
# right tail
beam_number[we:we + n_elem_tail] = 4
x[wn:wn + n_node_tail - 1] = x[end_of_fin_node]
y[wn:wn + n_node_tail - 1] = np.linspace(0.0, span_tail, n_node_tail)[1:]
z[wn:wn + n_node_tail - 1] = z[end_of_fin_node]
for ielem in range(n_elem_tail):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [-1.0, 0.0, 0.0]
conn[we, 0] = end_of_fin_node
elem_stiffness[we:we + n_elem_tail] = 2
elem_mass[we:we + n_elem_tail] = 2
boundary_conditions[wn + n_node_tail - 2] = -1
we += n_elem_tail
wn += n_node_tail - 1
# left tail
beam_number[we:we + n_elem_tail] = 5
x[wn:wn + n_node_tail - 1] = x[end_of_fin_node]
y[wn:wn + n_node_tail - 1] = np.linspace(0.0, -span_tail, n_node_tail)[1:]
z[wn:wn + n_node_tail - 1] = z[end_of_fin_node]
for ielem in range(n_elem_tail):
conn[we + ielem, :] = ((np.ones((3,)) * (we + ielem) * (n_node_elem - 1)) +
[0, 2, 1])
for inode in range(n_node_elem):
frame_of_reference_delta[we + ielem, inode, :] = [1.0, 0.0, 0.0]
conn[we, 0] = end_of_fin_node
elem_stiffness[we:we + n_elem_tail] = 2
elem_mass[we:we + n_elem_tail] = 2
boundary_conditions[wn + n_node_tail - 2] = -1
we += n_elem_tail
wn += n_node_tail - 1
with h5.File(route + '/' + case_name + '.fem.h5', 'a') as h5file:
coordinates = h5file.create_dataset('coordinates', data=np.column_stack((x, y, z)))
conectivities = h5file.create_dataset('connectivities', data=conn)
num_nodes_elem_handle = h5file.create_dataset(
'num_node_elem', data=n_node_elem)
num_nodes_handle = h5file.create_dataset(
'num_node', data=n_node)
num_elem_handle = h5file.create_dataset(
'num_elem', data=n_elem)
stiffness_db_handle = h5file.create_dataset(
'stiffness_db', data=stiffness)
stiffness_handle = h5file.create_dataset(
'elem_stiffness', data=elem_stiffness)
mass_db_handle = h5file.create_dataset(
'mass_db', data=mass)
mass_handle = h5file.create_dataset(
'elem_mass', data=elem_mass)
frame_of_reference_delta_handle = h5file.create_dataset(
'frame_of_reference_delta', data=frame_of_reference_delta)
structural_twist_handle = h5file.create_dataset(
'structural_twist', data=structural_twist)
bocos_handle = h5file.create_dataset(
'boundary_conditions', data=boundary_conditions)
beam_handle = h5file.create_dataset(
'beam_number', data=beam_number)
app_forces_handle = h5file.create_dataset(
'app_forces', data=app_forces)
lumped_mass_nodes_handle = h5file.create_dataset(
'lumped_mass_nodes', data=lumped_mass_nodes)
lumped_mass_handle = h5file.create_dataset(
'lumped_mass', data=lumped_mass)
lumped_mass_inertia_handle = h5file.create_dataset(
'lumped_mass_inertia', data=lumped_mass_inertia)
lumped_mass_position_handle = h5file.create_dataset(
'lumped_mass_position', data=lumped_mass_position)
def generate_aero_file():
global x, y, z
# control surfaces
n_control_surfaces = 2
control_surface = np.zeros((n_elem, n_node_elem), dtype=int) - 1
control_surface_type = np.zeros((n_control_surfaces,), dtype=int)
control_surface_deflection = np.zeros((n_control_surfaces,))
control_surface_chord = np.zeros((n_control_surfaces,), dtype=int)
control_surface_hinge_coord = np.zeros((n_control_surfaces,), dtype=float)
# control surface type 0 = static
# control surface type 1 = dynamic
control_surface_type[0] = 0
control_surface_deflection[0] = cs_deflection
control_surface_chord[0] = m
control_surface_hinge_coord[0] = -0.25 # nondimensional wrt elastic axis (+ towards the trailing edge)
control_surface_type[1] = 0
control_surface_deflection[1] = rudder_static_deflection
control_surface_chord[1] = 1
control_surface_hinge_coord[1] = -0. # nondimensional wrt elastic axis (+ towards the trailing edge)
we = 0
wn = 0
# right wing (surface 0, beam 0)
i_surf = 0
airfoil_distribution[we:we + n_elem_main, :] = 0
surface_distribution[we:we + n_elem_main] = i_surf
surface_m[i_surf] = m
aero_node[wn:wn + n_node_main] = True
temp_chord = np.linspace(chord_main, chord_main, n_node_main)
temp_sweep = np.linspace(0.0, 0 * np.pi / 180, n_node_main)
node_counter = 0
for i_elem in range(we, we + n_elem_main):
for i_local_node in range(n_node_elem):
if not i_local_node == 0:
node_counter += 1
chord[i_elem, i_local_node] = temp_chord[node_counter]
elastic_axis[i_elem, i_local_node] = ea_main
sweep[i_elem, i_local_node] = temp_sweep[node_counter]
we += n_elem_main
wn += n_node_main
# left wing (surface 1, beam 1)
i_surf = 1
airfoil_distribution[we:we + n_elem_main, :] = 0
# airfoil_distribution[wn:wn + n_node_main - 1] = 0
surface_distribution[we:we + n_elem_main] = i_surf
surface_m[i_surf] = m
aero_node[wn:wn + n_node_main - 1] = True
# chord[wn:wn + num_node_main - 1] = np.linspace(main_chord, main_tip_chord, num_node_main)[1:]
# chord[wn:wn + num_node_main - 1] = main_chord
# elastic_axis[wn:wn + num_node_main - 1] = main_ea
temp_chord = np.linspace(chord_main, chord_main, n_node_main)
node_counter = 0
for i_elem in range(we, we + n_elem_main):
for i_local_node in range(n_node_elem):
if not i_local_node == 0:
node_counter += 1
chord[i_elem, i_local_node] = temp_chord[node_counter]
elastic_axis[i_elem, i_local_node] = ea_main
sweep[i_elem, i_local_node] = -temp_sweep[node_counter]
we += n_elem_main
wn += n_node_main - 1
we += n_elem_fuselage
wn += n_node_fuselage - 1 - 1
#
# # fin (surface 2, beam 3)
i_surf = 2
airfoil_distribution[we:we + n_elem_fin, :] = 1
# airfoil_distribution[wn:wn + n_node_fin] = 0
surface_distribution[we:we + n_elem_fin] = i_surf
surface_m[i_surf] = m
aero_node[wn:wn + n_node_fin] = True
# chord[wn:wn + num_node_fin] = fin_chord
for i_elem in range(we, we + n_elem_fin):
for i_local_node in range(n_node_elem):
chord[i_elem, i_local_node] = chord_fin
elastic_axis[i_elem, i_local_node] = ea_fin
control_surface[i_elem, i_local_node] = 1
# twist[end_of_fuselage_node] = 0
# twist[wn:] = 0
# elastic_axis[wn:wn + num_node_main] = fin_ea
we += n_elem_fin
wn += n_node_fin - 1
#
# # # right tail (surface 3, beam 4)
i_surf = 3
airfoil_distribution[we:we + n_elem_tail, :] = 2
# airfoil_distribution[wn:wn + n_node_tail] = 0
surface_distribution[we:we + n_elem_tail] = i_surf
surface_m[i_surf] = m
# XXX not very elegant
aero_node[wn:] = True
# chord[wn:wn + num_node_tail] = tail_chord
# elastic_axis[wn:wn + num_node_main] = tail_ea
for i_elem in range(we, we + n_elem_tail):
for i_local_node in range(n_node_elem):
twist[i_elem, i_local_node] = -0
for i_elem in range(we, we + n_elem_tail):
for i_local_node in range(n_node_elem):
chord[i_elem, i_local_node] = chord_tail
elastic_axis[i_elem, i_local_node] = ea_tail
control_surface[i_elem, i_local_node] = 0
we += n_elem_tail
wn += n_node_tail
#
# # left tail (surface 4, beam 5)
i_surf = 4
airfoil_distribution[we:we + n_elem_tail, :] = 2
# airfoil_distribution[wn:wn + n_node_tail - 1] = 0
surface_distribution[we:we + n_elem_tail] = i_surf
surface_m[i_surf] = m
aero_node[wn:wn + n_node_tail - 1] = True
# chord[wn:wn + num_node_tail] = tail_chord
# elastic_axis[wn:wn + num_node_main] = tail_ea
# twist[we:we + num_elem_tail] = -tail_twist
for i_elem in range(we, we + n_elem_tail):
for i_local_node in range(n_node_elem):
twist[i_elem, i_local_node] = -0
for i_elem in range(we, we + n_elem_tail):
for i_local_node in range(n_node_elem):
chord[i_elem, i_local_node] = chord_tail
elastic_axis[i_elem, i_local_node] = ea_tail
control_surface[i_elem, i_local_node] = 0
we += n_elem_tail
wn += n_node_tail
with h5.File(route + '/' + case_name + '.aero.h5', 'a') as h5file:
airfoils_group = h5file.create_group('airfoils')
# add one airfoil
naca_airfoil_main = airfoils_group.create_dataset('0', data=np.column_stack(
generate_naca_camber(P=0, M=0)))
naca_airfoil_tail = airfoils_group.create_dataset('1', data=np.column_stack(
generate_naca_camber(P=0, M=0)))
naca_airfoil_fin = airfoils_group.create_dataset('2', data=np.column_stack(
generate_naca_camber(P=0, M=0)))
# chord
chord_input = h5file.create_dataset('chord', data=chord)
dim_attr = chord_input.attrs['units'] = 'm'
# twist
twist_input = h5file.create_dataset('twist', data=twist)
dim_attr = twist_input.attrs['units'] = 'rad'
# sweep
sweep_input = h5file.create_dataset('sweep', data=sweep)
dim_attr = sweep_input.attrs['units'] = 'rad'
# airfoil distribution
airfoil_distribution_input = h5file.create_dataset('airfoil_distribution', data=airfoil_distribution)
surface_distribution_input = h5file.create_dataset('surface_distribution', data=surface_distribution)
surface_m_input = h5file.create_dataset('surface_m', data=surface_m)
m_distribution_input = h5file.create_dataset('m_distribution', data=m_distribution.encode('ascii', 'ignore'))
aero_node_input = h5file.create_dataset('aero_node', data=aero_node)
elastic_axis_input = h5file.create_dataset('elastic_axis', data=elastic_axis)
control_surface_input = h5file.create_dataset('control_surface', data=control_surface)
control_surface_deflection_input = h5file.create_dataset('control_surface_deflection',
data=control_surface_deflection)
control_surface_chord_input = h5file.create_dataset('control_surface_chord', data=control_surface_chord)
control_surface_hinge_coord_input = h5file.create_dataset('control_surface_hinge_coord',
data=control_surface_hinge_coord)
control_surface_types_input = h5file.create_dataset('control_surface_type', data=control_surface_type)
def generate_naca_camber(M=0, P=0):
mm = M * 1e-2
p = P * 1e-1
def naca(x, mm, p):
if x < 1e-6:
return 0.0
elif x < p:
return mm / (p * p) * (2 * p * x - x * x)
elif x > p and x < 1 + 1e-6:
return mm / ((1 - p) * (1 - p)) * (1 - 2 * p + 2 * p * x - x * x)
x_vec = np.linspace(0, 1, 1000)
y_vec = np.array([naca(x, mm, p) for x in x_vec])
return x_vec, y_vec
def generate_solver_file():
file_name = route + '/' + case_name + '.sharpy'
settings = dict()
settings['SHARPy'] = {'case': case_name,
'route': route,
'flow': flow,
'write_screen': 'off',
'write_log': 'off',
'log_folder': route + '/output/',
'log_file': case_name + '.log'}
settings['NonLinearStatic'] = {'print_info': 'off',
'max_iterations': 150,
'num_load_steps': 1,
'delta_curved': 1e-1,
'min_delta': tolerance,
'gravity_on': gravity,
'gravity': 9.81}
settings['StaticUvlm'] = {'print_info': 'on',
'horseshoe': 'off',
'num_cores': num_cores,
'n_rollup': 0,
'rollup_dt': dt,
'rollup_aic_refresh': 1,
'rollup_tolerance': 1e-4,
'velocity_field_generator': 'SteadyVelocityField',
'velocity_field_input': {'u_inf': u_inf,
'u_inf_direction': [1., 0, 0]},
'rho': rho}
settings['StaticCoupled'] = {'print_info': 'off',
'structural_solver': 'NonLinearStatic',
'structural_solver_settings': settings['NonLinearStatic'],
'aero_solver': 'StaticUvlm',
'aero_solver_settings': settings['StaticUvlm'],
'max_iter': 100,
'n_load_steps': n_step,
'tolerance': fsi_tolerance,
'relaxation_factor': structural_relaxation_factor}
settings['StaticTrim'] = {'solver': 'StaticCoupled',
'solver_settings': settings['StaticCoupled'],
'initial_alpha': alpha,
'initial_deflection': cs_deflection,
'initial_thrust': thrust}
settings['NonLinearDynamicCoupledStep'] = {'print_info': 'off',
'max_iterations': 950,
'delta_curved': 1e-1,
'min_delta': tolerance,
'newmark_damp': 5e-3,
'gravity_on': gravity,
'gravity': 9.81,
'num_steps': n_tstep,
'dt': dt,
'initial_velocity': u_inf}
relative_motion = 'off'
settings['StepUvlm'] = {'print_info': 'off',
'num_cores': num_cores,
'convection_scheme': 2,
'gamma_dot_filtering': 7,
'velocity_field_generator': 'GustVelocityField',
'velocity_field_input': {'u_inf': int(not free_flight) * u_inf,
'u_inf_direction': [1., 0, 0],
'gust_shape': '1-cos',
'gust_parameters': {'gust_length': gust_length,
'gust_intensity': gust_intensity * u_inf},
'offset': gust_offset,
'relative_motion': relative_motion},
'rho': rho,
'n_time_steps': n_tstep,
'dt': dt}
settings['BeamLoads'] = {'csv_output': True}
solver = 'NonLinearDynamicCoupledStep'
settings['DynamicCoupled'] = {'structural_solver': solver,
'structural_solver_settings': settings[solver],
'aero_solver': 'StepUvlm',
'aero_solver_settings': settings['StepUvlm'],
'fsi_substeps': 200,
'fsi_tolerance': fsi_tolerance,
'relaxation_factor': relaxation_factor,
'minimum_steps': 1,
'relaxation_steps': 150,
'final_relaxation_factor': 0.5,
'n_time_steps': n_tstep,
'dt': dt,
'include_unsteady_force_contribution': 'on',}
settings['BeamLoader'] = {'unsteady': 'on',
'orientation': algebra.euler2quat(np.array([roll,
alpha,
beta]))}
settings['AerogridLoader'] = {'unsteady': 'on',
'aligned_grid': 'on',
'mstar': int(20 / tstep_factor),
'freestream_dir': ['1', '0', '0'],
'wake_shape_generator': 'StraightWake',
'wake_shape_generator_input': {'u_inf': u_inf,
'u_inf_direction': ['1', '0', '0'],
'dt': dt}}
import configobj
config = configobj.ConfigObj()
config.filename = file_name
for k, v in settings.items():
config[k] = v
config.write()
clean_test_files()
generate_fem()
generate_aero_file()
generate_solver_file()