16 #ifndef dealii_differentiation_ad_ad_helpers_h
17 #define dealii_differentiation_ad_ad_helpers_h
21 #if defined(DEAL_II_WITH_ADOLC) || defined(DEAL_II_TRILINOS_WITH_SACADO)
171 typename ScalarType =
double>
246 print(std::ostream &stream)
const;
269 std::ostream & stream)
const;
298 const bool ensure_persistent_setting =
true);
347 const bool clear_registered_tapes =
true);
455 const bool overwrite_tape =
false,
456 const bool keep_independent_values =
true);
622 const bool read_mode);
839 typename ScalarType =
double>
929 template <
typename VectorType>
933 const std::vector<::types::global_dof_index> &local_dof_indices);
955 const std::vector<ad_type> &
1004 template <
typename VectorType>
1008 const std::vector<::types::global_dof_index> &local_dof_indices);
1223 typename ScalarType =
double>
1538 typename ScalarType =
double>
1659 template <
int dim,
typename ExtractorType>
1679 static const unsigned int rank = 0;
1684 template <
typename NumberType>
1689 "The number of components doesn't match that of the corresponding tensor type.");
1692 "The rank doesn't match that of the corresponding tensor type.");
1701 template <
typename NumberType>
1707 template <
typename NumberType>
1713 static inline unsigned int
1729 (void)unrolled_index;
1740 template <
typename IndexType =
unsigned int,
int rank_in>
1743 const unsigned int column_offset)
1746 (void)table_indices;
1747 (void)column_offset;
1769 static const unsigned int rank = 1;
1774 template <
typename NumberType>
1779 "The number of components doesn't match that of the corresponding tensor type.");
1782 "The rank doesn't match that of the corresponding tensor type.");
1787 template <
typename NumberType>
1793 template <
typename NumberType>
1799 static inline unsigned int
1816 (void)unrolled_index;
1824 template <
int rank_in>
1827 const unsigned int column_offset)
1848 template <
typename IndexType =
unsigned int,
int rank_in>
1851 const unsigned int column_offset)
1855 "Cannot extract more table indices than the input table has!");
1857 return TensorType::component_to_unrolled_index(
1858 table_index_view(table_indices, column_offset));
1880 static const unsigned int rank = 1;
1885 template <
typename NumberType>
1891 template <
typename NumberType>
1897 template <
typename NumberType>
1903 static inline unsigned int
1919 (void)unrolled_index;
1927 template <
int rank_in>
1930 const unsigned int column_offset)
1951 template <
typename IndexType =
unsigned int,
int rank_in>
1954 const unsigned int column_offset)
1958 "Cannot extract more table indices than the input table has!");
1960 return TensorType::component_to_unrolled_index(
1961 table_index_view(table_indices, column_offset));
1988 template <
typename NumberType>
1994 template <
typename NumberType>
2000 template <
typename NumberType>
2006 static inline unsigned int
2022 (void)unrolled_index;
2030 template <
int rank_in>
2033 const unsigned int column_offset)
2038 table_indices[column_offset + 1]);
2056 template <
typename IndexType =
unsigned int,
int rank_in>
2059 const unsigned int column_offset)
2063 "Cannot extract more table indices than the input table has!");
2065 return TensorType::component_to_unrolled_index(
2066 table_index_view(table_indices, column_offset));
2093 template <
typename NumberType>
2099 template <
typename NumberType>
2105 template <
typename NumberType>
2111 static inline unsigned int
2130 return table_indices[0] != table_indices[1];
2137 template <
int rank_in>
2140 const unsigned int column_offset)
2145 table_indices[column_offset + 1]);
2163 template <
typename IndexType =
unsigned int,
int rank_in>
2166 const unsigned int column_offset)
2170 "Cannot extract more table indices than the input table has!");
2172 return TensorType::component_to_unrolled_index(
2173 table_index_view(table_indices, column_offset));
2183 template <
int dim,
typename NumberType,
typename ExtractorType>
2192 ExtractorType>::template tensor_type<NumberType>;
2205 template <
typename ExtractorType_Row,
typename ExtractorType_Col>
2219 template <
int rank,
int dim,
typename NumberType>
2243 template <
int ,
int dim,
typename NumberType>
2268 template <
int ,
int dim,
typename NumberType>
2292 template <
int ,
int dim,
typename NumberType>
2308 typename ExtractorType_Row,
2309 typename ExtractorType_Col>
2334 template <
int dim,
typename NumberType,
typename ExtractorType_Field>
2350 typename ExtractorType_Field,
2351 typename ExtractorType_Derivative>
2354 ExtractorType_Field,
2355 ExtractorType_Derivative>;
2364 typename IndexType =
unsigned int,
2365 typename ExtractorType>
2366 std::vector<IndexType>
2368 const bool ignore_symmetries =
true)
2370 (void)ignore_symmetries;
2373 const IndexType comp_first =
2376 std::iota(indices.begin(), indices.end(), comp_first);
2389 template <
int dim,
typename IndexType =
unsigned int>
2390 std::vector<IndexType>
2393 const bool ignore_symmetries =
true)
2398 if (ignore_symmetries ==
true)
2400 const IndexType comp_first =
2402 extractor_symm_tensor);
2404 std::iota(indices.begin(), indices.end(), comp_first);
2412 std::vector<IndexType> indices =
2413 extract_field_component_indices<dim>(extractor_tensor,
true);
2417 for (
unsigned int i = 0; i < indices.size(); ++i)
2423 const IndexType sti_new_index =
2426 indices[i] = sti_new_index;
2439 template <
typename TensorType,
typename NumberType>
2442 const unsigned int unrolled_index,
2447 t[TensorType::unrolled_to_component_indices(unrolled_index)] =
value;
2455 template <
int dim,
typename NumberType>
2457 const unsigned int unrolled_index,
2461 (void)unrolled_index;
2471 template <
typename NumberType>
2474 const unsigned int unrolled_index,
2478 (void)unrolled_index;
2488 template <
int dim,
typename NumberType>
2490 const unsigned int unrolled_index_row,
2491 const unsigned int unrolled_index_col,
2498 SubTensorType::n_independent_components);
2500 SubTensorType::n_independent_components);
2502 SubTensorType::unrolled_to_component_indices(unrolled_index_row);
2504 SubTensorType::unrolled_to_component_indices(unrolled_index_col);
2505 t[indices_row[0]][indices_row[1]][indices_col[0]][indices_col[1]] =
2516 typename NumberType,
2517 template <
int,
int,
typename>
class TensorType>
2520 const unsigned int unrolled_index)
2524 return t[TensorType<rank, dim, NumberType>::
2525 unrolled_to_component_indices(unrolled_index)];
2534 typename NumberType,
2535 template <
int,
int,
typename>
class TensorType>
2538 const unsigned int unrolled_index)
2541 (void)unrolled_index;
2551 template <
typename NumberType>
2552 inline const NumberType &
2556 (void)unrolled_index;
2567 typename NumberType,
2568 template <
int,
int,
typename>
class TensorType>
2571 const unsigned int unrolled_index)
2575 return t[TensorType<rank, dim, NumberType>::
2576 unrolled_to_component_indices(unrolled_index)];
2585 typename NumberType,
2586 template <
int,
int,
typename>
class TensorType>
2588 const unsigned int index)
2601 template <
typename NumberType>
2639 typename ScalarType =
double>
2641 :
public HelperBase<ADNumberTypeCode, ScalarType>
2706 const bool clear_registered_tapes =
true)
override;
2754 template <
typename ValueType,
typename ExtractorType>
2757 const ExtractorType &extractor);
2779 const std::vector<ad_type> &
2807 template <
typename ExtractorType>
2809 ExtractorType>::template tensor_type<ad_type>
2865 template <
typename ValueType,
typename ExtractorType>
2868 const ExtractorType &extractor);
2890 const bool symmetric_component,
3086 typename ScalarType =
double>
3223 template <
typename ExtractorType_Row>
3227 const ExtractorType_Row & extractor_row);
3267 template <
typename ExtractorType_Row,
typename ExtractorType_Col>
3271 ExtractorType_Col>::type
3273 const ExtractorType_Row & extractor_row,
3274 const ExtractorType_Col & extractor_col);
3479 typename ScalarType =
double>
3565 template <
typename ValueType,
typename ExtractorType>
3568 const ExtractorType &extractor);
3611 template <
typename ExtractorType_Row>
3615 const ExtractorType_Row & extractor_row);
3664 template <
typename ExtractorType_Row,
typename ExtractorType_Col>
3668 ExtractorType_Col>::type
3670 const ExtractorType_Row & extractor_row,
3671 const ExtractorType_Col & extractor_col);
3735 template <enum AD::NumberTypes ADNumberTypeCode,
typename ScalarType>
3736 template <
typename VectorType>
3740 const std::vector<::types::global_dof_index> &local_dof_indices)
3747 local_dof_indices.size() == this->n_independent_variables(),
3749 "Degree of freedom index vector size does not match number of independent variables"));
3750 for (
unsigned int i = 0; i < this->n_independent_variables(); ++i)
3752 Assert(this->registered_independent_variable_values[i] ==
false,
3753 ExcMessage(
"Independent variables already registered."));
3755 set_dof_values(values, local_dof_indices);
3760 template <enum AD::NumberTypes ADNumberTypeCode,
typename ScalarType>
3761 template <
typename VectorType>
3765 const std::vector<::types::global_dof_index> &local_dof_indices)
3767 Assert(local_dof_indices.size() == this->n_independent_variables(),
3769 "Vector size does not match number of independent variables"));
3770 for (
unsigned int i = 0; i < this->n_independent_variables(); ++i)
3772 i, values[local_dof_indices[i]]);
3783 typename ScalarType>
3784 template <
typename ValueType,
typename ExtractorType>
3788 const ExtractorType &extractor)
3795 const std::vector<unsigned int> index_set(
3796 internal::extract_field_component_indices<dim>(extractor));
3797 for (
const unsigned int index : index_set)
3800 this->registered_independent_variable_values[index] ==
false,
3802 "Overlapping indices for independent variables. "
3803 "One or more indices associated with the field that "
3804 "is being registered as an independent variable have "
3805 "already been associated with another field. This suggests "
3806 "that the component offsets used to construct their counterpart "
3807 "extractors are incompatible with one another. Make sure that "
3808 "the first component for each extractor properly takes into "
3809 "account the dimensionality of the preceding fields."));
3812 set_independent_variable(
value, extractor);
3819 typename ScalarType>
3820 template <
typename ValueType,
typename ExtractorType>
3824 const ExtractorType &extractor)
3826 const std::vector<unsigned int> index_set(
3827 internal::extract_field_component_indices<dim>(extractor));
3828 for (
unsigned int i = 0; i < index_set.size(); ++i)
3830 set_sensitivity_value(
3832 internal::Extractor<dim, ExtractorType>::symmetric_component(i),
3841 typename ScalarType>
3842 template <
typename ExtractorType>
3843 typename internal::Extractor<dim, ExtractorType>::template tensor_type<
3844 typename PointLevelFunctionsBase<dim, ADNumberTypeCode, ScalarType>::
3849 if (ADNumberTraits<ad_type>::is_taped ==
true)
3851 Assert(this->active_tape_index() !=
3858 this->finalize_sensitive_independent_variables();
3859 Assert(this->independent_variables.size() ==
3860 this->n_independent_variables(),
3862 this->n_independent_variables()));
3864 const std::vector<unsigned int> index_set(
3865 internal::extract_field_component_indices<dim>(extractor));
3866 typename internal::Extractor<dim,
3867 ExtractorType>::template tensor_type<ad_type>
3870 for (
unsigned int i = 0; i < index_set.size(); ++i)
3872 const unsigned int index = index_set[i];
3874 Assert(this->registered_independent_variable_values[index] ==
true,
3877 this->independent_variables[index];
3891 typename ScalarType>
3892 template <
typename ExtractorType_Row>
3893 typename internal::ScalarFieldGradient<
3896 ExtractorType_Row>::type
3899 const ExtractorType_Row & extractor_row)
3908 const std::vector<unsigned int> row_index_set(
3909 internal::extract_field_component_indices<dim>(extractor_row));
3910 Assert(out.n_independent_components == row_index_set.size(),
3911 ExcMessage(
"Not all tensor components have been extracted!"));
3912 for (
unsigned int r = 0; r < row_index_set.size(); ++r)
3922 typename ScalarType>
3923 template <
typename ExtractorType_Row,
typename ExtractorType_Col>
3924 typename internal::ScalarFieldHessian<
3928 ExtractorType_Col>::type
3931 const ExtractorType_Row & extractor_row,
3932 const ExtractorType_Col & extractor_col)
3934 using InternalHessian = internal::ScalarFieldHessian<dim,
3938 using InternalExtractorRow = internal::Extractor<dim, ExtractorType_Row>;
3939 using InternalExtractorCol = internal::Extractor<dim, ExtractorType_Col>;
3940 using HessianType =
typename InternalHessian::type;
3956 const std::vector<unsigned int> row_index_set(
3957 internal::extract_field_component_indices<dim>(
3958 extractor_row,
false ));
3959 const std::vector<unsigned int> col_index_set(
3960 internal::extract_field_component_indices<dim>(
3961 extractor_col,
false ));
3963 for (
unsigned int index = 0;
3964 index < HessianType::n_independent_components;
3968 HessianType::unrolled_to_component_indices(index);
3969 const unsigned int r =
3970 InternalExtractorRow::local_component(ti_out, 0);
3971 const unsigned int c =
3972 InternalExtractorCol::local_component(ti_out,
3973 InternalExtractorRow::rank);
3976 out, index, hessian[row_index_set[r]][col_index_set[c]]);
3990 typename ScalarType>
3991 template <
typename ValueType,
typename ExtractorType>
3995 const ExtractorType &extractor)
3997 const std::vector<unsigned int> index_set(
3998 internal::extract_field_component_indices<dim>(extractor));
3999 for (
unsigned int i = 0; i < index_set.size(); ++i)
4001 Assert(this->registered_marked_dependent_variables[index_set[i]] ==
4003 ExcMessage(
"Overlapping indices for dependent variables."));
4013 typename ScalarType>
4014 template <
typename ExtractorType_Row>
4018 ExtractorType_Row>::type
4020 const Vector<scalar_type> &values,
4021 const ExtractorType_Row & extractor_row)
4025 typename internal::VectorFieldValue<dim, scalar_type, ExtractorType_Row>::
4030 const std::vector<unsigned int> row_index_set(
4031 internal::extract_field_component_indices<dim>(extractor_row));
4032 Assert(out.n_independent_components == row_index_set.size(),
4033 ExcMessage(
"Not all tensor components have been extracted!"));
4034 for (
unsigned int r = 0; r < row_index_set.size(); ++r)
4044 typename ScalarType>
4045 template <
typename ExtractorType_Row,
typename ExtractorType_Col>
4050 ExtractorType_Col>::type
4053 const ExtractorType_Row & extractor_row,
4054 const ExtractorType_Col & extractor_col)
4060 using InternalExtractorRow = internal::Extractor<dim, ExtractorType_Row>;
4061 using InternalExtractorCol = internal::Extractor<dim, ExtractorType_Col>;
4062 using JacobianType =
typename InternalJacobian::type;
4078 const std::vector<unsigned int> row_index_set(
4079 internal::extract_field_component_indices<dim>(
4080 extractor_row,
false ));
4081 const std::vector<unsigned int> col_index_set(
4082 internal::extract_field_component_indices<dim>(
4083 extractor_col,
false ));
4085 for (
unsigned int index = 0;
4086 index < JacobianType::n_independent_components;
4090 JacobianType::unrolled_to_component_indices(index);
4091 const unsigned int r =
4092 InternalExtractorRow::local_component(ti_out, 0);
4093 const unsigned int c =
4094 InternalExtractorCol::local_component(ti_out,
4095 InternalExtractorRow::rank);
4098 out, index, jacobian[row_index_set[r]][col_index_set[c]]);
4114 #endif // defined(DEAL_II_WITH_ADOLC) || defined(DEAL_II_TRILINOS_WITH_SACADO)
4116 #endif // dealii_differentiation_ad_ad_helpers_h