-
Notifications
You must be signed in to change notification settings - Fork 9
/
MaterialProperty.templates.hh
1139 lines (1068 loc) · 46.4 KB
/
MaterialProperty.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 MATERIAL_PROPERTY_TEMPLATES_HH
#define MATERIAL_PROPERTY_TEMPLATES_HH
#include <MaterialProperty.hh>
#include <deal.II/base/point.h>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/types.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/mapping.h>
#include <deal.II/grid/manifold.h>
#include <deal.II/hp/fe_values.h>
#include <deal.II/hp/q_collection.h>
#include <deal.II/lac/la_parallel_vector.h>
#include <boost/algorithm/string/split.hpp>
#include <boost/optional.hpp>
#include <Kokkos_Core_fwd.hpp>
#include <algorithm>
#include <type_traits>
namespace adamantine
{
namespace internal
{
template <int dim>
void compute_average(
unsigned int const n_q_points, unsigned int const dofs_per_cell,
dealii::DoFHandler<dim> const &mp_dof_handler,
dealii::DoFHandler<dim> const &temperature_dof_handler,
dealii::hp::FEValues<dim> &hp_fe_values,
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host> const
&temperature,
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host>
&temperature_average)
{
std::vector<dealii::types::global_dof_index> mp_dof_indices(1);
std::vector<dealii::types::global_dof_index> enth_dof_indices(dofs_per_cell);
auto mp_cell = mp_dof_handler.begin_active();
auto mp_end_cell = mp_dof_handler.end();
auto enth_cell = temperature_dof_handler.begin_active();
for (; mp_cell != mp_end_cell; ++enth_cell, ++mp_cell)
{
ASSERT(mp_cell->is_locally_owned() == enth_cell->is_locally_owned(),
"Internal Error");
if ((mp_cell->is_locally_owned()) && (enth_cell->active_fe_index() == 0))
{
hp_fe_values.reinit(enth_cell);
dealii::FEValues<dim> const &fe_values =
hp_fe_values.get_present_fe_values();
mp_cell->get_dof_indices(mp_dof_indices);
dealii::types::global_dof_index const mp_dof_index = mp_dof_indices[0];
enth_cell->get_dof_indices(enth_dof_indices);
double volume = 0.;
for (unsigned int q = 0; q < n_q_points; ++q)
for (unsigned int i = 0; i < dofs_per_cell; ++i)
{
volume += fe_values.shape_value(i, q) * fe_values.JxW(q);
temperature_average[mp_dof_index] +=
fe_values.shape_value(i, q) * temperature[enth_dof_indices[i]] *
fe_values.JxW(q);
}
temperature_average[mp_dof_index] /= volume;
}
}
}
template <typename ViewType,
std::enable_if_t<
std::is_same_v<typename ViewType::memory_space,
typename dealii::MemorySpace::Host::kokkos_space>,
int> = 0>
double get_value(ViewType &view, unsigned int i, unsigned int j)
{
return view(i, j);
}
template <int dim>
void compute_average(
unsigned int const n_q_points, unsigned int const dofs_per_cell,
dealii::DoFHandler<dim> const &mp_dof_handler,
dealii::DoFHandler<dim> const &temperature_dof_handler,
dealii::hp::FEValues<dim> &hp_fe_values,
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Default> const
&temperature,
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Default>
&temperature_average)
{
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host>
temperature_host(temperature.get_partitioner());
temperature_host.import(temperature, dealii::VectorOperation::insert);
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host>
temperature_average_host(temperature_average.get_partitioner());
temperature_average_host = 0.;
compute_average(n_q_points, dofs_per_cell, mp_dof_handler,
temperature_dof_handler, hp_fe_values, temperature_host,
temperature_average_host);
temperature_average.import(temperature_average_host,
dealii::VectorOperation::insert);
}
template <typename ViewType,
std::enable_if_t<
!std::is_same_v<typename ViewType::memory_space,
typename dealii::MemorySpace::Host::kokkos_space>,
int> = 0>
double get_value(ViewType &view, unsigned int i, unsigned int j)
{
auto view_host =
Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace{}, view);
return view_host(i, j);
}
} // namespace internal
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
MaterialProperty(
MPI_Comm const &communicator,
dealii::parallel::distributed::Triangulation<dim> const &tria,
boost::property_tree::ptree const &database)
: _communicator(communicator), _fe(0), _mp_dof_handler(tria)
{
// Because deal.II cannot easily attach data to a cell. We store the state
// of the material in distributed::Vector. This allows to use deal.II to
// compute the new state after refinement of the mesh. However, this
// requires to use another DoFHandler.
reinit_dofs();
// Set the material state to the state defined in the geometry.
set_initial_state();
// Fill the _properties map
fill_properties(database);
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
MaterialProperty(MaterialProperty<dim, p_order, MaterialStates,
MemorySpaceType> const &other)
: _use_table(other._use_table),
_state_property_tables(other._state_property_tables),
_state_property_polynomials(other._state_property_polynomials),
_properties(other._properties), _state(other._state),
_property_values(other._property_values), _fe(0)
{
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
double
MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::get_cell_value(
typename dealii::Triangulation<dim>::active_cell_iterator const &cell,
StateProperty prop) const
{
unsigned int property = static_cast<unsigned int>(prop);
auto const mp_dof_index = get_dof_index(cell);
// FIXME this is extremely slow on CUDA but this function should not exist in
// the first place
return internal::get_value(_property_values, property, mp_dof_index);
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
double
MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::get_cell_value(
typename dealii::Triangulation<dim>::active_cell_iterator const &cell,
Property prop) const
{
dealii::types::material_id material_id = cell->material_id();
unsigned int property = static_cast<unsigned int>(prop);
// FIXME this is extremely slow on CUDA but this function should not exist in
// the first place
return internal::get_value(_properties, material_id, property);
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
double MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
get_mechanical_property(
typename dealii::Triangulation<dim>::active_cell_iterator const &cell,
StateProperty prop) const
{
unsigned int property =
static_cast<unsigned int>(prop) - g_n_thermal_state_properties;
ASSERT(property < g_n_mechanical_state_properties,
"Unknown mechanical property requested.");
return _mechanical_properties_host(cell->material_id(), property);
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
double MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
get_state_ratio(
typename dealii::Triangulation<dim>::active_cell_iterator const &cell,
typename MaterialStates::State material_state) const
{
auto const mp_dof_index = get_dof_index(cell);
auto const mat_state = static_cast<unsigned int>(material_state);
// FIXME this is extremely slow on CUDA but this function should not exist in
// the first place
return internal::get_value(_state, mat_state, mp_dof_index);
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates,
MemorySpaceType>::reinit_dofs()
{
_mp_dof_handler.distribute_dofs(_fe);
// Initialize _dofs_map
_dofs_map.clear();
unsigned int i = 0;
std::vector<dealii::types::global_dof_index> mp_dof(1);
for (auto cell :
dealii::filter_iterators(_mp_dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell()))
{
cell->get_dof_indices(mp_dof);
_dofs_map[mp_dof[0]] = i;
++i;
}
_state = Kokkos::View<double **, typename MemorySpaceType::kokkos_space>(
"state", MaterialStates::n_material_states, _dofs_map.size());
#ifdef ADAMANTINE_DEBUG
if constexpr (std::is_same_v<MemorySpaceType, dealii::MemorySpace::Host>)
{
Kokkos::parallel_for(
"adamantine::set_state_nan",
Kokkos::MDRangePolicy<Kokkos::DefaultHostExecutionSpace,
Kokkos::Rank<2>>(
{{0, 0}}, {{MaterialStates::n_material_states,
static_cast<int>(_dofs_map.size())}}),
KOKKOS_CLASS_LAMBDA(int i, int j) {
_state(i, j) = std::numeric_limits<double>::signaling_NaN();
});
}
#endif
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::update(
dealii::DoFHandler<dim> const &temperature_dof_handler,
dealii::LA::distributed::Vector<double, MemorySpaceType> const &temperature)
{
auto temperature_average =
compute_average_temperature(temperature_dof_handler, temperature);
// Set View to zero in purpose
_property_values =
Kokkos::View<double **, typename MemorySpaceType::kokkos_space>(
"property_values", g_n_thermal_state_properties, _dofs_map.size());
std::vector<dealii::types::global_dof_index> mp_dofs_vec;
std::vector<dealii::types::material_id> material_ids_vec;
for (auto cell :
dealii::filter_iterators(_mp_dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell()))
{
std::vector<dealii::types::global_dof_index> mp_dof(1);
cell->get_dof_indices(mp_dof);
mp_dofs_vec.push_back(_dofs_map.at(mp_dof[0]));
material_ids_vec.push_back(cell->material_id());
}
unsigned int const material_ids_size = material_ids_vec.size();
Kokkos::View<dealii::types::material_id *, Kokkos::HostSpace> material_ids(
material_ids_vec.data(), material_ids_size);
auto material_ids_mirror = Kokkos::create_mirror_view_and_copy(
typename MemorySpaceType::kokkos_space{}, material_ids);
Kokkos::View<dealii::types::global_dof_index *, Kokkos::HostSpace> mp_dofs(
mp_dofs_vec.data(), mp_dofs_vec.size());
auto mp_dofs_mirror = Kokkos::create_mirror_view_and_copy(
typename MemorySpaceType::kokkos_space{}, mp_dofs);
double *temperature_average_local = temperature_average.get_values();
using ExecutionSpace = std::conditional_t<
std::is_same_v<MemorySpaceType, dealii::MemorySpace::Host>,
Kokkos::DefaultHostExecutionSpace, Kokkos::DefaultExecutionSpace>;
Kokkos::parallel_for(
"adamantine::update_material_properties",
Kokkos::RangePolicy<ExecutionSpace>(0, material_ids_size),
KOKKOS_CLASS_LAMBDA(int i) {
unsigned int constexpr solid =
static_cast<unsigned int>(MaterialStates::State::solid);
unsigned int constexpr prop_solidus =
static_cast<unsigned int>(Property::solidus);
unsigned int constexpr prop_liquidus =
static_cast<unsigned int>(Property::liquidus);
// Set solid_ratio to one, since that's the value when MaterialStates is
// Solid
double solid_ratio = 1.;
double liquid_ratio = 0.;
dealii::types::material_id material_id = material_ids_mirror(i);
double const solidus = _properties(material_id, prop_solidus);
double const liquidus = _properties(material_id, prop_liquidus);
unsigned int const dof = mp_dofs_mirror(i);
if constexpr (!std::is_same_v<MaterialStates, Solid>)
{
unsigned int constexpr liquid =
static_cast<unsigned int>(MaterialStates::State::liquid);
// First determine the ratio of liquid.
if (temperature_average_local[dof] < solidus)
liquid_ratio = 0.;
else if (temperature_average_local[dof] > liquidus)
liquid_ratio = 1.;
else
liquid_ratio = (temperature_average_local[dof] - solidus) /
(liquidus - solidus);
if constexpr (std::is_same_v<MaterialStates, SolidLiquid>)
{
solid_ratio = 1. - liquid_ratio;
}
else if constexpr (std::is_same_v<MaterialStates, SolidLiquidPowder>)
{
unsigned int constexpr powder =
static_cast<unsigned int>(MaterialStates::State::powder);
// Because the powder can only become liquid, the solid can only
// become liquid, and the liquid can only become solid, the ratio of
// powder can only decrease.
double powder_ratio =
Kokkos::min(1. - liquid_ratio, _state(powder, dof));
solid_ratio = 1. - liquid_ratio - powder_ratio;
// Update _state
_state(powder, dof) = powder_ratio;
}
_state(liquid, dof) = liquid_ratio;
}
// Update _state
_state(solid, dof) = solid_ratio;
if (_use_table)
{
for (unsigned int property = 0;
property < g_n_thermal_state_properties; ++property)
{
for (unsigned int material_state = 0;
material_state < MaterialStates::n_material_states;
++material_state)
{
_property_values(property, dof) +=
_state(material_state, dof) *
compute_property_from_table(
_state_property_tables, material_id, material_state,
property, temperature_average_local[dof]);
}
}
}
else
{
for (unsigned int property = 0;
property < g_n_thermal_state_properties; ++property)
{
for (unsigned int material_state = 0;
material_state < MaterialStates::n_material_states;
++material_state)
{
for (unsigned int i = 0; i <= p_order; ++i)
{
_property_values(property, dof) +=
_state(material_state, dof) *
_state_property_polynomials(material_id, material_state,
property, i) *
std::pow(temperature_average_local[dof], i);
}
}
}
}
// If we are in the mushy state, i.e., part liquid part solid, we need
// to modify the rho C_p to take into account the latent heat.
if constexpr (!std::is_same_v<MaterialStates, Solid>)
{
if ((liquid_ratio > 0.) && (liquid_ratio < 1.))
{
unsigned int const specific_heat_prop =
static_cast<unsigned int>(StateProperty::specific_heat);
unsigned int const latent_heat_prop =
static_cast<unsigned int>(Property::latent_heat);
for (unsigned int material_state = 0;
material_state < MaterialStates::n_material_states;
++material_state)
{
_property_values(specific_heat_prop, dof) +=
_state(material_state, dof) *
_properties(material_id, latent_heat_prop) /
(liquidus - solidus);
}
}
}
// The radiation heat transfer coefficient is not a real material
// property but it is derived from other material properties: h_rad =
// emissitivity * stefan-boltzmann constant * (T + T_infty) (T^2 +
// T^2_infty).
unsigned int const emissivity_prop =
static_cast<unsigned int>(StateProperty::emissivity);
unsigned int const radiation_heat_transfer_coef_prop =
static_cast<unsigned int>(
StateProperty::radiation_heat_transfer_coef);
unsigned int const radiation_temperature_infty_prop =
static_cast<unsigned int>(Property::radiation_temperature_infty);
double const T = temperature_average_local[dof];
double const T_infty =
_properties(material_id, radiation_temperature_infty_prop);
double const emissivity = _property_values(emissivity_prop, dof);
_property_values(radiation_heat_transfer_coef_prop, dof) =
emissivity * Constant::stefan_boltzmann * (T + T_infty) *
(T * T + T_infty * T_infty);
});
}
// TODO When we can get rid of this function, we can remove
// StateProperty::radiation_heat_transfer_coef
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
update_boundary_material_properties(
dealii::DoFHandler<dim> const &temperature_dof_handler,
dealii::LA::distributed::Vector<double, MemorySpaceType> const
&temperature)
{
auto temperature_average =
compute_average_temperature(temperature_dof_handler, temperature);
// Initialize the View to zero in purpose
_property_values =
Kokkos::View<double **, typename MemorySpaceType::kokkos_space>(
"property_values", g_n_thermal_state_properties, _dofs_map.size());
std::vector<dealii::types::global_dof_index> mp_dof(1);
// We don't need to loop over all the active cells. We only need to loop over
// the cells at the boundary and at the interface with FE_Nothing. However, to
// do this we need to use the temperature_dof_handler instead of the
// _mp_dof_handler.
for (auto cell :
dealii::filter_iterators(_mp_dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell()))
{
dealii::types::material_id material_id = cell->material_id();
cell->get_dof_indices(mp_dof);
unsigned int const dof = _dofs_map.at(mp_dof[0]);
if (_use_table)
{
// We only care about properties that are used to compute the boundary
// condition. So we start at 3.
for (unsigned int property = 3; property < g_n_thermal_state_properties;
++property)
{
for (unsigned int material_state = 0;
material_state < MaterialStates::n_material_states;
++material_state)
{
_property_values(property, dof) +=
_state(material_state, dof) *
compute_property_from_table(
_state_property_tables, material_id, material_state, property,
temperature_average.local_element(dof));
}
}
}
else
{
// We only care about properties that are used to compute the boundary
// condition. So we start at 3.
for (unsigned int property = 3; property < g_n_thermal_state_properties;
++property)
{
for (unsigned int material_state = 0;
material_state < MaterialStates::n_material_states;
++material_state)
{
for (unsigned int i = 0; i <= p_order; ++i)
{
_property_values(property, dof) +=
_state(material_state, dof) *
_state_property_polynomials(material_id, material_state,
property, i) *
std::pow(temperature_average.local_element(dof), i);
}
}
}
}
// The radiation heat transfer coefficient is not a real material property
// but it is derived from other material properties:
// h_rad = emissitivity * stefan-boltzmann constant * (T + T_infty) (T^2 +
// T^2_infty).
unsigned int const emissivity_prop =
static_cast<unsigned int>(StateProperty::emissivity);
unsigned int const radiation_heat_transfer_coef_prop =
static_cast<unsigned int>(StateProperty::radiation_heat_transfer_coef);
unsigned int const radiation_temperature_infty_prop =
static_cast<unsigned int>(Property::radiation_temperature_infty);
double const T = temperature_average.local_element(dof);
double const T_infty =
_properties(material_id, radiation_temperature_infty_prop);
double const emissivity = _property_values(emissivity_prop, dof);
_property_values(radiation_heat_transfer_coef_prop, dof) =
emissivity * Constant::stefan_boltzmann * (T + T_infty) *
(T * T + T_infty * T_infty);
}
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::set_state(
[[maybe_unused]] dealii::Table<2, dealii::VectorizedArray<double>> const
&liquid_ratio,
[[maybe_unused]] dealii::Table<2, dealii::VectorizedArray<double>> const
&powder_ratio,
[[maybe_unused]] std::map<typename dealii::DoFHandler<dim>::cell_iterator,
std::pair<unsigned int, unsigned int>>
&cell_it_to_mf_cell_map,
[[maybe_unused]] dealii::DoFHandler<dim> const &dof_handler)
{
if constexpr (std::is_same_v<MaterialStates, Solid>)
{
// When there is only Solid, we know we can set all of _state to one.
Kokkos::deep_copy(_state, 1.);
}
else
{
auto constexpr solid_state =
static_cast<unsigned int>(MaterialStates::State::solid);
auto constexpr liquid_state =
static_cast<unsigned int>(MaterialStates::State::liquid);
std::vector<dealii::types::global_dof_index> mp_dof(1.);
if constexpr (std::is_same_v<MaterialStates, SolidLiquid>)
{
for (auto const &cell : dealii::filter_iterators(
dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell()))
{
typename dealii::Triangulation<dim>::active_cell_iterator cell_tria(
cell);
auto mp_dof_index = get_dof_index(cell_tria);
auto const &mf_cell_vector = cell_it_to_mf_cell_map[cell];
unsigned int const n_q_points =
dof_handler.get_fe().tensor_degree() + 1;
double liquid_ratio_sum = 0.;
for (unsigned int q = 0; q < n_q_points; ++q)
{
liquid_ratio_sum +=
liquid_ratio(mf_cell_vector.first, q)[mf_cell_vector.second];
}
_state(liquid_state, mp_dof_index) = liquid_ratio_sum / n_q_points;
_state(solid_state, mp_dof_index) =
1. - _state(liquid_state, mp_dof_index);
}
}
else if constexpr (std::is_same_v<MaterialStates, SolidLiquidPowder>)
{
auto constexpr powder_state =
static_cast<unsigned int>(MaterialStates::State::powder);
for (auto const &cell : dealii::filter_iterators(
dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell()))
{
typename dealii::Triangulation<dim>::active_cell_iterator cell_tria(
cell);
auto mp_dof_index = get_dof_index(cell_tria);
auto const &mf_cell_vector = cell_it_to_mf_cell_map[cell];
unsigned int const n_q_points =
dof_handler.get_fe().tensor_degree() + 1;
double liquid_ratio_sum = 0.;
double powder_ratio_sum = 0.;
for (unsigned int q = 0; q < n_q_points; ++q)
{
liquid_ratio_sum +=
liquid_ratio(mf_cell_vector.first, q)[mf_cell_vector.second];
powder_ratio_sum +=
powder_ratio(mf_cell_vector.first, q)[mf_cell_vector.second];
}
_state(liquid_state, mp_dof_index) = liquid_ratio_sum / n_q_points;
_state(powder_state, mp_dof_index) = powder_ratio_sum / n_q_points;
_state(solid_state, mp_dof_index) = 1. -
_state(liquid_state, mp_dof_index) -
_state(powder_state, mp_dof_index);
}
}
}
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
set_state_device(
Kokkos::View<double *, typename MemorySpaceType::kokkos_space>
liquid_ratio,
Kokkos::View<double *, typename MemorySpaceType::kokkos_space>
powder_ratio,
std::map<typename dealii::DoFHandler<dim>::cell_iterator,
std::vector<unsigned int>> const &_cell_it_to_mf_pos,
dealii::DoFHandler<dim> const &dof_handler)
{
// Create a mapping between the matrix free dofs and material property dofs
unsigned int const n_q_points = dof_handler.get_fe().tensor_degree() + 1;
Kokkos::View<unsigned int **, dealii::MemorySpace::Default::kokkos_space>
mapping(Kokkos::view_alloc("mapping", Kokkos::WithoutInitializing),
_state.extent(1), n_q_points);
auto mapping_host =
Kokkos::create_mirror_view(Kokkos::WithoutInitializing, mapping);
Kokkos::View<dealii::types::global_dof_index *,
dealii::MemorySpace::Host::kokkos_space>
mp_dof_host("mp_dof_host", _state.extent(1));
// We only loop over the part of the domain which has material, i.e., not over
// FE_Nothing cell. This is because _cell_it_to_mf_pos does not exist for
// FE_Nothing cells. However, we have set the state of the material on the
// entire domain. This is not a problem since that state is unchanged and does
// not need to be updated.
unsigned int cell_i = 0;
for (auto const &cell : dealii::filter_iterators(
dof_handler.active_cell_iterators(),
dealii::IteratorFilters::ActiveFEIndexEqualTo(0, true)))
{
typename dealii::Triangulation<dim>::active_cell_iterator cell_tria(cell);
auto mp_dof_index = get_dof_index(cell_tria);
auto const &mf_cell_vector = _cell_it_to_mf_pos.at(cell);
for (unsigned int q = 0; q < n_q_points; ++q)
{
mapping_host(cell_i, q) = mf_cell_vector[q];
}
mp_dof_host(cell_i) = mp_dof_index;
++cell_i;
}
Kokkos::deep_copy(mapping, mapping_host);
Kokkos::View<dealii::types::global_dof_index *,
dealii::MemorySpace::Default::kokkos_space>
mp_dof(Kokkos::view_alloc("mp_dof", Kokkos::WithoutInitializing),
mp_dof_host.extent(0));
Kokkos::deep_copy(mp_dof, mp_dof_host);
if constexpr (std::is_same_v<MaterialStates, Solid>)
{
// When there is only Solid, we can just set _state to one.
Kokkos::deep_copy(_state, 1.);
}
else
{
auto const solid_state =
static_cast<unsigned int>(MaterialStates::State::solid);
auto const liquid_state =
static_cast<unsigned int>(MaterialStates::State::liquid);
if constexpr (std::is_same_v<MaterialStates, SolidLiquid>)
{
using ExecutionSpace = std::conditional_t<
std::is_same_v<MemorySpaceType, dealii::MemorySpace::Host>,
Kokkos::DefaultHostExecutionSpace, Kokkos::DefaultExecutionSpace>;
Kokkos::parallel_for(
"adamantine::set_state_device",
Kokkos::RangePolicy<ExecutionSpace>(0, cell_i),
KOKKOS_CLASS_LAMBDA(int i) {
double liquid_ratio_sum = 0.;
for (unsigned int q = 0; q < n_q_points; ++q)
{
liquid_ratio_sum += liquid_ratio(mapping(i, q));
}
_state(liquid_state, mp_dof(i)) = liquid_ratio_sum / n_q_points;
_state(solid_state, mp_dof(i)) =
1. - _state(liquid_state, mp_dof(i));
});
}
else if constexpr (std::is_same_v<MaterialStates, SolidLiquidPowder>)
{
auto const powder_state =
static_cast<unsigned int>(MaterialStates::State::powder);
using ExecutionSpace = std::conditional_t<
std::is_same_v<MemorySpaceType, dealii::MemorySpace::Host>,
Kokkos::DefaultHostExecutionSpace, Kokkos::DefaultExecutionSpace>;
Kokkos::parallel_for(
"adamantine::set_state_device",
Kokkos::RangePolicy<ExecutionSpace>(0, cell_i),
KOKKOS_CLASS_LAMBDA(int i) {
double liquid_ratio_sum = 0.;
double powder_ratio_sum = 0.;
for (unsigned int q = 0; q < n_q_points; ++q)
{
liquid_ratio_sum += liquid_ratio(mapping(i, q));
powder_ratio_sum += powder_ratio(mapping(i, q));
}
_state(liquid_state, mp_dof(i)) = liquid_ratio_sum / n_q_points;
_state(powder_state, mp_dof(i)) = powder_ratio_sum / n_q_points;
_state(solid_state, mp_dof(i)) = 1. -
_state(liquid_state, mp_dof(i)) -
_state(powder_state, mp_dof(i));
});
}
}
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
set_cell_state(
std::vector<std::array<double, MaterialStates::n_material_states>> const
&cell_state)
{
auto state_host =
Kokkos::create_mirror_view(Kokkos::WithoutInitializing, _state);
for (unsigned int i = 0; i < cell_state.size(); ++i)
{
for (unsigned int j = 0; j < MaterialStates::n_material_states; ++j)
{
state_host(j, i) = cell_state[i][j];
}
}
Kokkos::deep_copy(_state, state_host);
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates,
MemorySpaceType>::set_initial_state()
{
// Set the material state to the one defined by the user_index
std::vector<dealii::types::global_dof_index> mp_dofs_vec;
std::vector<unsigned int> user_indices_vec;
for (auto cell :
dealii::filter_iterators(_mp_dof_handler.active_cell_iterators(),
dealii::IteratorFilters::LocallyOwnedCell()))
{
std::vector<dealii::types::global_dof_index> mp_dof(1);
cell->get_dof_indices(mp_dof);
mp_dofs_vec.push_back(_dofs_map.at(mp_dof[0]));
user_indices_vec.push_back(cell->user_index());
}
typename MemorySpaceType::kokkos_space memory_space;
Kokkos::View<dealii::types::global_dof_index *, Kokkos::HostSpace>
mp_dofs_host(mp_dofs_vec.data(), mp_dofs_vec.size());
auto mp_dofs =
Kokkos::create_mirror_view_and_copy(memory_space, mp_dofs_host);
Kokkos::View<unsigned int *, Kokkos::HostSpace> user_indices_host(
user_indices_vec.data(), user_indices_vec.size());
auto user_indices =
Kokkos::create_mirror_view_and_copy(memory_space, user_indices_host);
Kokkos::deep_copy(_state, 0.);
using ExecutionSpace = std::conditional_t<
std::is_same_v<MemorySpaceType, dealii::MemorySpace::Host>,
Kokkos::DefaultHostExecutionSpace, Kokkos::DefaultExecutionSpace>;
Kokkos::parallel_for(
"adamantine::set_initial_state",
Kokkos::RangePolicy<ExecutionSpace>(0, user_indices.extent(0)),
KOKKOS_CLASS_LAMBDA(int i) { _state(user_indices(i), mp_dofs(i)) = 1.; });
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
void MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>::
fill_properties(boost::property_tree::ptree const &database)
{
// PropertyTreeInput materials.property_format
std::string property_format = database.get<std::string>("property_format");
_use_table = (property_format == "table");
// PropertyTreeInput materials.n_materials
unsigned int const n_materials = database.get<unsigned int>("n_materials");
// Find all the material_ids being used.
std::vector<dealii::types::material_id> material_ids;
for (dealii::types::material_id id = 0;
id < dealii::numbers::invalid_material_id; ++id)
{
if (database.count("material_" + std::to_string(id)) != 0)
material_ids.push_back(id);
if (material_ids.size() == n_materials)
break;
}
// When using the polynomial format we allocate one contiguous block of
// memory. Thus, the largest material_id should be as small as possible
unsigned int const n_material_ids =
*std::max_element(material_ids.begin(), material_ids.end()) + 1;
_properties = Kokkos::View<double *[g_n_properties],
typename MemorySpaceType::kokkos_space>(
Kokkos::view_alloc("properties", Kokkos::WithoutInitializing),
n_material_ids);
auto properties_host =
Kokkos::create_mirror_view(Kokkos::WithoutInitializing, _properties);
if (_use_table)
{
// View is initialized to zero in purpose
_state_property_tables =
Kokkos::View<double * [MaterialStates::n_material_states]
[g_n_thermal_state_properties][table_size][2],
typename MemorySpaceType::kokkos_space>(
"state_property_tables", n_material_ids);
// Mechanical properties only exist for the solid state. View is initialized
// to zero in purpose.
_mechanical_properties_tables_host =
Kokkos::View<double *[g_n_mechanical_state_properties][table_size][2],
typename dealii::MemorySpace::Host::kokkos_space>(
"mechanical_properties_tables_host", n_material_ids);
}
else
{
// View is initialized to zero in purpose
_state_property_polynomials =
Kokkos::View<double * [MaterialStates::n_material_states]
[g_n_thermal_state_properties][p_order + 1],
typename MemorySpaceType::kokkos_space>(
"state_property_polynomials", n_material_ids);
// Mechanical properties only exist for the solid state. View is initialized
// to zero in purpose
_mechanical_properties_polynomials_host =
Kokkos::View<double *[g_n_mechanical_state_properties][p_order + 1],
typename dealii::MemorySpace::Host::kokkos_space>(
"mechanical_properties_polynomials_host", n_material_ids);
}
auto state_property_tables_host = Kokkos::create_mirror_view_and_copy(
Kokkos::DefaultHostExecutionSpace{}, _state_property_tables);
auto state_property_polynomials_host = Kokkos::create_mirror_view_and_copy(
Kokkos::DefaultHostExecutionSpace{}, _state_property_polynomials);
for (auto const material_id : material_ids)
{
// Get the material property tree.
boost::property_tree::ptree const &material_database =
database.get_child("material_" + std::to_string(material_id));
// For each material, loop over the possible states.
for (unsigned int state = 0; state < MaterialStates::n_material_states;
++state)
{
// The state may or may not exist for the material.
boost::optional<boost::property_tree::ptree const &> state_database =
material_database.get_child_optional(material_state_names[state]);
if (state_database)
{
// For each state, loop over the possible properties.
for (unsigned int p = 0; p < g_n_state_properties; ++p)
{
// The property may or may not exist for that state
boost::optional<std::string> const property =
state_database.get().get_optional<std::string>(
state_property_names[p]);
// If the property exists, put it in the map. If the property does not
// exist, we have a nullptr.
if (property)
{
// Remove blank spaces
std::string property_string = property.get();
property_string.erase(
std::remove_if(property_string.begin(), property_string.end(),
[](unsigned char x) { return std::isspace(x); }),
property_string.end());
if (_use_table)
{
std::vector<std::string> parsed_property;
boost::split(parsed_property, property_string,
[](char c) { return c == ';'; });
unsigned int const parsed_property_size = parsed_property.size();
ASSERT_THROW(parsed_property_size <= table_size,
"Too many coefficients, increase the table size");
for (unsigned int i = 0; i < parsed_property_size; ++i)
{
std::vector<std::string> t_v;
boost::split(t_v, parsed_property[i],
[](char c) { return c == ','; });
ASSERT(t_v.size() == 2, "Error reading material property.");
if (p < g_n_thermal_state_properties)
{
state_property_tables_host(material_id, state, p, i, 0) =
std::stod(t_v[0]);
state_property_tables_host(material_id, state, p, i, 1) =
std::stod(t_v[1]);
}
else
{
if (state ==
static_cast<unsigned int>(MaterialStates::State::solid))
{
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties, i, 0) =
std::stod(t_v[0]);
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties, i, 1) =
std::stod(t_v[1]);
}
}
}
// fill the rest with the last value
for (unsigned int i = parsed_property_size; i < table_size; ++i)
{
if (p < g_n_thermal_state_properties)
{
state_property_tables_host(material_id, state, p, i, 0) =
state_property_tables_host(material_id, state, p, i - 1,
0);
state_property_tables_host(material_id, state, p, i, 1) =
state_property_tables_host(material_id, state, p, i - 1,
1);
}
else
{
if (state ==
static_cast<unsigned int>(MaterialStates::State::solid))
{
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties, i, 0) =
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties,
i - 1, 0);
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties, i, 1) =
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties,
i - 1, 1);
}
}
}
}
else
{
std::vector<std::string> parsed_property;
boost::split(parsed_property, property_string,
[](char c) { return c == ','; });
unsigned int const parsed_property_size = parsed_property.size();
ASSERT_THROW(
parsed_property_size <= p_order + 1,
"Too many coefficients, increase the polynomial order");
for (unsigned int i = 0; i < parsed_property_size; ++i)
{
if (p < g_n_thermal_state_properties)
{
state_property_polynomials_host(material_id, state, p, i) =
std::stod(parsed_property[i]);
}
else if (state == static_cast<unsigned int>(
MaterialStates::State::solid))
{
_mechanical_properties_polynomials_host(
material_id, p - g_n_thermal_state_properties, i) =
std::stod(parsed_property[i]);
}
}
}
}
else if (state_property_names[p] == "elastic_limit" &&
state ==
static_cast<unsigned int>(MaterialStates::State::solid))
{
// If the elastic limit is not provided, we solve a purely elastic
// problem. We set the elastic limit to infinity.
double infinity = std::numeric_limits<double>::infinity();
if (_use_table)
{
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties, 0, 0) =
infinity;
_mechanical_properties_tables_host(
material_id, p - g_n_thermal_state_properties, 0, 1) =
infinity;
}
else
{
_mechanical_properties_polynomials_host(
material_id, p - g_n_thermal_state_properties, 0) = infinity;
}
}
}
}
}
// Check for the properties that are associated to a material but that
// are independent of an individual state. These properties are duplicated
// for every state.
for (unsigned int p = 0; p < g_n_properties; ++p)