15#ifndef dealii_solver_cg_h
16#define dealii_solver_cg_h
42 template <
typename,
typename>
176template <
typename VectorType = Vector<
double>>
215 template <
typename MatrixType,
typename PreconditionerType>
219 void solve(const MatrixType &A,
222 const PreconditionerType &preconditioner);
230 boost::signals2::connection
231 connect_coefficients_slot(
232 const
std::function<
void(typename VectorType::value_type,
233 typename VectorType::value_type)> &slot);
241 boost::signals2::connection
242 connect_condition_number_slot(const
std::function<
void(
double)> &slot,
243 const
bool every_iteration = false);
251 boost::signals2::connection
252 connect_eigenvalues_slot(
253 const
std::function<
void(const
std::vector<
double> &)> &slot,
254 const
bool every_iteration = false);
263 print_vectors(const
unsigned int step,
266 const VectorType &d) const;
274 compute_eigs_and_cond(
275 const
std::vector<typename VectorType::value_type> &diagonal,
276 const
std::vector<typename VectorType::value_type> &offdiagonal,
277 const
boost::signals2::signal<
void(const
std::vector<
double> &)>
279 const
boost::signals2::signal<
void(
double)> &cond_signal);
289 boost::signals2::signal<
void(typename VectorType::value_type,
290 typename VectorType::value_type)>
297 boost::signals2::signal<
void(
double)> condition_number_signal;
303 boost::signals2::signal<
void(
double)> all_condition_numbers_signal;
309 boost::signals2::signal<
void(const
std::vector<
double> &)> eigenvalues_signal;
315 boost::signals2::signal<
void(const
std::vector<
double> &)>
316 all_eigenvalues_signal;
326 bool determine_beta_by_flexible_formula;
356template <typename VectorType =
Vector<
double>>
398template <
typename VectorType>
402 const AdditionalData &data)
404 , additional_data(data)
405 , determine_beta_by_flexible_formula(false)
410template <
typename VectorType>
414 , additional_data(data)
415 , determine_beta_by_flexible_formula(false)
420template <
typename VectorType>
425 const VectorType &)
const
430template <
typename VectorType>
433 const std::vector<typename VectorType::value_type> &diagonal,
434 const std::vector<typename VectorType::value_type> &offdiagonal,
435 const boost::signals2::signal<
void(
const std::vector<double> &)>
437 const boost::signals2::signal<
void(
double)> &cond_signal)
440 if (!cond_signal.empty() || !eigenvalues_signal.empty())
444 for (size_type i = 0; i <
diagonal.size(); ++i)
448 T(i, i + 1) = offdiagonal[i];
450 T.compute_eigenvalues();
454 auto condition_number = T.eigenvalue(T.n() - 1) / T.eigenvalue(0);
457 cond_signal(
std::abs(condition_number));
461 if (!eigenvalues_signal.empty())
464 for (
unsigned int j = 0; j < T.n(); ++j)
487 template <
typename VectorType,
489 typename PreconditionerType>
490 struct IterationWorkerBase
492 using Number =
typename VectorType::value_type;
495 const PreconditionerType &preconditioner;
515 Number r_dot_preconditioner_dot_r;
518 double residual_norm;
519 Number previous_alpha;
521 IterationWorkerBase(
const MatrixType &A,
522 const PreconditionerType &preconditioner,
527 , preconditioner(preconditioner)
538 , r_dot_preconditioner_dot_r(Number())
542 , previous_alpha(Number())
546 startup(
const VectorType &b)
566 residual_norm = r.l2_norm();
574 template <
typename VectorType,
576 typename PreconditionerType,
578 struct IterationWorker
579 :
public IterationWorkerBase<VectorType, MatrixType, PreconditionerType>
582 IterationWorkerBase<VectorType, MatrixType, PreconditionerType>;
584 IterationWorker(
const MatrixType &A,
585 const PreconditionerType &preconditioner,
589 : BaseClass(A, preconditioner, flexible, memory, x)
593 using BaseClass::alpha;
594 using BaseClass::beta;
596 using BaseClass::preconditioner;
598 using BaseClass::r_dot_preconditioner_dot_r;
599 using BaseClass::residual_norm;
605 do_iteration(
const unsigned int iteration_index)
607 using Number =
typename VectorType::value_type;
609 const Number previous_r_dot_preconditioner_dot_r =
610 r_dot_preconditioner_dot_r;
612 if (std::is_same_v<PreconditionerType, PreconditionIdentity> ==
false)
614 preconditioner.vmult(v, r);
615 r_dot_preconditioner_dot_r = r * v;
618 r_dot_preconditioner_dot_r = residual_norm * residual_norm;
620 const VectorType &direction =
621 std::is_same_v<PreconditionerType, PreconditionIdentity> ? r : v;
623 if (iteration_index > 1)
628 r_dot_preconditioner_dot_r / previous_r_dot_preconditioner_dot_r;
630 beta -= (r * z) / previous_r_dot_preconditioner_dot_r;
631 p.sadd(beta, 1., direction);
634 p.equ(1., direction);
641 const Number p_dot_A_dot_p = p * v;
644 this->previous_alpha = alpha;
645 alpha = r_dot_preconditioner_dot_r / p_dot_A_dot_p;
652 finalize_after_convergence(
const unsigned int)
664 template <
typename MatrixType,
typename VectorType>
665 using vmult_functions_t =
decltype(std::declval<const MatrixType>().vmult(
666 std::declval<VectorType &>(),
667 std::declval<const VectorType &>(),
669 const std::function<
void(
const unsigned int,
const unsigned int)> &>(),
670 std::declval<
const std::function<
void(
const unsigned int,
671 const unsigned int)> &>()));
673 template <
typename MatrixType,
typename VectorType>
674 constexpr bool has_vmult_functions =
675 is_supported_operation<vmult_functions_t, MatrixType, VectorType>;
680 template <
typename PreconditionerType>
681 using apply_to_subrange_t =
682 decltype(std::declval<const PreconditionerType>()
683 .apply_to_subrange(0U, 0U,
nullptr,
nullptr));
685 template <
typename PreconditionerType>
686 constexpr bool has_apply_to_subrange =
687 is_supported_operation<apply_to_subrange_t, PreconditionerType>;
692 template <
typename PreconditionerType>
694 decltype(std::declval<const PreconditionerType>().apply(0U, 0.0));
696 template <
typename PreconditionerType>
697 constexpr bool has_apply =
698 is_supported_operation<apply_t, PreconditionerType>;
704 template <
typename VectorType,
706 typename PreconditionerType>
707 struct IterationWorker<
711 std::enable_if_t<has_vmult_functions<MatrixType, VectorType> &&
712 (has_apply_to_subrange<PreconditionerType> ||
713 has_apply<PreconditionerType>)&&std::
714 is_same_v<VectorType,
715 LinearAlgebra::distributed::Vector<
716 typename VectorType::value_type,
719 :
public IterationWorkerBase<VectorType, MatrixType, PreconditionerType>
721 using Number =
typename VectorType::value_type;
723 Number next_r_dot_preconditioner_dot_r;
724 Number previous_beta;
726 IterationWorker(
const MatrixType &A,
727 const PreconditionerType &preconditioner,
731 : IterationWorkerBase<VectorType, MatrixType, PreconditionerType>(
737 , next_r_dot_preconditioner_dot_r(0.)
744 do_iteration(
const unsigned int iteration_index)
746 if (iteration_index > 1)
748 previous_beta = this->beta;
749 this->beta = next_r_dot_preconditioner_dot_r /
750 this->r_dot_preconditioner_dot_r;
753 std::array<VectorizedArray<Number>, 7> vectorized_sums = {};
758 [&](
const unsigned int begin,
const unsigned int end) {
759 operation_before_loop(iteration_index, begin, end);
761 [&](
const unsigned int begin,
const unsigned int end) {
762 operation_after_loop(begin, end, vectorized_sums);
765 std::array<Number, 7> scalar_sums;
766 for (
unsigned int i = 0; i < 7; ++i)
767 scalar_sums[i] = vectorized_sums[i][0];
768 for (
unsigned int l = 1; l < VectorizedArray<Number>::size(); ++
l)
769 for (
unsigned int i = 0; i < 7; ++i)
770 scalar_sums[i] += vectorized_sums[i][l];
774 this->r.get_mpi_communicator(),
777 this->r_dot_preconditioner_dot_r = scalar_sums[6];
779 const Number p_dot_A_dot_p = scalar_sums[0];
782 this->previous_alpha = this->alpha;
783 this->alpha = this->r_dot_preconditioner_dot_r / p_dot_A_dot_p;
789 this->alpha * (-2. * scalar_sums[2] + this->alpha * scalar_sums[1])));
791 next_r_dot_preconditioner_dot_r =
std::abs(
792 this->r_dot_preconditioner_dot_r +
793 this->alpha * (-2. * scalar_sums[4] + this->alpha * scalar_sums[5]));
798 template <
typename U =
void>
799 std::enable_if_t<has_apply<PreconditionerType>, U>
800 operation_before_loop(
const unsigned int iteration_index,
801 const unsigned int start_range,
802 const unsigned int end_range)
const
804 Number *x = this->x.begin();
805 Number *r = this->r.begin();
806 Number *p = this->p.begin();
807 Number *v = this->v.begin();
808 const Number alpha = this->alpha;
809 const Number beta = this->beta;
811 const unsigned int end_regular =
812 start_range + (end_range - start_range) / n_lanes * n_lanes;
813 if (iteration_index == 1)
818 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
823 for (
unsigned int l = 0;
l < n_lanes; ++
l)
824 pj[l] = this->preconditioner.apply(j + l, rj[l]);
829 for (
unsigned int j = end_regular; j < end_range; ++j)
831 p[j] = this->preconditioner.apply(j, r[j]);
835 else if (iteration_index % 2 == 0 && beta != Number())
837 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
845 for (
unsigned int l = 0;
l < n_lanes; ++
l)
846 prec_rj[l] = this->preconditioner.apply(j + l, rj[l]);
848 pj = beta * pj + prec_rj;
853 for (
unsigned int j = end_regular; j < end_range; ++j)
855 r[j] -= alpha * v[j];
856 p[j] = beta * p[j] + this->preconditioner.apply(j, r[j]);
860 else if (iteration_index % 2 == 0 && beta == Number())
865 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
877 for (
unsigned int l = 0;
l < n_lanes; ++
l)
878 prec_rj[l] = this->preconditioner.apply(j + l, rj[l]);
879 prec_rj.
store(p + j);
883 for (
unsigned int j = end_regular; j < end_range; ++j)
885 r[j] -= alpha * v[j];
886 x[j] += alpha * p[j];
887 p[j] = this->preconditioner.apply(j, r[j]);
893 const Number alpha_plus_previous_alpha_over_beta =
894 alpha + this->previous_alpha / this->previous_beta;
895 const Number previous_alpha_over_beta =
896 this->previous_alpha / this->previous_beta;
897 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
902 xj += alpha_plus_previous_alpha_over_beta * pj;
906 for (
unsigned int l = 0;
l < n_lanes; ++
l)
908 prec_rj[
l] = this->preconditioner.apply(j + l, rj[l]);
909 prec_vj[
l] = this->preconditioner.apply(j + l, vj[l]);
911 xj -= previous_alpha_over_beta * prec_rj;
915 prec_rj -= alpha * prec_vj;
916 pj = beta * pj + prec_rj;
921 for (
unsigned int j = end_regular; j < end_range; ++j)
923 x[j] += alpha_plus_previous_alpha_over_beta * p[j];
924 x[j] -= previous_alpha_over_beta *
925 this->preconditioner.apply(j, r[j]);
926 r[j] -= alpha * v[j];
927 p[j] = beta * p[j] + this->preconditioner.apply(j, r[j]);
935 template <
typename U =
void>
936 std::enable_if_t<has_apply<PreconditionerType>, U>
937 operation_after_loop(
938 const unsigned int start_range,
939 const unsigned int end_range,
942 const Number *r = this->r.begin();
943 const Number *p = this->p.begin();
944 const Number *v = this->v.begin();
945 std::array<VectorizedArray<Number>, 7> my_sums = {};
947 const unsigned int end_regular =
948 start_range + (end_range - start_range) / n_lanes * n_lanes;
949 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
956 for (
unsigned int l = 0;
l < n_lanes; ++
l)
958 prec_vj[
l] = this->preconditioner.apply(j + l, vj[l]);
959 prec_rj[
l] = this->preconditioner.apply(j + l, rj[l]);
961 my_sums[0] += pj * vj;
962 my_sums[1] += vj * vj;
963 my_sums[2] += rj * vj;
964 my_sums[3] += rj * rj;
965 my_sums[4] += rj * prec_vj;
966 my_sums[5] += vj * prec_vj;
967 my_sums[6] += rj * prec_rj;
969 for (
unsigned int j = end_regular; j < end_range; ++j)
971 const Number prec_v = this->preconditioner.apply(j, v[j]);
972 const Number prec_r = this->preconditioner.apply(j, r[j]);
973 my_sums[0][0] += p[j] * v[j];
974 my_sums[1][0] += v[j] * v[j];
975 my_sums[2][0] += r[j] * v[j];
976 my_sums[3][0] += r[j] * r[j];
977 my_sums[4][0] += r[j] * prec_v;
978 my_sums[5][0] += v[j] * prec_v;
979 my_sums[6][0] += r[j] * prec_r;
981 for (
unsigned int i = 0; i < vectorized_sums.size(); ++i)
982 vectorized_sums[i] += my_sums[i];
987 template <
typename U =
void>
988 std::enable_if_t<has_apply<PreconditionerType>, U>
989 finalize_after_convergence(
const unsigned int iteration_index)
991 if (iteration_index % 2 == 1 || this->beta == Number())
992 this->x.add(this->alpha, this->p);
995 using Number =
typename VectorType::value_type;
996 const unsigned int end_range = this->x.locally_owned_size();
998 Number *x = this->x.begin();
999 Number *r = this->r.begin();
1000 Number *p = this->p.begin();
1006 const Number alpha_plus_previous_alpha_over_beta =
1007 this->alpha + this->previous_alpha / this->beta;
1008 const Number previous_alpha_over_beta =
1009 this->previous_alpha / this->beta;
1012 for (
unsigned int j = 0; j < end_range; ++j)
1014 x[j] += alpha_plus_previous_alpha_over_beta * p[j] -
1015 previous_alpha_over_beta *
1016 this->preconditioner.apply(j, r[j]);
1024 template <
typename U =
void>
1025 std::enable_if_t<!has_apply<PreconditionerType>, U>
1026 operation_before_loop(
const unsigned int iteration_index,
1027 const unsigned int start_range,
1028 const unsigned int end_range)
const
1030 Number *x = this->x.begin() + start_range;
1031 Number *r = this->r.begin() + start_range;
1032 Number *p = this->p.begin() + start_range;
1033 Number *v = this->v.begin() + start_range;
1034 const Number alpha = this->alpha;
1035 const Number beta = this->beta;
1036 constexpr unsigned int grain_size = 128;
1037 std::array<Number, grain_size> prec_r;
1038 if (iteration_index == 1)
1040 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1042 const unsigned int length =
std::min(grain_size, end_range - j);
1043 this->preconditioner.apply_to_subrange(j,
1048 for (
unsigned int i = 0; i < length; ++i)
1058 else if (iteration_index % 2 == 0 && beta != Number())
1060 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1062 const unsigned int length =
std::min(grain_size, end_range - j);
1064 for (
unsigned int i = 0; i < length; ++i)
1065 r[i] -= this->alpha * v[i];
1066 this->preconditioner.apply_to_subrange(j,
1071 for (
unsigned int i = 0; i < length; ++i)
1073 p[i] = this->beta * p[i] + prec_r[i];
1081 else if (iteration_index % 2 == 0 && beta == Number())
1086 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1088 const unsigned int length =
std::min(grain_size, end_range - j);
1090 for (
unsigned int i = 0; i < length; ++i)
1091 r[i] -= this->alpha * v[i];
1092 this->preconditioner.apply_to_subrange(j,
1097 for (
unsigned int i = 0; i < length; ++i)
1099 x[i] += this->alpha * p[i];
1110 const Number alpha_plus_previous_alpha_over_beta =
1111 this->alpha + this->previous_alpha / this->previous_beta;
1112 const Number previous_alpha_over_beta =
1113 this->previous_alpha / this->previous_beta;
1114 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1116 const unsigned int length =
std::min(grain_size, end_range - j);
1117 this->preconditioner.apply_to_subrange(j,
1122 for (
unsigned int i = 0; i < length; ++i)
1124 x[i] += alpha_plus_previous_alpha_over_beta * p[i] -
1125 previous_alpha_over_beta * prec_r[i];
1126 r[i] -= this->alpha * v[i];
1128 this->preconditioner.apply_to_subrange(j,
1133 for (
unsigned int i = 0; i < length; ++i)
1135 p[i] = this->beta * p[i] + prec_r[i];
1149 template <
typename U =
void>
1150 std::enable_if_t<!has_apply<PreconditionerType>, U>
1151 operation_after_loop(
1152 const unsigned int start_range,
1153 const unsigned int end_range,
1156 const Number *r = this->r.begin();
1157 const Number *p = this->p.begin();
1158 const Number *v = this->v.begin();
1159 std::array<VectorizedArray<Number>, 7> my_sums = {};
1160 constexpr unsigned int grain_size = 128;
1163 const unsigned int end_regular =
1164 start_range + (end_range - start_range) / grain_size * grain_size;
1165 std::array<Number, grain_size> prec_r;
1166 std::array<Number, grain_size> prec_v;
1167 for (
unsigned int j = start_range; j < end_regular; j += grain_size)
1169 this->preconditioner.apply_to_subrange(j,
1173 this->preconditioner.apply_to_subrange(j,
1178 for (
unsigned int i = 0; i < grain_size;
1184 prec_rj.
load(prec_r.data() + i);
1185 prec_vj.
load(prec_v.data() + i);
1187 my_sums[0] += pj * vj;
1188 my_sums[1] += vj * vj;
1189 my_sums[2] += rj * vj;
1190 my_sums[3] += rj * rj;
1191 my_sums[4] += rj * prec_vj;
1192 my_sums[5] += vj * prec_vj;
1193 my_sums[6] += rj * prec_rj;
1196 const unsigned int length = end_range - end_regular;
1198 this->preconditioner.apply_to_subrange(end_regular,
1199 end_regular + length,
1202 this->preconditioner.apply_to_subrange(end_regular,
1203 end_regular + length,
1206 for (
unsigned int j = end_regular; j < end_range; ++j)
1208 my_sums[0][0] += p[j] * v[j];
1209 my_sums[1][0] += v[j] * v[j];
1210 my_sums[2][0] += r[j] * v[j];
1211 my_sums[3][0] += r[j] * r[j];
1212 my_sums[4][0] += r[j] * prec_v[j - end_regular];
1213 my_sums[5][0] += v[j] * prec_v[j - end_regular];
1214 my_sums[6][0] += r[j] * prec_r[j - end_regular];
1216 for (
unsigned int i = 0; i < vectorized_sums.size(); ++i)
1217 vectorized_sums[i] += my_sums[i];
1223 template <
typename U =
void>
1224 std::enable_if_t<!has_apply<PreconditionerType>, U>
1225 finalize_after_convergence(
const unsigned int iteration_index)
1227 if (iteration_index % 2 == 1 || this->beta == Number())
1228 this->x.add(this->alpha, this->p);
1231 const unsigned int end_range = this->x.locally_owned_size();
1233 Number *x = this->x.begin();
1234 Number *r = this->r.begin();
1235 Number *p = this->p.begin();
1236 const Number alpha_plus_previous_alpha_over_beta =
1237 this->alpha + this->previous_alpha / this->beta;
1238 const Number previous_alpha_over_beta =
1239 this->previous_alpha / this->beta;
1241 constexpr unsigned int grain_size = 128;
1242 std::array<Number, grain_size> prec_r;
1243 for (
unsigned int j = 0; j < end_range; j += grain_size)
1245 const unsigned int length =
std::min(grain_size, end_range - j);
1246 this->preconditioner.apply_to_subrange(j,
1251 for (
unsigned int i = 0; i < length; ++i)
1252 x[i] += alpha_plus_previous_alpha_over_beta * p[i] -
1253 previous_alpha_over_beta * prec_r[i];
1267template <
typename VectorType>
1269template <
typename MatrixType,
typename PreconditionerType>
1275 const VectorType &b,
1276 const PreconditionerType &preconditioner)
1278 using number =
typename VectorType::value_type;
1285 const bool do_eigenvalues =
1286 !condition_number_signal.empty() || !all_condition_numbers_signal.empty() ||
1287 !eigenvalues_signal.empty() || !all_eigenvalues_signal.empty();
1290 std::vector<typename VectorType::value_type>
diagonal;
1291 std::vector<typename VectorType::value_type> offdiagonal;
1293 typename VectorType::value_type eigen_beta_alpha = 0;
1297 internal::SolverCG::
1298 IterationWorker<VectorType, MatrixType, PreconditionerType>
1300 A, preconditioner, determine_beta_by_flexible_formula, this->memory, x);
1304 solver_state = this->iteration_status(0, worker.residual_norm, x);
1312 worker.do_iteration(it);
1314 print_vectors(it, x, worker.r, worker.p);
1318 this->coefficients_signal(worker.previous_alpha, worker.beta);
1323 diagonal.push_back(number(1.) / worker.previous_alpha +
1325 eigen_beta_alpha = worker.beta / worker.previous_alpha;
1326 offdiagonal.push_back(
std::sqrt(worker.beta) /
1327 worker.previous_alpha);
1329 compute_eigs_and_cond(diagonal,
1331 all_eigenvalues_signal,
1332 all_condition_numbers_signal);
1335 solver_state = this->iteration_status(it, worker.residual_norm, x);
1338 worker.finalize_after_convergence(it);
1340 compute_eigs_and_cond(diagonal,
1343 condition_number_signal);
1351template <
typename VectorType>
1354 const std::function<
void(
typename VectorType::value_type,
1355 typename VectorType::value_type)> &slot)
1357 return coefficients_signal.connect(slot);
1362template <
typename VectorType>
1365 const std::function<
void(
double)> &slot,
1366 const bool every_iteration)
1368 if (every_iteration)
1370 return all_condition_numbers_signal.connect(slot);
1374 return condition_number_signal.connect(slot);
1380template <
typename VectorType>
1383 const std::function<
void(
const std::vector<double> &)> &slot,
1384 const bool every_iteration)
1386 if (every_iteration)
1388 return all_eigenvalues_signal.connect(slot);
1392 return eigenvalues_signal.connect(slot);
1398template <
typename VectorType>
1402 const AdditionalData &)
1405 this->determine_beta_by_flexible_formula =
true;
1410template <
typename VectorType>
1413 const AdditionalData &)
1416 this->determine_beta_by_flexible_formula =
true;
SolverCG(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
boost::signals2::connection connect_condition_number_slot(const std::function< void(double)> &slot, const bool every_iteration=false)
boost::signals2::connection connect_eigenvalues_slot(const std::function< void(const std::vector< double > &)> &slot, const bool every_iteration=false)
boost::signals2::connection connect_coefficients_slot(const std::function< void(typename VectorType::value_type, typename VectorType::value_type)> &slot)
virtual ~SolverCG() override=default
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
static void compute_eigs_and_cond(const std::vector< typename VectorType::value_type > &diagonal, const std::vector< typename VectorType::value_type > &offdiagonal, const boost::signals2::signal< void(const std::vector< double > &)> &eigenvalues_signal, const boost::signals2::signal< void(double)> &cond_signal)
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
SolverFlexibleCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
SolverFlexibleCG(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
static constexpr std::size_t size()
void store(OtherNumber *ptr) const
void load(const OtherNumber *ptr)
#define DEAL_II_OPENMP_SIMD_PRAGMA
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcDivideByZero()
#define AssertThrow(cond, exc)
@ diagonal
Matrix is diagonal.
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)