-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathhelmholtz.cc
1898 lines (1598 loc) · 77.7 KB
/
helmholtz.cc
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
/* ---------------------------------------------------------------------
*
* Copyright (C) 2021 by Wolfgang Bangerth and SAATI Co.
*
* ---------------------------------------------------------------------
*/
/* notes from wolfgang:
You can always get the latest version from the github repository:
https://urlsand.esvalabs.com/?u=https%3A%2F%2Fgithub.com%2Fbangerth%2Fhelmholtz&e=a9a925b0&h=f79f86f5&f=y&p=y
Click on helmholtz.cc
https://urlsand.esvalabs.com/?u=https%3A%2F%2Fgithub.com%2Fbangerth%2Fhelmholtz%2Fblob%2Fmaster%2Fhelmholtz.cc&e=a9a925b0&h=3338f99f&f=y&p=y
and then either click on "Raw" to see the actual file (which is also here:
https://urlsand.esvalabs.com/?u=https%3A%2F%2Fraw.githubusercontent.com%2Fbangerth%2Fhelmholtz%2Fmaster%2Fhelmholtz.cc&e=a9a925b0&h=a2d05443&f=y&p=y
) or two symbols to the right to copy it into your OS buffer for pasting into
an editor.
*/
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/function.h>
#include <deal.II/base/function_lib.h>
#include <deal.II/base/thread_management.h>
#include <deal.II/base/parameter_handler.h>
#include <deal.II/base/timer.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/vector.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
#include <deal.II/lac/solver_cg.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/sparse_direct.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/reference_cell.h>
#include <deal.II/grid/grid_in.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>
#include <deal.II/fe/fe_simplex_p.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/numerics/vector_tools.h>
#include <deal.II/numerics/matrix_tools.h>
#include <deal.II/numerics/data_out.h>
#include <deal.II/numerics/data_out_faces.h>
#include <deal.II/numerics/vector_tools_point_value.h>
#include <deal.II/numerics/vector_tools_point_gradient.h>
#include <deal.II/fe/fe_interface_values.h>
#include <deal.II/meshworker/mesh_loop.h>
#include <fstream>
#include <future>
#include <thread>
#include <iostream>
#include <sstream>
#include <cmath>
#include <cstdio>
#include <complex>
#include <memory>
#include <regex>
#ifdef _MSC_BUILD
#define dir_sep_d "\\"
#else
#define dir_sep_d "/"
#endif
// i'm having some errors writing to files. i'm testing to see if it has anything to do with the
// directory separators. apparenlty unix uses / and windows uses \.
// I don't think this is the problem, but i'm leaving this in until i fix it just in case it's related to this.
// these defines are also put into ares.exe, so they can't change unless Ares chages too.
#define dii_termination_signal_filename_base_d "dii_termination_signal.txt"
#define dii_output_log_filename_base_d "dii_output_log.txt"
#define dii_frequency_response_filename_base_d "dii_frequency_response.txt"
#define dii_simulation_end_filename_base_d "dii_finished_signal.txt"
// "simulation end" is also the "success" signal filname
#define dii_biharmonic_prm_filename_base_d "dii_fea_2d_mesh_solver.prm"
#define dii_visualization_dir_name_base_d "visualization"
// added for helmholtz 3D volume solver
#define dii_solver_failure_signal_base_d "dll_solver_failure_signal.txt"
#define dii_frequency_response_filename_csv_base_d "dii_frequency_response.csv"
#define dii_helmholtz_prm_filename_base_d "dii_helmholtz.prm"
#define dii_port_areas_base_d "dii_port_areas.txt"
std::string dii_termination_signal_filename_st = dii_termination_signal_filename_base_d;
std::string dii_output_log_filename_st = dii_output_log_filename_base_d;
std::string dii_frequency_response_filename_st = dii_frequency_response_filename_base_d;
std::string dii_simulation_end_filename_st = dii_simulation_end_filename_base_d;
std::string dii_visualization_dir_name_st = dii_visualization_dir_name_base_d;
std::string dii_solver_failure_signal_st = dii_solver_failure_signal_base_d;
std::string dii_frequency_response_filename_csv_st = dii_frequency_response_filename_csv_base_d;
std::string dii_helmholtz_prm_filename_st = dii_helmholtz_prm_filename_base_d;
std::string dii_port_areas_base_st = dii_port_areas_base_d;
std::string instance_folder;
std::string output_file_prefix;
int index_solver_base = 0; // added to the solver index when generating individual frequency output files
std::ofstream logger;
namespace TransmissionProblem
{
using namespace dealii;
using ScalarType = std::complex<double>;
// The following namespace defines material parameters. We use SI units.
namespace MaterialParameters
{
std::vector<std::pair<double,unsigned int>> frequencies;
std::unique_ptr<Functions::InterpolatedTensorProductGridData<1>> density_real;
std::unique_ptr<Functions::InterpolatedTensorProductGridData<1>> density_imag;
std::unique_ptr<Functions::InterpolatedTensorProductGridData<1>> bulk_modulus_real;
std::unique_ptr<Functions::InterpolatedTensorProductGridData<1>> bulk_modulus_imag;
}
std::string mesh_file_name;
double geometry_conversion_factor;
unsigned int fe_degree = 2;
int n_mesh_refinement_steps = 5;
unsigned int n_threads = 0;
std::vector<Point<3>> evaluation_points;
bool mesh_summary_only = false;
void
declare_parameters (ParameterHandler &prm)
{
prm.declare_entry ("Mesh file name", "./mesh.msh",
Patterns::FileName(),
"The name of the file from which to read the mesh. "
"The extents of the geometry so described are scaled "
"by the 'Geometry conversion factor to SI units' parameter.");
prm.declare_entry ("Geometry conversion factor to meters", "1",
Patterns::Double(),
"A conversion factor from whatever length units are used "
"for the geometry description and for the description of "
"arbitrary evaluation points, to meters. For example, if the mesh and "
"the evaluation points are given in multiples of inches, "
"then this factor should be set to 0.0254 because one inch "
"equals 0.0254 meters = 25.4 mm.");
prm.declare_entry ("Material properties file name",
"./material_properties.txt",
Patterns::FileName(Patterns::FileName::input),
"The name of the file from which to read the mechanical material response.");
prm.declare_entry ("Evaluation points", "",
Patterns::List (Patterns::List(Patterns::Double(),3,3,","),
0, Patterns::List::max_int_value, ";"),
"A list of points at which the program evaluates both the pressure "
"and the (volumetric) velocity. Each point is specified by x,y,z "
"coordinates using the same units as used for the mesh (and then "
"scaled by the value of the 'Geometry conversion factor to meters' "
"parameter). Points are separated by semicolons.");
prm.declare_entry ("Frequencies", "linear_spacing(100,10000,100)",
Patterns::Anything(),
"A description of the frequencies to compute. See "
"the readme.md file for a description of the format "
"for this entry.");
prm.declare_entry ("Number of mesh refinement steps", "5",
Patterns::Integer(-100,10),
"The number of global mesh refinement steps applied "
"to the coarse mesh if positive or zero. If negative, "
"then it denotes the number of mesh points per wave length "
"as described in readme.md.");
prm.declare_entry ("Finite element polynomial degree", "2",
Patterns::Integer(1,5),
"The polynomial degree to be used for the finite element.");
prm.declare_entry ("Number of threads", "0",
Patterns::Integer(0),
"The number of threads this program may use at the same time. "
"Threads are used to compute the frequency response for "
"different frequencies at the same time since these are "
"independent computations. A value of zero means that the "
"program may use as many threads as it pleases, whereas a "
"positive number limits how many threads (and consequently "
"CPU cores) the program will use.");
prm.declare_entry ("Mesh summary only", "false",
Patterns::Bool(),
"If set to `true', the program stops after reading the mesh "
"file and only output the summary information about the mesh "
"that would ordinarily be computed as part of the simulation "
"process. Specifically, this includes information about "
"port boundary indicators, port areas, and similar.");
}
void
read_parameters (ParameterHandler &prm)
{
// First read parameter values from the input file 'helmholtz.prm'
prm.parse_input (instance_folder + dir_sep_d + dii_helmholtz_prm_filename_st);
// Start with geometry things: The mesh, the scaling factor, evaluation points
mesh_file_name = prm.get ("Mesh file name");
geometry_conversion_factor = prm.get_double ("Geometry conversion factor to meters");
{
const auto points = Utilities::split_string_list (prm.get("Evaluation points"), ';');
for (const auto &point : points)
{
const auto coordinates = Utilities::split_string_list (point, ',');
AssertDimension (coordinates.size(), 3);
Point<3> p;
for (unsigned int d=0; d<3; ++d)
p[d] = Utilities::string_to_double(coordinates[d]);
p *= geometry_conversion_factor;
evaluation_points.push_back (p);
}
}
// Next up: Material parameters.
const std::string material_properties_file_name
= prm.get ("Material properties file name");
logger << "debug Material properties file name = " << material_properties_file_name << "\r\n";
{
std::ifstream material_input (instance_folder + dir_sep_d + material_properties_file_name);
AssertThrow (material_input,
ExcMessage("Can't read material properties from file <"
+ material_properties_file_name
+ ">."));
// discard the first line as a comment
{
std::string dummy;
std::getline(material_input, dummy);
}
std::vector<double> frequencies;
std::vector<double> rho_real;
std::vector<double> rho_imag;
std::vector<double> K_real;
std::vector<double> K_imag;
// Read until we find the end of the file
while (true)
{
double f, rr, ri, Kr, Ki;
material_input >> f >> rr >> ri >> Kr >> Ki;
if (material_input.good())
{
frequencies.push_back (f);
rho_real.push_back (rr);
rho_imag.push_back (ri);
K_real.push_back (Kr);
K_imag.push_back (Ki);
}
else
break;
}
logger << "INFO Material parameters file contains data for "
<< frequencies.size()
<< " frequencies ranging from "
<< frequencies.front()
<< " to "
<< frequencies.back()
<< "Hz."
<< std::endl;
// Now convert these arrays into Table objects as desired by the
// InterpolatedTensorProductGridData class
Table<1,double> rr (frequencies.size(), rho_real.begin());
Table<1,double> ri (frequencies.size(), rho_imag.begin());
Table<1,double> Kr (frequencies.size(), K_real.begin());
Table<1,double> Ki (frequencies.size(), K_imag.begin());
// And finally use these to initialize the
// InterpolatedTensorProductGridData objects themselves
const std::array<std::vector<double>,1> f({{frequencies}});
MaterialParameters::density_real
= std::make_unique<Functions::InterpolatedTensorProductGridData<1>>(f, rr);
MaterialParameters::density_imag
= std::make_unique<Functions::InterpolatedTensorProductGridData<1>>(f, ri);
MaterialParameters::bulk_modulus_real
= std::make_unique<Functions::InterpolatedTensorProductGridData<1>>(f, Kr);
MaterialParameters::bulk_modulus_imag
= std::make_unique<Functions::InterpolatedTensorProductGridData<1>>(f, Ki);
}
// Read and parse the entry that determines which frequencies to compute.
// Recall that the format is one of the following:
// - linear_spacing(min,max,n_steps)
// - exp_spacing(min,max,n_steps)
// - list(...)
const std::string frequency_descriptor = prm.get ("Frequencies");
if (frequency_descriptor.find ("linear_spacing") == 0)
{
// Get the rest of the string, and eat any space at the start and end
const std::string parenthesized_expr
= Utilities::trim (frequency_descriptor.substr
(std::string("linear_spacing").size(),
std::string::npos));
AssertThrow (parenthesized_expr.size() >= 2
&&
parenthesized_expr.front() == '('
&&
parenthesized_expr.back() == ')',
ExcMessage ("Wrong format for 'linear_spacing'."));
// Then get the interior part, again trim spaces, and split at
// commas
const std::vector<std::string> min_max_steps
= Utilities::split_string_list
(Utilities::trim (parenthesized_expr.substr
(1,
parenthesized_expr.size() - 2)),
',');
AssertThrow (min_max_steps.size() == 3,
ExcMessage ("Wrong format for 'linear_spacing'."));
const double min_omega = Utilities::string_to_double(min_max_steps[0])
* 2 * numbers::PI;
const double max_omega = Utilities::string_to_double(min_max_steps[1])
* 2 * numbers::PI;
const unsigned int n_frequencies = Utilities::string_to_int(min_max_steps[2]);
const double delta_omega = (max_omega - min_omega)
/ (n_frequencies-1)
* (1.-1e-12);
unsigned int index = 0;
for (double omega = min_omega;
omega <= max_omega;
omega += delta_omega, ++index)
MaterialParameters::frequencies.emplace_back (omega,index);
}
else if (frequency_descriptor.find ("exp_spacing") == 0)
{
// Get the rest of the string, and eat any space at the start and end
const std::string parenthesized_expr
= Utilities::trim (frequency_descriptor.substr
(std::string("exp_spacing").size(),
std::string::npos));
AssertThrow (parenthesized_expr.size() >= 2
&&
parenthesized_expr.front() == '('
&&
parenthesized_expr.back() == ')',
ExcMessage ("Wrong format for 'exp_spacing'."));
// Then get the interior part, again trim spaces, and split at
// commas
const std::vector<std::string> min_max_steps
= Utilities::split_string_list
(Utilities::trim (parenthesized_expr.substr
(1,
parenthesized_expr.size() - 2)),
',');
AssertThrow (min_max_steps.size() == 3,
ExcMessage ("Wrong format for 'exp_spacing'."));
const double log_min_omega = std::log(Utilities::string_to_double(min_max_steps[0])
* 2 * numbers::PI);
const double log_max_omega = std::log(Utilities::string_to_double(min_max_steps[1])
* 2 * numbers::PI);
const unsigned int n_frequencies = Utilities::string_to_int(min_max_steps[2]);
const double delta_log_omega = (log_max_omega - log_min_omega)
/ (n_frequencies - 1)
* (1.-1e-12);
unsigned int index = 0;
for (double log_omega = log_min_omega;
log_omega <= log_max_omega;
log_omega += delta_log_omega, ++index)
MaterialParameters::frequencies.emplace_back (std::exp(log_omega),index);
}
else if (frequency_descriptor.find ("list") == 0)
{
// Get the rest of the string, and eat any space at the start and end
const std::string parenthesized_expr
= Utilities::trim (frequency_descriptor.substr
(std::string("list").size(),
std::string::npos));
AssertThrow (parenthesized_expr.size() >= 2
&&
parenthesized_expr.front() == '('
&&
parenthesized_expr.back() == ')',
ExcMessage ("Wrong format for 'list' frequency spacing."));
// Then get the interior part, again trim spaces, and split at
// commas
const std::vector<double> freq_list =
Utilities::string_to_double
(Utilities::split_string_list
(Utilities::trim (parenthesized_expr.substr
(1,
parenthesized_expr.size() - 2)),
','));
AssertThrow (freq_list.size() >= 1,
ExcMessage ("Wrong format for 'list' frequency spacing."));
for (unsigned int i=0; i<freq_list.size(); ++i)
MaterialParameters::frequencies.emplace_back (freq_list[i], i);
// Because MaterialParameters::frequencies stores angular
// frequencies, we need to multiply by 2*pi
for (auto &f : MaterialParameters::frequencies)
f.first *= 2 * numbers::PI;
}
else
AssertThrow (false,
ExcMessage ("The format for the description of the frequencies to "
"be solved for, namely <"
+ frequency_descriptor + ">, did not match any of "
"the recognized formats."));
fe_degree = prm.get_integer ("Finite element polynomial degree");
n_mesh_refinement_steps = prm.get_integer ("Number of mesh refinement steps");
n_threads = prm.get_integer ("Number of threads");
mesh_summary_only = prm.get_bool ("Mesh summary only");
}
// A data structure that is used to collect the results of the computations
// for one frequency. The main class fills this for a given frequency
// in various places of its member functions, and at the end puts it into
// a global map.
struct OutputData
{
// For each source port (first index), store integrated pressures
// and velocities at all other ports (second index)
FullMatrix<ScalarType> P, U;
// For each source port (first index), store pressure and velocity
// at all evaluation points (second index)
Table<2,std::complex<double>> evaluation_point_pressures;
Table<2,Tensor<1,3,std::complex<double>>> evaluation_point_velocities;
std::vector<std::string> visualization_file_names;
};
TimerOutput timer_output = TimerOutput (logger, TimerOutput::summary,
TimerOutput::wall_times);
// Create a file that indicates that program execution failed for
// some reason. Then exit the program with an error code.
void create_failure_file_and_exit ()
{
{
std::ofstream failure_signal (instance_folder + dir_sep_d +
output_file_prefix +
dii_solver_failure_signal_st);
failure_signal << "FAILURE" << std::endl;
}
std::exit(1);
}
// Check whether an external program has left a signal that
// indicates that the current program run should terminate without
// computing any further frequency responses. This is done by
// placing the word "STOP" into the file "termination_signal" in the
// current directory.
//
// Once detected, we delete the file again and terminate the
// program.
void check_for_termination_signal()
{
// Try and see whether we can open the file at all. If we can't,
// then no termination signal has been sent. If so, return 'true',
// but before that set a flag that ensures we don't have to do the
// expensive test with the file in any further calls. (We'll try
// to abort the program below, but this may block for a bit
// because we need to wait for the lock that guards access to the
// output file.)
std::ifstream in(instance_folder + dir_sep_d +
output_file_prefix +
dii_termination_signal_filename_st);
if (!in)
return;
// OK, the file exists, but does it contain the right content?
std::string line;
std::getline(in, line);
if (line == "STOP")
{
// Close the file handle and remove the file.
in.close();
std::remove ((instance_folder + dir_sep_d +
output_file_prefix +
dii_termination_signal_filename_st).c_str());
logger << "INFO *** Terminating program upon request." << std::endl;
create_failure_file_and_exit();
}
// The file exists, but it has the wrong content (or no content so
// far). This means no termination. In the best of all cases, we
// will have caught the driver program having created but not
// written to the file. The next time we check, we might find the
// file in the correct state.
}
/**
* Provide a standardized way of writing complex numbers. Specifically,
* write the number `p` in the form `a+bj`. If `field_width` is a positive
* number, then both the real and imaginary parts are given this many
* characters to occupy (excluding the character 'j'), where the real
* part is right-aligned and the imaginary part left-aligned. If
* `field_width` is zero, no alignment happens.
*/
void write_complex_number (const std::complex<double> &p,
const unsigned int field_width,
std::ostream &out)
{
if (field_width > 0)
{
out << std::setw(field_width) << std::right << std::real(p)
<< (std::imag(p) >= 0 ? '+' : '-');
std::ostringstream s;
s << std::fabs(std::imag(p)) << "*j";
out << std::setw(field_width+1) << std::left << s.str();
}
else
out << std::real(p)
<< (std::imag(p) >= 0 ? '+' : '-')
<< std::fabs(std::imag(p))
<< "*j";
}
// @sect3{The main class}
//
// The following is the principal class of this tutorial program. It has
// the structure of many of the other tutorial programs and there should
// really be nothing particularly surprising about its contents or
// the constructor that follows it.
template <int dim>
class HelmholtzProblem
{
public:
HelmholtzProblem(const double omega, const unsigned int frequency_number);
void run();
private:
void make_grid();
void setup_system();
void assemble_system(const unsigned int current_source_port);
void solve();
void postprocess(const unsigned int current_source_port);
void output_results(const unsigned int current_source_port);
void output_evaluated_data ();
// The frequency that this instance of the class is supposed to
// solve for, its number in the list, along with the density and
// wave speed to be used for this frequency.
const double omega;
const unsigned int frequency_number;
const std::complex<double> density;
const std::complex<double> wave_speed;
Triangulation<dim> triangulation;
std::vector<types::boundary_id> port_boundary_ids;
std::vector<double> port_areas;
std::unique_ptr<Mapping<dim>> mapping;
Quadrature<dim> cell_quadrature;
Quadrature<dim-1> face_quadrature;
std::unique_ptr<FiniteElement<dim>> fe;
DoFHandler<dim> dof_handler;
SparsityPattern sparsity_pattern;
SparseMatrix<ScalarType> system_matrix;
Vector<ScalarType> solution;
Vector<ScalarType> system_rhs;
OutputData output_data;
};
template <int dim>
HelmholtzProblem<dim>::HelmholtzProblem(const double omega,
const unsigned int frequency_number)
: omega (omega)
, frequency_number(frequency_number)
, density (MaterialParameters::density_real->value(Point<1>(omega)),
MaterialParameters::density_imag->value(Point<1>(omega)))
, wave_speed (std::sqrt(std::complex<double>(MaterialParameters::bulk_modulus_real->value(Point<1>(omega)),
MaterialParameters::bulk_modulus_imag->value(Point<1>(omega)))
/ density))
, dof_handler(triangulation)
{
logger << "INFO omega=" << omega << std::endl;
AssertThrow (std::real(1./wave_speed)>0
&&
std::imag(1./wave_speed)<=0,
ExcMessage("You need to choose bulk moduli and densities "
"that lead to wave numberss for which the real "
"part is positive and the imaginary part is negative or zero. "
"\n\n"
"The values you have chosen resulted in a wave speed of "
+ std::to_string(std::real(1./wave_speed))
+ std::to_string(std::imag(1./wave_speed))
+"*j if the frequency is omega=1."));
}
// Next up are the functions that create the initial mesh (a once refined
// unit square) and set up the constraints, vectors, and matrices on
// each mesh. Again, both of these are essentially unchanged from many
// previous tutorial programs.
template <int dim>
void HelmholtzProblem<dim>::make_grid()
{
std::unique_ptr<TimerOutput::Scope> timer_section = (n_threads==1 ? std::make_unique<TimerOutput::Scope>(timer_output, "Make grid") : nullptr);
// Read in the file in question. Do this on only one thread at a
// time to avoid file access collisions
{
static std::mutex mutex;
std::lock_guard<std::mutex> lock(mutex);
static std::unique_ptr<Triangulation<dim>> mesh_from_file;
static std::once_flag read_mesh_file;
std::call_once (read_mesh_file,
[this]()
{
std::unique_ptr<TimerOutput::Scope> timer_section = (n_threads==1 ? std::make_unique<TimerOutput::Scope>(timer_output, "Read mesh file") : nullptr);
mesh_from_file = std::make_unique<Triangulation<dim>>();
std::ifstream input (instance_folder + "/" + mesh_file_name);
AssertThrow (input,
ExcMessage ("The file <" + instance_folder + "/" + mesh_file_name +
"> can not be read from when "
"trying to load a mesh."));
// Determine what format we want to read the mesh in: .mphtxt =>
// COMSOL; .msh => GMSH; .inp => ABAQUS
try
{
GridIn<dim> grid_in;
grid_in.attach_triangulation (*mesh_from_file);
if (std::regex_match(mesh_file_name,
std::regex(".*\\.mphtxt", std::regex_constants::basic)))
{
logger << "INFO Reading mesh file <" << mesh_file_name
<< "> in COMSOL .mphtxt format" << std::endl;
grid_in.read_comsol_mphtxt (input);
}
else if (std::regex_match(mesh_file_name,
std::regex(".*\\.msh", std::regex_constants::basic)))
{
logger << "INFO Reading mesh file <" << mesh_file_name
<< "> in GMSH .msh format" << std::endl;
grid_in.read_msh (input);
}
else if (std::regex_match(mesh_file_name,
std::regex(".*\\.inp", std::regex_constants::basic)))
{
logger << "INFO Reading mesh file <" << mesh_file_name
<< "> in ABAQUS .inp format" << std::endl;
grid_in.read_abaqus (input);
}
else
AssertThrow (false,
ExcMessage ("The file ending for the mesh file <"
+ mesh_file_name +
"> is not supported."));
}
catch (const ExcIO &)
{
AssertThrow (false,
ExcMessage ("Couldn't read from mesh file <"
+ instance_folder + "/" + mesh_file_name
+ ">."));
}
});
logger << "INFO Cloning the mesh" << std::endl;
triangulation.copy_triangulation (*mesh_from_file);
}
logger << "INFO The starting mesh has " << triangulation.n_active_cells() << " cells" << std::endl;
// Scale the triangulation by the geometry factor
GridTools::scale (geometry_conversion_factor, triangulation);
// Now implement the heuristic for mesh refinement described in
// readme.md: If positive, just do a number of global refinement
// steps. If negative, interpret it as the number of mesh points
// per wave length.
if (n_mesh_refinement_steps >= 0)
triangulation.refine_global(n_mesh_refinement_steps);
else
{
const int N = -n_mesh_refinement_steps;
std::complex<double> k = omega / wave_speed;
// We want a mesh size that can has N mesh points per wave
// length of the oscillatory part (which is 2*pi/kr) with the
// real part of k, but that can also resolve the exponential
// decay due to the imaginary part of k with N points. So the
// length scale L we need to resolve is as follows (with the
// additional modification that it cannot be larger than the
// size of the domain).
//
// If we have a lossless medium (imag(k)=0), then the
// attenuation lengthscale is infinite.
const double L = std::min ({2*numbers::PI/(std::real(k)),
(imag(k) != 0 ?
-1./(std::imag(k))
:
std::numeric_limits<double>::infinity()),
GridTools::diameter (triangulation)});
AssertThrow (L>0,
ExcMessage ("You are using material parameters for the current "
"frequency (omega=" +
std::to_string(omega) +
") for which one of the length scales used "
"for mesh refinement is negative. In particular, "
"the wave length is 2*pi/real(k)=" +
std::to_string(2*numbers::PI/(std::real(k))) +
" and the attenuation length scale is -1/imag(k)=" +
std::to_string(-1./(std::imag(k))) +
"."));
// To have N mesh points, we need a delta_x that divides by
// (N/fe_degree), given that that there are delta_x/fe_degree
// grid points in each cell:
const double delta_x = L/N*fe_degree;
while (GridTools::maximal_cell_diameter(triangulation,
triangulation.get_reference_cells()[0].template get_default_linear_mapping<dim>())
>= delta_x)
triangulation.refine_global();
}
logger << "INFO The refined mesh has " << triangulation.n_active_cells() << " cells" << std::endl;
// Depending on what reference cell the triangulation uses, pick
// the right finite element
Assert (triangulation.get_reference_cells().size() == 1,
ExcMessage("We don't know how to deal with this kind of mesh."));
if (triangulation.get_reference_cells()[0] == ReferenceCells::Tetrahedron)
{
fe = std::make_unique<FE_SimplexP<dim>>(TransmissionProblem::fe_degree);
mapping = ReferenceCells::Tetrahedron.get_default_mapping<dim,dim>(1);
}
else if (triangulation.get_reference_cells()[0] == ReferenceCells::Hexahedron)
{
fe = std::make_unique<FE_Q<dim>>(TransmissionProblem::fe_degree);
mapping = ReferenceCells::Hexahedron.get_default_mapping<dim,dim>(1);
}
else
Assert (false, ExcMessage("We don't know how to deal with this kind of mesh."));
cell_quadrature = triangulation.get_reference_cells()[0]
.template get_gauss_type_quadrature<dim>(TransmissionProblem::fe_degree+1);
face_quadrature = triangulation.get_reference_cells()[0].face_reference_cell(0)
.template get_gauss_type_quadrature<dim-1>(TransmissionProblem::fe_degree+1);
// Finally output the mesh with ports correctly colored. This has
// to be done only once, so guard everything appropriately.
static std::once_flag output_boundary_visualization;
std::call_once (output_boundary_visualization,
[this]()
{
class BoundaryIds : public DataPostprocessorScalar<dim>
{
public:
BoundaryIds()
: DataPostprocessorScalar<dim>("boundary_id", update_quadrature_points)
{}
virtual void evaluate_scalar_field(
const DataPostprocessorInputs::Scalar<dim> &inputs,
std::vector<Vector<double>> &computed_quantities) const override
{
AssertDimension(computed_quantities.size(),
inputs.solution_values.size());
const typename DoFHandler<dim>::active_cell_iterator
cell = inputs.template get_cell<dim>();
const unsigned int face = inputs.get_face_number();
for (auto &output : computed_quantities)
{
AssertDimension(output.size(), 1);
output(0) = cell->face(face)->boundary_id();
}
}
};
BoundaryIds boundary_ids;
DataOutFaces<dim> data_out_faces;
std::unique_ptr<FiniteElement<dim>> fe;
Assert (triangulation.get_reference_cells().size() == 1,
ExcMessage("We don't know how to deal with this kind of mesh."));
if (triangulation.get_reference_cells()[0] == ReferenceCells::Tetrahedron)
fe = std::make_unique<FE_SimplexP<dim>>(TransmissionProblem::fe_degree);
else if (triangulation.get_reference_cells()[0] == ReferenceCells::Hexahedron)
fe = std::make_unique<FE_Q<dim>>(TransmissionProblem::fe_degree);
else
Assert (false, ExcMessage("We don't know how to deal with this kind of mesh."));
DoFHandler<dim> dummy_dof_handler(triangulation);
dummy_dof_handler.distribute_dofs(*fe);
Vector<double> dummy_solution (dummy_dof_handler.n_dofs());
data_out_faces.attach_dof_handler(dummy_dof_handler);
data_out_faces.add_data_vector(dummy_solution, boundary_ids);
data_out_faces.build_patches();
const std::string file_name = instance_folder + dir_sep_d +
output_file_prefix +
dii_visualization_dir_name_st + dir_sep_d +
"surface.vtu";
std::ofstream out(file_name);
AssertThrow (out,
ExcMessage ("The file <" + file_name +
"> can not be written to when trying to write "
"surface visualization data."));
try
{
data_out_faces.write_vtu(out);
}
catch (const ExcIO &)
{
AssertThrow (false,
ExcMessage ("Couldn't write to surface visualization output file <"
+ file_name
+ ">."));
}
});
// Figure out what boundary ids we have that describe ports. We
// take these as all of those boundary ids that are non-zero
port_boundary_ids = triangulation.get_boundary_ids();
static std::once_flag output_port_ids;
std::call_once (output_port_ids,
[this]()
{
logger << "INFO Found boundary ids ";
for (const auto &id : port_boundary_ids)
logger << id << ' ';
logger << std::endl;
});
port_boundary_ids.erase (std::find(port_boundary_ids.begin(),
port_boundary_ids.end(),
types::boundary_id(0)));
AssertThrow (port_boundary_ids.size() > 0,
ExcMessage ("The geometry read in from the mesh file has no "
"parts of the boundary marked with boundary indicators "
"other than zero. This means that it has no parts of "
"the boundary specifically marked as part of a 'port'. "
"We don't know what to compute then."));
// Now also correctly size the matrices we compute for each
// frequency
output_data.P.reinit (port_boundary_ids.size(),
port_boundary_ids.size());
output_data.U.reinit (port_boundary_ids.size(),
port_boundary_ids.size());
output_data.evaluation_point_pressures.reinit (port_boundary_ids.size(),
evaluation_points.size());
output_data.evaluation_point_velocities.reinit (port_boundary_ids.size(),
evaluation_points.size());
// As a final step, compute the areas of the various ports so we
// can later normalize when computing average pressures and
// velocities:
port_areas.resize (port_boundary_ids.size(), 0.);
FEFaceValues<dim> fe_face_values(*mapping,
*fe,
face_quadrature,
update_JxW_values);
for (const auto &cell : triangulation.active_cell_iterators())
for (const auto &face : cell->face_iterators())
if (face->at_boundary())
{
const auto boundary_id = face->boundary_id();
Assert ((boundary_id == 0) ||
(std::find(port_boundary_ids.begin(),
port_boundary_ids.end(),
face->boundary_id()) !=
port_boundary_ids.end()),
ExcInternalError());
if (boundary_id != 0)
{
const unsigned int this_port = (std::find(port_boundary_ids.begin(),
port_boundary_ids.end(),
boundary_id)
- port_boundary_ids.begin());
fe_face_values.reinit(cell, face);
for (const unsigned int q_point : fe_face_values.quadrature_point_indices())
port_areas[this_port] += fe_face_values.JxW(q_point);
}
}
// Output the port areas, but again do so only once
static std::once_flag output_port_areas;
std::call_once (output_port_areas,
[this]()
{
const std::string file_name = instance_folder + dir_sep_d +
output_file_prefix +
dii_port_areas_base_st;
std::ofstream port_area_output(file_name);
AssertThrow (port_area_output,
ExcMessage ("The file <" + file_name +
"> can not be written to when trying to write "
"port area data."));
for (const types::boundary_id b_id : port_boundary_ids)
{
const unsigned int this_port = (std::find(port_boundary_ids.begin(),
port_boundary_ids.end(),
b_id)
- port_boundary_ids.begin());
port_area_output << b_id << ' ' << port_areas[this_port] << std::endl;
}
});
}
template <int dim>
void HelmholtzProblem<dim>::setup_system()
{
std::unique_ptr<TimerOutput::Scope> timer_section = (n_threads==1 ? std::make_unique<TimerOutput::Scope>(timer_output, "Set up system") : nullptr);
dof_handler.distribute_dofs(*fe);
logger << "INFO The mesh has " << dof_handler.n_dofs() << " unknowns" << std::endl;
DynamicSparsityPattern c_sparsity(dof_handler.n_dofs());
DoFTools::make_sparsity_pattern (dof_handler, c_sparsity);