/
ThermalPhysics.templates.hh
1268 lines (1162 loc) · 48.5 KB
/
ThermalPhysics.templates.hh
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) 2016 - 2024, the adamantine authors.
*
* This file is subject to the Modified BSD License and may not be distributed
* without copyright and license information. Please refer to the file LICENSE
* for the text and further information on this license.
*/
#ifndef THERMAL_PHYSICS_TEMPLATES_HH
#define THERMAL_PHYSICS_TEMPLATES_HH
#include <CubeHeatSource.hh>
#include <ElectronBeamHeatSource.hh>
#include <GoldakHeatSource.hh>
#include <ThermalOperator.hh>
#include <ThermalOperatorDevice.hh>
#include <ThermalPhysics.hh>
#include <Timer.hh>
#include <deal.II/base/geometry_info.h>
#include <deal.II/base/index_set.h>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/distributed/cell_data_transfer.templates.h>
#include <deal.II/distributed/solution_transfer.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/fe/fe_nothing.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/grid/filtered_iterator.h>
#include <deal.II/hp/fe_values.h>
#include <deal.II/hp/q_collection.h>
#include <deal.II/lac/precondition.h>
#include <deal.II/lac/read_write_vector.h>
#include <deal.II/lac/solver_gmres.h>
#include <deal.II/lac/vector_operation.h>
#ifdef ADAMANTINE_WITH_CALIPER
#include <caliper/cali.h>
#endif
#include <algorithm>
#include <memory>
namespace adamantine
{
namespace
{
template <int dim, int fe_degree, typename MemorySpaceType,
std::enable_if_t<
std::is_same<MemorySpaceType, dealii::MemorySpace::Host>::value,
int> = 0>
dealii::LA::distributed::Vector<double, MemorySpaceType>
evaluate_thermal_physics_impl(
std::shared_ptr<ThermalOperatorBase<dim, MemorySpaceType>> thermal_operator,
double const t, double const current_source_height,
dealii::LA::distributed::Vector<double, MemorySpaceType> const &y,
std::vector<Timer> &timers)
{
timers[evol_time_eval_th_ph].start();
thermal_operator->set_time_and_source_height(t, current_source_height);
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host> value(
y.get_partitioner());
value = 0.;
// Apply the Thermal Operator.
thermal_operator->vmult_add(value, y);
// Multiply by the inverse of the mass matrix.
value.scale(*thermal_operator->get_inverse_mass_matrix());
timers[evol_time_eval_th_ph].stop();
return value;
}
template <int dim, int fe_degree, typename MemorySpaceType,
std::enable_if_t<
std::is_same<MemorySpaceType, dealii::MemorySpace::Host>::value,
int> = 0>
void init_dof_vector(
double const value,
dealii::LinearAlgebra::distributed::Vector<double, MemorySpaceType> &vector)
{
unsigned int const local_size = vector.locally_owned_size();
for (unsigned int i = 0; i < local_size; ++i)
vector.local_element(i) = value;
}
template <int dim, bool use_table, int p_order, int fe_degree,
typename MaterialStates, typename MemorySpaceType,
std::enable_if_t<std::is_same<MemorySpaceType,
dealii::MemorySpace::Default>::value,
int> = 0>
dealii::LA::distributed::Vector<double, MemorySpaceType>
evaluate_thermal_physics_impl(
std::shared_ptr<ThermalOperatorBase<dim, MemorySpaceType>> const
&thermal_operator,
dealii::hp::FECollection<dim> const &fe_collection, double const t,
dealii::DoFHandler<dim> const &dof_handler,
std::vector<std::shared_ptr<HeatSource<dim>>> const &heat_sources,
double current_source_height, BoundaryType boundary_type,
MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>
&material_properties,
dealii::AffineConstraints<double> const &affine_constraints,
dealii::LA::distributed::Vector<double, MemorySpaceType> const &y,
std::vector<Timer> &timers)
{
auto thermal_operator_dev = std::dynamic_pointer_cast<ThermalOperatorDevice<
dim, use_table, p_order, fe_degree, MaterialStates, MemorySpaceType>>(
thermal_operator);
timers[evol_time_update_bound_mat_prop].start();
thermal_operator_dev->update_boundary_material_properties(y);
timers[evol_time_update_bound_mat_prop].stop();
timers[evol_time_eval_th_ph].start();
dealii::LA::distributed::Vector<double, MemorySpaceType> value_dev(
y.get_partitioner());
// Apply the Thermal Operator.
thermal_operator_dev->vmult(value_dev, y);
// Compute the source term.
// TODO do this on the GPU
for (auto &beam : heat_sources)
beam->update_time(t);
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host> source(
y.get_partitioner());
source = 0.;
// Compute inv_rho_cp at the cell level on the host. We would not need to do
// this if everything was done on the GPU.
thermal_operator_dev->update_inv_rho_cp_cell();
dealii::hp::QCollection<dim> source_q_collection;
source_q_collection.push_back(dealii::QGauss<dim>(fe_degree + 1));
source_q_collection.push_back(dealii::QGauss<dim>(1));
dealii::hp::FEValues<dim> hp_fe_values(fe_collection, source_q_collection,
dealii::update_quadrature_points |
dealii::update_values |
dealii::update_JxW_values);
unsigned int const dofs_per_cell = fe_collection.max_dofs_per_cell();
unsigned int const n_q_points = source_q_collection.max_n_quadrature_points();
std::vector<dealii::types::global_dof_index> local_dof_indices(dofs_per_cell);
dealii::QGauss<dim - 1> face_quadrature(fe_degree + 1);
dealii::FEFaceValues<dim> fe_face_values(
fe_collection[0], face_quadrature,
dealii::update_values | dealii::update_quadrature_points |
dealii::update_JxW_values);
unsigned int const n_face_q_points = face_quadrature.size();
dealii::Vector<double> cell_source(dofs_per_cell);
// Loop over the locally owned cells with an active FE index of zero
for (auto const &cell : dealii::filter_iterators(
dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell(),
dealii::IteratorFilters::ActiveFEIndexEqualTo(0)))
{
cell_source = 0.;
hp_fe_values.reinit(cell);
dealii::FEValues<dim> const &fe_values =
hp_fe_values.get_present_fe_values();
for (unsigned int i = 0; i < dofs_per_cell; ++i)
{
for (unsigned int q = 0; q < n_q_points; ++q)
{
double const inv_rho_cp = thermal_operator_dev->get_inv_rho_cp(cell, q);
double quad_pt_source = 0.;
dealii::Point<dim> const &q_point = fe_values.quadrature_point(q);
for (auto &beam : heat_sources)
quad_pt_source += beam->value(q_point, current_source_height);
cell_source[i] += inv_rho_cp * quad_pt_source *
fe_values.shape_value(i, q) * fe_values.JxW(q);
}
}
// If we don't have a adiabatic boundary conditions, we need to add boundary
// conditions
if (!(boundary_type & BoundaryType::adiabatic))
{
for (unsigned int f = 0; f < dealii::GeometryInfo<dim>::faces_per_cell;
++f)
{
// We need to add the boundary conditions on the faces on the boundary
// but also on the faces at the interface with FE_Nothing
auto const &face = cell->face(f);
if ((face->at_boundary()) &&
((!face->at_boundary()) &&
(cell->neighbor(f)->active_fe_index() != 0)))
{
double conv_temperature_infty = 0.;
double conv_heat_transfer_coef = 0.;
double rad_temperature_infty = 0.;
double rad_heat_transfer_coef = 0.;
if (boundary_type & BoundaryType::convective)
{
conv_temperature_infty = material_properties.get_cell_value(
cell, Property::convection_temperature_infty);
conv_heat_transfer_coef = material_properties.get_cell_value(
cell, StateProperty::convection_heat_transfer_coef);
}
if (boundary_type & BoundaryType::radiative)
{
rad_temperature_infty = material_properties.get_cell_value(
cell, Property::radiation_temperature_infty);
rad_heat_transfer_coef = material_properties.get_cell_value(
cell, StateProperty::radiation_heat_transfer_coef);
}
fe_face_values.reinit(cell, face);
for (unsigned int i = 0; i < dofs_per_cell; ++i)
{
for (unsigned int q = 0; q < n_face_q_points; ++q)
{
double const inv_rho_cp =
thermal_operator_dev->get_inv_rho_cp(cell, q);
cell_source[i] +=
inv_rho_cp *
(conv_heat_transfer_coef * conv_temperature_infty +
rad_heat_transfer_coef * rad_temperature_infty) *
fe_face_values.shape_value(i, q) * fe_face_values.JxW(q);
}
}
}
}
}
cell->get_dof_indices(local_dof_indices);
affine_constraints.distribute_local_to_global(cell_source,
local_dof_indices, source);
}
source.compress(dealii::VectorOperation::add);
// Add source
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Default>
source_dev(source.get_partitioner());
source_dev.import(source, dealii::VectorOperation::insert);
value_dev += source_dev;
// Multiply by the inverse of the mass matrix.
value_dev.scale(*thermal_operator_dev->get_inverse_mass_matrix());
timers[evol_time_eval_th_ph].stop();
return value_dev;
}
template <int dim, int fe_degree, typename MemorySpaceType,
std::enable_if_t<std::is_same<MemorySpaceType,
dealii::MemorySpace::Default>::value,
int> = 0>
void init_dof_vector(
double const value,
dealii::LinearAlgebra::distributed::Vector<double, MemorySpaceType> &vector)
{
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host>
vector_host(vector.get_partitioner());
unsigned int const local_size = vector_host.locally_owned_size();
for (unsigned int i = 0; i < local_size; ++i)
vector_host.local_element(i) = value;
vector.import(vector_host, dealii::VectorOperation::insert);
}
} // namespace
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::
ThermalPhysics(MPI_Comm const &communicator,
boost::property_tree::ptree const &database,
Geometry<dim> &geometry,
MaterialProperty<dim, p_order, MaterialStates,
MemorySpaceType> &material_properties)
: _boundary_type(BoundaryType::invalid), _geometry(geometry),
_dof_handler(_geometry.get_triangulation()),
_cell_weights(
_dof_handler,
dealii::parallel::CellWeights<dim>::ndofs_weighting({1, 1})),
_material_properties(material_properties)
{
// Create the FECollection
_fe_collection.push_back(dealii::FE_Q<dim>(fe_degree));
_fe_collection.push_back(dealii::FE_Nothing<dim>());
// Create the QCollection
_q_collection.push_back(QuadratureType(fe_degree + 1));
_q_collection.push_back(QuadratureType(fe_degree + 1));
// Create the heat sources
boost::property_tree::ptree const &source_database =
database.get_child("sources");
// PropertyTreeInput sources.n_beams
unsigned int const n_beams = source_database.get<unsigned int>("n_beams");
_heat_sources.resize(n_beams);
for (unsigned int i = 0; i < n_beams; ++i)
{
// PropertyTreeInput sources.beam_X.type
boost::property_tree::ptree const &beam_database =
source_database.get_child("beam_" + std::to_string(i));
std::string type = beam_database.get<std::string>("type");
if (type == "goldak")
{
_heat_sources[i] = std::make_shared<GoldakHeatSource<dim>>(beam_database);
}
else if (type == "electron_beam")
{
_heat_sources[i] =
std::make_shared<ElectronBeamHeatSource<dim>>(beam_database);
}
else if (type == "cube")
{
_heat_sources[i] = std::make_shared<CubeHeatSource<dim>>(beam_database);
}
else
{
ASSERT_THROW(false, "Error: Beam type '" +
beam_database.get<std::string>("type") +
"' not recognized.");
}
}
// Create the boundary condition type
// PropertyTreeInput boundary.type
std::string boundary_type_str = database.get<std::string>("boundary.type");
// Parse the string
size_t pos_str = 0;
std::string boundary;
std::string delimiter = ",";
auto parse_boundary_type = [&](std::string const &boundary)
{
if (boundary == "adiabatic")
{
_boundary_type = BoundaryType::adiabatic;
}
else
{
if (boundary == "radiative")
{
_boundary_type |= BoundaryType::radiative;
}
else if (boundary == "convective")
{
_boundary_type |= BoundaryType::convective;
}
else
{
ASSERT_THROW(false, "Unknown boundary type.");
}
}
};
while ((pos_str = boundary_type_str.find(delimiter)) != std::string::npos)
{
boundary = boundary_type_str.substr(0, pos_str);
parse_boundary_type(boundary);
boundary_type_str.erase(0, pos_str + delimiter.length());
}
parse_boundary_type(boundary_type_str);
// Create the thermal operator
if (std::is_same<MemorySpaceType, dealii::MemorySpace::Host>::value)
{
if (_material_properties.properties_use_table())
{
_thermal_operator =
std::make_shared<ThermalOperator<dim, true, p_order, fe_degree,
MaterialStates, MemorySpaceType>>(
communicator, _boundary_type, _material_properties,
_heat_sources);
}
else
{
_thermal_operator =
std::make_shared<ThermalOperator<dim, false, p_order, fe_degree,
MaterialStates, MemorySpaceType>>(
communicator, _boundary_type, _material_properties,
_heat_sources);
}
}
else
{
if (_material_properties.properties_use_table())
{
_thermal_operator = std::make_shared<ThermalOperatorDevice<
dim, true, p_order, fe_degree, MaterialStates, MemorySpaceType>>(
communicator, _boundary_type, _material_properties);
}
else
{
_thermal_operator = std::make_shared<ThermalOperatorDevice<
dim, false, p_order, fe_degree, MaterialStates, MemorySpaceType>>(
communicator, _boundary_type, _material_properties);
}
}
// Create the time stepping scheme
boost::property_tree::ptree const &time_stepping_database =
database.get_child("time_stepping");
// PropertyTreeInput time_stepping.method
std::string method = time_stepping_database.get<std::string>("method");
std::transform(method.begin(), method.end(), method.begin(),
[](unsigned char c) { return std::tolower(c); });
if (method.compare("forward_euler") == 0)
_time_stepping =
std::make_unique<dealii::TimeStepping::ExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::FORWARD_EULER);
else if (method.compare("rk_third_order") == 0)
_time_stepping =
std::make_unique<dealii::TimeStepping::ExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::RK_THIRD_ORDER);
else if (method.compare("rk_fourth_order") == 0)
_time_stepping =
std::make_unique<dealii::TimeStepping::ExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::RK_CLASSIC_FOURTH_ORDER);
else if (method.compare("heun_euler") == 0)
{
_time_stepping = std::make_unique<
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::HEUN_EULER);
_embedded_method = true;
}
else if (method.compare("bogacki_shampine") == 0)
{
_time_stepping = std::make_unique<
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::BOGACKI_SHAMPINE);
_embedded_method = true;
}
else if (method.compare("dopri") == 0)
{
_time_stepping = std::make_unique<
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::DOPRI);
_embedded_method = true;
}
else if (method.compare("fehlberg") == 0)
{
_time_stepping = std::make_unique<
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::FEHLBERG);
_embedded_method = true;
}
else if (method.compare("cash_karp") == 0)
{
_time_stepping = std::make_unique<
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::CASH_KARP);
_embedded_method = true;
}
else if (method.compare("backward_euler") == 0)
{
_time_stepping =
std::make_unique<dealii::TimeStepping::ImplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::BACKWARD_EULER);
_implicit_method = true;
}
else if (method.compare("implicit_midpoint") == 0)
{
_time_stepping =
std::make_unique<dealii::TimeStepping::ImplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::IMPLICIT_MIDPOINT);
_implicit_method = true;
}
else if (method.compare("crank_nicolson") == 0)
{
_time_stepping =
std::make_unique<dealii::TimeStepping::ImplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::CRANK_NICOLSON);
_implicit_method = true;
}
else if (method.compare("sdirk2") == 0)
{
_time_stepping =
std::make_unique<dealii::TimeStepping::ImplicitRungeKutta<LA_Vector>>(
dealii::TimeStepping::SDIRK_TWO_STAGES);
_implicit_method = true;
}
if (_embedded_method == true)
{
// PropertyTreeInput time_steppping.coarsening_parameter
double coarsen_param =
time_stepping_database.get("coarsening_parameter", 1.2);
// PropertyTreeInput time_steppping.refining_parameter
double refine_param = time_stepping_database.get("refining_parameter", 0.8);
// PropertyTreeInput time_stepping.min_time_step
double min_delta = time_stepping_database.get("min_time_step", 1e-14);
// PropertyTreeInput time_stepping.max_time_step
double max_delta = time_stepping_database.get("max_time_step", 1e100);
// PropertyTreeInput time_stepping.refining_tolerance
double refine_tol = time_stepping_database.get("refining_tolerance", 1e-8);
// PropertyTreeInput time_stepping.coarsening_tolerance
double coarsen_tol =
time_stepping_database.get("coarsening_tolerance", 1e-12);
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector> *embedded_rk =
static_cast<
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector> *>(
_time_stepping.get());
embedded_rk->set_time_adaptation_parameters(coarsen_param, refine_param,
min_delta, max_delta,
refine_tol, coarsen_tol);
}
// If the time stepping scheme is implicit, set the parameters for the solver
// and create the implicit operator.
if (_implicit_method == true)
{
// PropertyTreeInput time_stepping.max_iteration
_max_iter = time_stepping_database.get("max_iteration", 1000);
// PropertyTreeInput time_stepping.tolerance
_tolerance = time_stepping_database.get("tolerance", 1e-12);
// PropertyTreeInput time_stepping.right_preconditioning
_right_preconditioning =
time_stepping_database.get("right_preconditioning", false);
// PropertyTreeInput time_stepping.n_tmp_vectors
_max_n_tmp_vectors = time_stepping_database.get("n_tmp_vectors", 30);
// PropertyTreeInput time_stepping.newton_max_iteration
unsigned int newton_max_iter =
time_stepping_database.get("newton_max_iteration", 100);
// PropertyTreeInput time_stepping.newton_tolerance
double newton_tolerance =
time_stepping_database.get("newton_tolerance", 1e-6);
dealii::TimeStepping::ImplicitRungeKutta<LA_Vector> *implicit_rk =
static_cast<dealii::TimeStepping::ImplicitRungeKutta<LA_Vector> *>(
_time_stepping.get());
implicit_rk->set_newton_solver_parameters(newton_max_iter,
newton_tolerance);
// PropertyTreeInput time_stepping.jfnk
bool jfnk = time_stepping_database.get("jfnk", false);
_implicit_operator = std::make_unique<ImplicitOperator<MemorySpaceType>>(
_thermal_operator, jfnk);
}
// Set material on part of the domain
// PropertyTreeInput geometry.material_height
double const material_height = database.get("geometry.material_height", 1e9);
for (auto const &cell :
dealii::filter_iterators(_dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell()))
{
// If the center of the cell is below material_height, it contains material
// otherwise it does not.
if (cell->center()[axis<dim>::z] < material_height)
{
cell->set_active_fe_index(0);
// Set material deposition cos and sin. We arbitrarily choose cos = 1 and
// sin = 0
_deposition_cos.push_back(1.);
_deposition_sin.push_back(0.);
// Set the initial material as non-melted
_has_melted.push_back(false);
}
else
cell->set_active_fe_index(1);
}
// Set the initial height of the heat source. Right now this is just the
// maximum heat source height, which can lead to unexpected behavior for
// different sources with different heights.
double temp_height = std::numeric_limits<double>::lowest();
for (auto const &source : _heat_sources)
{
temp_height = std::max(temp_height, source->get_current_height(0.0));
}
_current_source_height = temp_height;
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::setup()
{
setup_dofs();
update_material_deposition_orientation();
compute_inverse_mass_matrix();
get_state_from_material_properties();
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::setup_dofs()
{
_dof_handler.distribute_dofs(_fe_collection);
dealii::IndexSet locally_relevant_dofs;
dealii::DoFTools::extract_locally_relevant_dofs(_dof_handler,
locally_relevant_dofs);
_affine_constraints.clear();
_affine_constraints.reinit(locally_relevant_dofs);
dealii::DoFTools::make_hanging_node_constraints(_dof_handler,
_affine_constraints);
_affine_constraints.close();
_thermal_operator->reinit(_dof_handler, _affine_constraints, _q_collection);
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::compute_inverse_mass_matrix()
{
_thermal_operator->compute_inverse_mass_matrix(_dof_handler,
_affine_constraints);
if (_implicit_method == true)
_implicit_operator->set_inverse_mass_matrix(
_thermal_operator->get_inverse_mass_matrix());
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::
mark_has_melted(
double const threshold_temperature,
dealii::LA::distributed::Vector<double, MemorySpaceType> &temperature)
{
temperature.update_ghost_values();
auto dofs_per_cell = _dof_handler.get_fe().dofs_per_cell;
dealii::hp::FEValues<dim> hp_fe_values(
_dof_handler.get_fe_collection(), _q_collection,
dealii::UpdateFlags::update_values |
dealii::UpdateFlags::update_JxW_values);
unsigned int const n_q_points = _q_collection.max_n_quadrature_points();
unsigned int cell_id = 0;
for (auto const &cell : dealii::filter_iterators(
_dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell(),
dealii::IteratorFilters::ActiveFEIndexEqualTo(0)))
{
if (!_has_melted[cell_id])
{
hp_fe_values.reinit(cell);
dealii::FEValues<dim> const &fe_values =
hp_fe_values.get_present_fe_values();
std::vector<dealii::types::global_dof_index> local_dof_indices(
fe_values.dofs_per_cell);
cell->get_dof_indices(local_dof_indices);
double cell_temperature = 0.0;
double cell_volume = 0.0;
for (unsigned int i = 0; i < dofs_per_cell; ++i)
{
for (unsigned int q = 0; q < n_q_points; ++q)
{
cell_temperature += fe_values.shape_value(i, q) *
temperature(local_dof_indices[i]) *
fe_values.JxW(q);
cell_volume += fe_values.shape_value(i, q) * fe_values.JxW(q);
}
}
cell_temperature /= cell_volume;
// Set the indicator that this cell has melted
if (cell_temperature > threshold_temperature)
{
_has_melted[cell_id] = true;
}
}
++cell_id;
}
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::
add_material(
std::vector<std::vector<
typename dealii::DoFHandler<dim>::active_cell_iterator>> const
&elements_to_activate,
std::vector<double> const &new_deposition_cos,
std::vector<double> const &new_deposition_sin,
std::vector<bool> &new_has_melted, unsigned int const activation_start,
unsigned int const activation_end,
double const new_material_temperature,
dealii::LA::distributed::Vector<double, MemorySpaceType> &solution)
{
#ifdef ADAMANTINE_WITH_CALIPER
CALI_CXX_MARK_FUNCTION;
#endif
// Update the material state from the ThermalOperator to MaterialProperty
// because, for now, we need to use state from MaterialProperty to perform the
// transfer to the refined mesh.
set_state_to_material_properties();
_thermal_operator->clear();
// The data on each cell is stored in the following order: solution, direction
// of deposition (cosine and sine), prior melting indictor, and state ratio.
std::vector<std::vector<double>> data_to_transfer;
unsigned int const n_dofs_per_cell = _dof_handler.get_fe().n_dofs_per_cell();
unsigned int const direction_data_size = 2;
unsigned int const phase_history_data_size = 1;
unsigned int constexpr n_material_states = MaterialStates::n_material_states;
unsigned int const data_size_per_cell =
n_dofs_per_cell + direction_data_size + phase_history_data_size +
n_material_states;
dealii::Vector<double> cell_solution(n_dofs_per_cell);
std::vector<double> dummy_cell_data(data_size_per_cell,
std::numeric_limits<double>::infinity());
solution.update_ghost_values();
// We need to move the solution on the host because we cannot use
// CellDataTransfer on the device.
dealii::IndexSet rw_index_set = solution.locally_owned_elements();
rw_index_set.add_indices(solution.get_partitioner()->ghost_indices());
dealii::LA::ReadWriteVector<double> rw_solution(rw_index_set);
rw_solution.import(solution, dealii::VectorOperation::insert);
auto state_host = Kokkos::create_mirror_view_and_copy(
Kokkos::HostSpace{}, _material_properties.get_state());
unsigned int locally_owned_cell_id = 0;
unsigned int activated_cell_id = 0;
unsigned int cell_id = 0;
std::map<typename dealii::DoFHandler<dim>::active_cell_iterator, int>
cell_to_id;
for (auto const &cell : _dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
if (cell->active_fe_index() == 0)
{
std::vector<double> cell_data(
direction_data_size + phase_history_data_size + n_material_states);
cell->get_dof_values(rw_solution, cell_solution);
cell_data.insert(cell_data.begin(), cell_solution.begin(),
cell_solution.end());
cell_data[n_dofs_per_cell] = _deposition_cos[activated_cell_id];
cell_data[n_dofs_per_cell + 1] = _deposition_sin[activated_cell_id];
if (_has_melted[activated_cell_id])
cell_data[n_dofs_per_cell + direction_data_size] = 1.0;
else
cell_data[n_dofs_per_cell + direction_data_size] = 0.0;
for (unsigned int i = 0; i < n_material_states; ++i)
cell_data[n_dofs_per_cell + direction_data_size +
phase_history_data_size + i] =
state_host(i, locally_owned_cell_id);
data_to_transfer.push_back(cell_data);
++activated_cell_id;
}
else
{
std::vector<double> cell_data = dummy_cell_data;
for (unsigned int i = 0; i < n_material_states; ++i)
cell_data[n_dofs_per_cell + direction_data_size +
phase_history_data_size + i] =
state_host(i, locally_owned_cell_id);
data_to_transfer.push_back(cell_data);
}
++locally_owned_cell_id;
}
else
{
data_to_transfer.push_back(dummy_cell_data);
}
cell_to_id[cell] = cell_id;
++cell_id;
}
// Activate elements by updating the fe_index
for (unsigned int i = activation_start; i < activation_end; ++i)
{
for (auto const &cell : elements_to_activate[i])
{
if (cell->active_fe_index() != 0)
{
cell->set_future_fe_index(0);
data_to_transfer[cell_to_id[cell]][n_dofs_per_cell] =
new_deposition_cos[i];
data_to_transfer[cell_to_id[cell]][n_dofs_per_cell + 1] =
new_deposition_sin[i];
if (data_to_transfer[cell_to_id[cell]]
[n_dofs_per_cell + direction_data_size] > 0.5)
new_has_melted[i] = true;
else
new_has_melted[i] = false;
}
}
}
dealii::parallel::distributed::Triangulation<dim> &triangulation =
dynamic_cast<dealii::parallel::distributed::Triangulation<dim> &>(
const_cast<dealii::Triangulation<dim> &>(
_dof_handler.get_triangulation()));
triangulation.prepare_coarsening_and_refinement();
dealii::parallel::distributed::CellDataTransfer<
dim, dim, std::vector<std::vector<double>>>
cell_data_trans(triangulation);
cell_data_trans.prepare_for_coarsening_and_refinement(data_to_transfer);
#ifdef ADAMANTINE_WITH_CALIPER
CALI_MARK_BEGIN("refine triangulation");
#endif
triangulation.execute_coarsening_and_refinement();
#ifdef ADAMANTINE_WITH_CALIPER
CALI_MARK_END("refine triangulation");
#endif
setup_dofs();
// Update MaterialProperty DoFHandler and resize the state vectors
_material_properties.reinit_dofs();
// Recompute the inverse of the mass matrix
compute_inverse_mass_matrix();
initialize_dof_vector(std::numeric_limits<double>::infinity(), solution);
rw_index_set = solution.locally_owned_elements();
rw_index_set.add_indices(solution.get_partitioner()->ghost_indices());
rw_solution.reinit(rw_index_set);
for (auto val : solution.locally_owned_elements())
rw_solution[val] = new_material_temperature;
// Unpack the material state and repopulate the material state
std::vector<std::vector<double>> transferred_data(
triangulation.n_active_cells(), std::vector<double>(data_size_per_cell));
cell_data_trans.unpack(transferred_data);
auto state = _material_properties.get_state();
state_host = Kokkos::create_mirror_view(state);
_deposition_cos.clear();
_deposition_sin.clear();
_has_melted.clear();
cell_id = 0;
locally_owned_cell_id = 0;
for (auto const &cell : _dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
if (transferred_data[cell_id][0] !=
std::numeric_limits<double>::infinity())
{
std::copy(transferred_data[cell_id].begin(),
transferred_data[cell_id].begin() + n_dofs_per_cell,
cell_solution.begin());
cell->set_dof_values(cell_solution, rw_solution);
}
if (cell->active_fe_index() == 0)
{
_deposition_cos.push_back(transferred_data[cell_id][n_dofs_per_cell]);
_deposition_sin.push_back(
transferred_data[cell_id][n_dofs_per_cell + 1]);
if (transferred_data[cell_id][n_dofs_per_cell + direction_data_size] >
0.5)
_has_melted.push_back(true);
else
_has_melted.push_back(false);
}
for (unsigned int i = 0; i < n_material_states; ++i)
{
state_host(i, locally_owned_cell_id) =
transferred_data[cell_id][n_dofs_per_cell + direction_data_size +
phase_history_data_size + i];
}
++locally_owned_cell_id;
}
++cell_id;
}
Kokkos::deep_copy(state, state_host);
get_state_from_material_properties();
_thermal_operator->set_material_deposition_orientation(_deposition_cos,
_deposition_sin);
// Communicate the results.
solution.import(rw_solution, dealii::VectorOperation::insert);
// Set the value to the newly create DoFs. Here we need to be careful with the
// hanging nodes. When there is a hanging node, the dofs at the vertices are
// "doubled": there is a dof associated to the coarse cell and a dof
// associated to the fine cell. The final value is decided by
// AffineConstraints. Thus, we need to make sure that the newly activated
// cells are at the same level than their neighbors.
rw_solution.reinit(solution.locally_owned_elements());
rw_solution.import(solution, dealii::VectorOperation::insert);
Kokkos::parallel_for("adamantine::set_new_material_temperature",
Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>(
0, rw_solution.locally_owned_size()),
[&](int i)
{
if (rw_solution.local_element(i) ==
std::numeric_limits<double>::infinity())
{
rw_solution.local_element(i) =
new_material_temperature;
}
});
solution.import(rw_solution, dealii::VectorOperation::insert);
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<
dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::update_physics_parameters(boost::property_tree::ptree const
&heat_source_database)
{
// Update the heat source from heat_source_database to reflect changes during
// the simulation (i.e. due to data assimilation)
unsigned int source_index = 0;
for (auto const &source : _heat_sources)
{
// PropertyTreeInput sources.beam_X
boost::property_tree::ptree const &beam_database =
heat_source_database.get_child("beam_" + std::to_string(source_index));
// PropertyTreeInput sources.beam_X.type
std::string type = beam_database.get<std::string>("type");
if (type == "goldak" || type == "electron_beam")
source->set_beam_properties(beam_database);
source_index++;
}
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
double ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::
evolve_one_time_step(
double t, double delta_t,
dealii::LA::distributed::Vector<double, MemorySpaceType> &solution,
std::vector<Timer> &timers)
{
// Update the height of the heat source. Right now this is just the
// maximum heat source height, which can lead to unexpected behavior for
// different sources with different heights.
double temp_height = std::numeric_limits<double>::lowest();
for (auto const &source : _heat_sources)
{
temp_height = std::max(temp_height, source->get_current_height(t));
}
_current_source_height = temp_height;
auto eval = [&](double const t, LA_Vector const &y)
{ return evaluate_thermal_physics(t, y, timers); };
auto id_m_Jinv = [&](double const t, double const tau, LA_Vector const &y)
{ return id_minus_tau_J_inverse(t, tau, y, timers); };
double time = _time_stepping->evolve_one_time_step(eval, id_m_Jinv, t,
delta_t, solution);
// If the method is embedded, get the next time step. Otherwise, just use the
// current time step.
if (_embedded_method == false)
_delta_t_guess = delta_t;
else
{
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector> *embedded_rk =
static_cast<
dealii::TimeStepping::EmbeddedExplicitRungeKutta<LA_Vector> *>(
_time_stepping.get());
_delta_t_guess = embedded_rk->get_status().delta_t_guess;
}
// Return the time at the end of the time step. This may be different than
// t+delta_t for embedded methods.
return time;
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::
initialize_dof_vector(
double const value,
dealii::LA::distributed::Vector<double, MemorySpaceType> &vector) const
{
// Resize the vector and initialize it to zero
_thermal_operator->initialize_dof_vector(vector);
if (value != 0.)
{
init_dof_vector<dim, fe_degree, MemorySpaceType>(value, vector);
}
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::get_state_from_material_properties()
{
_thermal_operator->get_state_from_material_properties();
}
template <int dim, int p_order, int fe_degree, typename MaterialStates,
typename MemorySpaceType, typename QuadratureType>
void ThermalPhysics<dim, p_order, fe_degree, MaterialStates, MemorySpaceType,
QuadratureType>::set_state_to_material_properties()
{
_thermal_operator->set_state_to_material_properties();
}