PRISMS-PF Manual
Loading...
Searching...
No Matches
user_constants.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/tensor.h>
7#include <deal.II/base/utilities.h>
8
9#include <boost/algorithm/string/predicate.hpp>
10#include <boost/range/algorithm_ext/erase.hpp>
11#include <boost/variant.hpp>
12
16
17#include <prismspf/config.h>
18
19#include <map>
20
22
26template <unsigned int dim>
28{
29public:
30 using InputVariant = boost::variant<double,
31 int,
32 bool,
33 dealii::Tensor<1, dim>,
34 dealii::Tensor<2, dim>,
35 dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>>;
36
41 construct_user_constant(std::vector<std::string> &model_constants_strings);
42
49 [[nodiscard]] double
50 get_model_constant_double(const std::string &constant_name) const;
51
58 [[nodiscard]] int
59 get_model_constant_int(const std::string &constant_name) const;
60
67 [[nodiscard]] bool
68 get_model_constant_bool(const std::string &constant_name) const;
69
76 [[nodiscard]] dealii::Tensor<1, dim>
77 get_model_constant_rank_1_tensor(const std::string &constant_name) const;
78
85 [[nodiscard]] dealii::Tensor<2, dim>
86 get_model_constant_rank_2_tensor(const std::string &constant_name) const;
87
94 [[nodiscard]] dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>
96
100 void
106
110 void
111 print() const;
112
113private:
117 unsigned int
118 compute_tensor_parentheses(const unsigned int &n_elements,
119 const std::vector<std::string> &tensor_elements);
120
124 void
125 remove_parentheses(std::vector<std::string> &tensor_elements);
126
130 dealii::Tensor<1, dim>
131 compute_rank_1_tensor_constant(const unsigned int &n_elements,
132 const std::vector<std::string> &tensor_elements);
133
137 dealii::Tensor<2, dim>
138 compute_rank_2_tensor_constant(const unsigned int &n_elements,
139 const std::vector<std::string> &tensor_elements);
140
145 primitive_model_constant(std::vector<std::string> &model_constants_strings);
146
147 [[nodiscard]] dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>
148 get_cij_tensor(std::vector<double> elastic_constants,
149 const std::string &elastic_const_symmetry) const;
150
151 [[nodiscard]] dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>
153 const std::vector<double> &constants) const;
154
158 std::map<std::string, InputVariant> model_constants;
159
163 class VariantPrinter : public boost::static_visitor<>
164 {
165 public:
166 void
167 operator()(double value) const
168 {
170 }
171
172 void
173 operator()(int value) const
174 {
176 }
177
178 void
179 operator()(bool value) const
180 {
181 ConditionalOStreams::pout_summary() << std::boolalpha << value;
182 }
183
184 void
185 operator()(const dealii::Tensor<1, dim> &value) const
186 {
187 ConditionalOStreams::pout_summary() << "Tensor<1, " << dim << ">: ";
188 for (unsigned int i = 0; i < dim; ++i)
189 {
190 ConditionalOStreams::pout_summary() << value[i] << ' ';
191 }
192 }
193
194 void
195 operator()(const dealii::Tensor<2, dim> &value) const
196 {
197 ConditionalOStreams::pout_summary() << "Tensor<2, " << dim << ">: ";
198 for (unsigned int i = 0; i < dim; ++i)
199 {
200 for (unsigned int j = 0; j < dim; ++j)
201 {
202 ConditionalOStreams::pout_summary() << value[i][j] << ' ';
203 }
204 }
205 }
206
207 template <unsigned int D = dim>
208 void
209 operator()(const dealii::Tensor<2, (2 * D) - 1 + (D / 3)> &value) const
210 requires((D != ((2 * D) - 1 + (D / 3))))
211 {
212 constexpr unsigned int dimension = (2 * D) - 1 + (D / 3);
213 ConditionalOStreams::pout_summary() << "Tensor<2, " << dimension << ">: ";
214 for (unsigned int i = 0; i < dimension; ++i)
215 {
216 for (unsigned int j = 0; j < dimension; ++j)
217 {
218 ConditionalOStreams::pout_summary() << value[i][j] << ' ';
219 }
220 }
221 }
222 };
223};
224
225template <unsigned int dim>
226inline double
228{
229 Assert(model_constants.find(constant_name) != model_constants.end(),
230 dealii::ExcMessage(
231 "Mismatch between constants in parameters.prm and CustomPDE.h. The constant "
232 "that you attempted to access was " +
233 constant_name + "."));
234
235 return boost::get<double>(model_constants.at(constant_name));
236}
237
238template <unsigned int dim>
239inline int
241{
242 Assert(model_constants.find(constant_name) != model_constants.end(),
243 dealii::ExcMessage(
244 "Mismatch between constants in parameters.prm and CustomPDE.h. The constant "
245 "that you attempted to access was " +
246 constant_name + "."));
247
248 return boost::get<int>(model_constants.at(constant_name));
249}
250
251template <unsigned int dim>
252inline bool
254{
255 Assert(model_constants.find(constant_name) != model_constants.end(),
256 dealii::ExcMessage(
257 "Mismatch between constants in parameters.prm and CustomPDE.h. The constant "
258 "that you attempted to access was " +
259 constant_name + "."));
260
261 return boost::get<bool>(model_constants.at(constant_name));
262}
263
264template <unsigned int dim>
265inline dealii::Tensor<1, dim>
267 const std::string &constant_name) const
268{
269 Assert(model_constants.find(constant_name) != model_constants.end(),
270 dealii::ExcMessage(
271 " Mismatch between constants in parameters.prm and "
272 "CustomPDE.h. The constant that you attempted to access was " +
273 constant_name + "."));
274
275 return boost::get<dealii::Tensor<1, dim>>(model_constants.at(constant_name));
276}
277
278template <unsigned int dim>
279inline dealii::Tensor<2, dim>
281 const std::string &constant_name) const
282{
283 Assert(model_constants.find(constant_name) != model_constants.end(),
284 dealii::ExcMessage(
285 "Mismatch between constants in parameters.prm and CustomPDE.h. The constant "
286 "that you attempted to access was " +
287 constant_name + "."));
288
289 return boost::get<dealii::Tensor<2, dim>>(model_constants.at(constant_name));
290}
291
292template <unsigned int dim>
293inline dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>
295 const std::string &constant_name) const
296{
297 Assert(model_constants.find(constant_name) != model_constants.end(),
298 dealii::ExcMessage(
299 "Mismatch between constants in parameters.prm and CustomPDE.h. The constant "
300 "that you attempted to access was " +
301 constant_name + "."));
302
303 return boost::get<dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>>(
304 model_constants.at(constant_name));
305}
306
307template <unsigned int dim>
308inline unsigned int
310 const unsigned int &n_elements,
311 const std::vector<std::string> &tensor_elements)
312{
313 unsigned int open_parentheses = 0;
314 unsigned int close_parentheses = 0;
315
316 for (unsigned int element = 0; element < n_elements; element++)
317 {
318 for (const char character : tensor_elements.at(element))
319 {
320 if (character == '(')
321 {
323 }
324 else if (character == ')')
325 {
327 }
328 }
329 }
330
332 {
333 AssertThrow(false,
334 dealii::ExcMessage("User-defined elastic constant list does not have "
335 "the same number of open and close parentheses."));
336 }
337
338 return open_parentheses;
339}
340
341template <unsigned int dim>
342inline void
344{
345 for (std::string &element : tensor_elements)
346 {
347 boost::range::remove_erase(element, '(');
348 boost::range::remove_erase(element, ')');
349 }
350}
351
352template <unsigned int dim>
353inline dealii::Tensor<1, dim>
355 const unsigned int &n_elements,
356 const std::vector<std::string> &tensor_elements)
357{
359 dealii::ExcMessage("The columns in user-defined constant tensors must be "
360 "equal to the maximum number of dimensions."));
361
362 dealii::Tensor<1, dim> temp;
363 for (unsigned int i = 0; i < dim; i++)
364 {
365 temp[i] = dealii::Utilities::string_to_double(tensor_elements.at(i));
366 }
367
368 return temp;
369}
370
371template <unsigned int dim>
372inline dealii::Tensor<2, dim>
374 const unsigned int &n_elements,
375 const std::vector<std::string> &tensor_elements)
376{
378 dealii::ExcMessage("User-defined constant tensor does not have the "
379 "correct number of elements, matrices must be 3x3."));
380
381 const unsigned int row_length = 3;
382
383 dealii::Tensor<2, dim> temp;
384 for (unsigned int i = 0; i < dim; i++)
385 {
386 for (unsigned int j = 0; j < dim; j++)
387 {
388 temp[i][j] =
389 dealii::Utilities::string_to_double(tensor_elements.at((i * row_length) + j));
390 }
391 }
392
393 return temp;
394}
395
396template <unsigned int dim>
399 std::vector<std::string> &model_constants_strings)
400{
401 // Ensure that the input includes a value and a type
403 model_constants_strings.size() > 1,
404 dealii::ExcMessage(
405 "At least two fields are required for user-defined variables (value and type)."));
406
407 std::vector<std::string> model_constants_type_strings =
408 dealii::Utilities::split_string_list(model_constants_strings.at(
409 model_constants_strings.size() - 1),
410 ' ');
411
412 if (model_constants_strings.size() == 2)
413 {
414 return primitive_model_constant(model_constants_strings);
415 }
416
417 if (boost::iequals(model_constants_type_strings.at(0), "tensor"))
418 {
419 const unsigned int n_elements = model_constants_strings.size() - 1;
420
421 const unsigned int open_parentheses =
422 compute_tensor_parentheses(n_elements, model_constants_strings);
423
425 FeatureNotImplemented("3rd rank tensors and above"));
426
427 remove_parentheses(model_constants_strings);
428
429 // Rank 1 tensor
430 if (open_parentheses == 1)
431 {
432 return compute_rank_1_tensor_constant(n_elements, model_constants_strings);
433 }
434 // Rank 2 tensor
435 return compute_rank_2_tensor_constant(n_elements, model_constants_strings);
436 }
437 if (boost::iequals(model_constants_type_strings.at(1), "elastic") &&
438 boost::iequals(model_constants_type_strings.at(2), "constants"))
439 {
440 const unsigned int n_elements = model_constants_strings.size() - 1;
441
442 remove_parentheses(model_constants_strings);
443
444 // Load in the elastic constants as a vector
445 std::vector<double> temp_elastic_constants(n_elements);
446 for (unsigned int i = 0; i < n_elements; i++)
447 {
449 dealii::Utilities::string_to_double(model_constants_strings.at(i));
450 }
451 const std::string &elastic_const_symmetry = model_constants_type_strings.at(0);
452 dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)> temp =
454 return temp;
455 }
456
457 AssertThrow(false,
458 dealii::ExcMessage(
459 "Only user-defined constant tensors may have multiple elements."));
460 return 0;
461}
462
463template <unsigned int dim>
466 std::vector<std::string> &model_constants_strings)
467{
468 std::vector<std::string> model_constants_type_strings =
469 dealii::Utilities::split_string_list(model_constants_strings.at(
470 model_constants_strings.size() - 1),
471 ' ');
472
473 if (boost::iequals(model_constants_type_strings.at(0), "double"))
474 {
475 return dealii::Utilities::string_to_double(model_constants_strings.at(0));
476 }
477 if (boost::iequals(model_constants_type_strings.at(0), "int"))
478 {
479 return dealii::Utilities::string_to_int(model_constants_strings.at(0));
480 }
481 if (boost::iequals(model_constants_type_strings.at(0), "bool"))
482 {
483 return boost::iequals(model_constants_strings.at(0), "true");
484 }
485
486 AssertThrow(false,
487 dealii::ExcMessage(
488 "The type for user-defined variables must be `double`, `int`, "
489 "`bool`, `tensor`, or `elastic constants`."));
490 return 0;
491}
492
493template <unsigned int dim>
494inline dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>
496 const std::string &elastic_const_symmetry) const
497{
498 // First set the material model
500 if (elastic_const_symmetry == "isotropic")
501 {
503 }
504 else if (elastic_const_symmetry == "transverse")
505 {
507 }
508 else if (elastic_const_symmetry == "orthotropic")
509 {
511 }
512 else if (elastic_const_symmetry == "anisotropic")
513 {
515 }
516 else
517 {
518 AssertThrow(false, dealii::ExcMessage("Invalid elasticity tensor type"));
519 }
520
521 // If the material model is anisotropic for a 2D calculation but the elastic
522 // constants are given for a 3D calculation, change the elastic constant
523 // vector to the 2D form
524 constexpr unsigned int max_number = 21;
525 if ((mat_model == Anisotropic) && (dim == 2) && elastic_constants.size() == max_number)
526 {
527 std::vector<double> elastic_constants_temp = elastic_constants;
528 elastic_constants.clear();
529 const std::vector<unsigned int> indices_2d = {0, 1, 5, 6, 10, 14};
530 std::transform(indices_2d.begin(),
531 indices_2d.end(),
532 std::back_inserter(elastic_constants),
533 [&elastic_constants_temp](unsigned int index)
534 {
535 return elastic_constants_temp.at(index);
536 });
537 }
538
539 return get_cij_matrix(mat_model, elastic_constants);
540}
541
542template <unsigned int dim>
543inline dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)>
545 const std::vector<double> &constants) const
546{
547 // Initialize compliance tensor
548 dealii::Tensor<2, (2 * dim) - 1 + (dim / 3)> compliance;
549
550 switch (dim)
551 {
552 case 1:
553 {
554 const int xx_dir = 0;
555
556 switch (model)
557 {
558 // For the 1D case, it make little sense to accept anything besides an
559 // isotropic elasticity tensor. One hiccup, is that if a user is debugging
560 // an application and switches to 1D, they will have to modify the
561 // elasticity constants to align with that. While more burdensome, there's
562 // less of a chance of producing spurious behavior.
563 case Isotropic:
564 {
565 const double modulus = constants.at(0);
566
568 break;
569 }
570 default:
571 AssertThrow(false,
572 dealii::ExcMessage(
573 "Invalid elasticity model type for 1D. We only accept "
574 "isotropic elasticity tensors."));
575 }
576 break;
577 }
578 case 2:
579 {
580 // The voigt indexing scheme for 2 dimensions
581 const int xx_dir = 0;
582 const int yy_dir = 1;
583 const int xy_dir = 2;
584
585 switch (model)
586 {
587 // Like the 1D case, it is nonsensical to have transverse or orthotropic
588 // compliance tensors, so we throw an error.
589 case Isotropic:
590 {
591 // For isotropic compliance tensors, we can simplify the computation to
592 // two parameters: $\lambda$ and $\mu$, where $\mu$ is the shear
593 // modulus. In cartesian coordinates,
594 // $$$
595 // c_{ijkl} = \lambda \delta_{ij} \delta_{kl} + \mu (\delta_{ik}
596 // \delta_{jl} + \delta_{il} \delta_{kj})
597 // $$$
598 const double modulus = constants.at(0);
599 const double poisson = constants.at(1);
600
601 const double shear_modulus = modulus / (2 * (1 + poisson));
602 const double lambda =
603 poisson * modulus / ((1 + poisson) * (1 - 2 * poisson));
604
606 lambda + 2 * shear_modulus;
609 break;
610 }
611 case Anisotropic:
612 {
613 // In the anisotropic case, every entry is specified (given the symmetry
614 // constraints). Also, ignore magic numbers because it is simpler to
615 // hardcode this.
616
617 // NOLINTBEGIN(cppcoreguidelines-avoid-magic-numbers,readability-magic-numbers)
618
623 constants.at(3);
625 constants.at(4);
627 constants.at(5);
628
629 // NOLINTEND(cppcoreguidelines-avoid-magic-numbers,readability-magic-numbers)
630 break;
631 }
632 default:
633 AssertThrow(false, dealii::ExcMessage("Invalid elasticity model type"));
634 }
635 break;
636 }
637 case 3:
638 {
639 const int xx_dir = 0;
640 const int yy_dir = 1;
641 const int zz_dir = 2;
642 const int yz_dir = 3;
643 const int xz_dir = 4;
644 const int xy_dir = 5;
645
646 switch (model)
647 {
648 case Isotropic:
649 {
650 // For isotropic compliance tensors, we can simplify the computation to
651 // two parameters: $\lambda$ and $\mu$, where $\mu$ is the shear
652 // modulus. In cartesian coordinates,
653 // $$$
654 // c_{ijkl} = \lambda \delta_{ij} \delta_{kl} + \mu (\delta_{ik}
655 // \delta_{jl} + \delta_{il} \delta_{kj})
656 // $$$
657 const double modulus = constants.at(0);
658 const double poisson = constants.at(1);
659
660 const double shear_modulus = modulus / (2 * (1 + poisson));
661 const double lambda =
662 poisson * modulus / ((1 + poisson) * (1 - 2 * poisson));
663
671 break;
672 }
673 case Anisotropic:
674 {
675 // In the anisotropic case, every entry is specified (given the symmetry
676 // constraints). Also, ignore magic numbers because it is simpler to
677 // hardcode this.
678
679 // NOLINTBEGIN(cppcoreguidelines-avoid-magic-numbers,readability-magic-numbers)
680
702
703 // NOLINTEND(cppcoreguidelines-avoid-magic-numbers,readability-magic-numbers)
704 break;
705 }
706 case Transverse:
707 // TODO (landinjm): implement
708 case Orthotropic:
709 // TODO (landinjm): implement
710 default:
711 AssertThrow(false, dealii::ExcMessage("Invalid elasticity model type"));
712 }
713 break;
714 }
715 default:
716 {
717 Assert(false, UnreachableCode());
718 }
719 }
720
721 return compliance;
722}
723
724template <unsigned int dim>
725void
727{
728 if (!model_constants.empty())
729 {
731 << "================================================\n"
732 << " User Constants\n"
733 << "================================================\n";
734
735 for (const auto &[constant_name, variant] : model_constants)
736 {
738 boost::apply_visitor(VariantPrinter(), variant);
740 }
741 ConditionalOStreams::pout_summary() << "\n" << std::flush;
742 }
743}
744
static dealii::ConditionalOStream & pout_summary()
Log output stream for writing a summary.log file.
Definition conditional_ostreams.cc:34
Class for printing of variant types. This is bad practice and should be fixed.
Definition user_constants.h:164
void operator()(bool value) const
Definition user_constants.h:179
void operator()(int value) const
Definition user_constants.h:173
void operator()(double value) const
Definition user_constants.h:167
void operator()(const dealii::Tensor< 2, dim > &value) const
Definition user_constants.h:195
void operator()(const dealii::Tensor< 2,(2 *D) - 1+(D/3)> &value) const
Definition user_constants.h:209
void operator()(const dealii::Tensor< 1, dim > &value) const
Definition user_constants.h:185
Class the stores and manages user-defined constants.
Definition user_constants.h:28
InputVariant construct_user_constant(std::vector< std::string > &model_constants_strings)
Assign the specified user constant to whatever type.
Definition user_constants.h:398
double get_model_constant_double(const std::string &constant_name) const
Retrieve the double from the model_constants that are defined from the parameters....
Definition user_constants.h:227
dealii::Tensor< 2,(2 *dim) - 1+(dim/3)> get_model_constant_elasticity_tensor(const std::string &constant_name) const
Retrieve the elasticity tensor from the model_constants that are defined from the parameters....
Definition user_constants.h:294
dealii::Tensor< 2,(2 *dim) - 1+(dim/3)> get_cij_matrix(const ElasticityModel &model, const std::vector< double > &constants) const
Definition user_constants.h:544
dealii::Tensor< 2, dim > get_model_constant_rank_2_tensor(const std::string &constant_name) const
Retrieve the rank 2 tensor from the model_constants that are defined from the parameters....
Definition user_constants.h:280
bool get_model_constant_bool(const std::string &constant_name) const
Retrieve the bool from the model_constants that are defined from the parameters.prm parser....
Definition user_constants.h:253
void add_user_constant(const std::string &constant_name, std::vector< std::string > &model_constants_strings)
Print all user-specified constants.
Definition user_constants.h:101
std::map< std::string, InputVariant > model_constants
List of user-defined constants.
Definition user_constants.h:158
InputVariant primitive_model_constant(std::vector< std::string > &model_constants_strings)
Assign the primitive user constants (e.g., int, double, bool).
Definition user_constants.h:465
boost::variant< double, int, bool, dealii::Tensor< 1, dim >, dealii::Tensor< 2, dim >, dealii::Tensor< 2,(2 *dim) - 1+(dim/3)> > InputVariant
Definition user_constants.h:35
dealii::Tensor< 2,(2 *dim) - 1+(dim/3)> get_cij_tensor(std::vector< double > elastic_constants, const std::string &elastic_const_symmetry) const
Definition user_constants.h:495
dealii::Tensor< 1, dim > compute_rank_1_tensor_constant(const unsigned int &n_elements, const std::vector< std::string > &tensor_elements)
Compute a 1st rank tensor from user inputs .
Definition user_constants.h:354
int get_model_constant_int(const std::string &constant_name) const
Retrieve the int from the model_constants that are defined from the parameters.prm parser....
Definition user_constants.h:240
dealii::Tensor< 2, dim > compute_rank_2_tensor_constant(const unsigned int &n_elements, const std::vector< std::string > &tensor_elements)
Compute a 2nd rank tensor from user inputs .
Definition user_constants.h:373
unsigned int compute_tensor_parentheses(const unsigned int &n_elements, const std::vector< std::string > &tensor_elements)
Compute the number of tensor rows.
Definition user_constants.h:309
void print() const
Print all user-specified constants.
Definition user_constants.h:726
dealii::Tensor< 1, dim > get_model_constant_rank_1_tensor(const std::string &constant_name) const
Retrieve the rank 1 tensor from the model_constants that are defined from the parameters....
Definition user_constants.h:266
void remove_parentheses(std::vector< std::string > &tensor_elements)
Remove and leading and trailing parentheses.
Definition user_constants.h:343
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
Definition conditional_ostreams.cc:20
ElasticityModel
Symmetry of elastic tensor.
Definition type_enums.h:40
@ Transverse
Definition type_enums.h:42
@ Anisotropic
Definition type_enums.h:44
@ Orthotropic
Definition type_enums.h:43
@ Isotropic
Definition type_enums.h:41