16#ifndef dealii_solver_cg_h
17#define dealii_solver_cg_h
42 template <
typename,
typename>
176template <
typename VectorType = Vector<
double>>
214 template <
typename MatrixType,
typename PreconditionerType>
218 const VectorType & b,
219 const PreconditionerType &preconditioner);
227 boost::signals2::connection
229 const std::function<
void(
typename VectorType::value_type,
230 typename VectorType::value_type)> &slot);
238 boost::signals2::connection
240 const bool every_iteration =
false);
248 boost::signals2::connection
250 const std::function<
void(
const std::vector<double> &)> &slot,
251 const bool every_iteration =
false);
261 const VectorType & x,
262 const VectorType & r,
263 const VectorType & d)
const;
272 const std::vector<typename VectorType::value_type> &diagonal,
273 const std::vector<typename VectorType::value_type> &offdiagonal,
274 const boost::signals2::signal<
void(
const std::vector<double> &)>
276 const boost::signals2::signal<
void(
double)> &cond_signal);
286 boost::signals2::signal<void(
typename VectorType::value_type,
287 typename VectorType::value_type)>
312 boost::signals2::signal<void(
const std::vector<double> &)>
353template <
typename VectorType = Vector<
double>>
394template <
typename VectorType>
397 const AdditionalData & data)
399 , additional_data(data)
400 , determine_beta_by_flexible_formula(false)
405template <
typename VectorType>
408 , additional_data(data)
409 , determine_beta_by_flexible_formula(false)
414template <
typename VectorType>
419 const VectorType &)
const
424template <
typename VectorType>
427 const std::vector<typename VectorType::value_type> &diagonal,
428 const std::vector<typename VectorType::value_type> &offdiagonal,
429 const boost::signals2::signal<
void(
const std::vector<double> &)>
430 & eigenvalues_signal,
431 const boost::signals2::signal<
void(
double)> &cond_signal)
434 if (!cond_signal.empty() || !eigenvalues_signal.empty())
438 for (size_type i = 0; i <
diagonal.size(); ++i)
442 T(i, i + 1) = offdiagonal[i];
444 T.compute_eigenvalues();
448 auto condition_number = T.eigenvalue(T.n() - 1) / T.eigenvalue(0);
451 cond_signal(
std::abs(condition_number));
455 if (!eigenvalues_signal.empty())
458 for (
unsigned int j = 0; j < T.n(); ++j)
481 template <
typename VectorType,
483 typename PreconditionerType>
484 struct IterationWorkerBase
486 using Number =
typename VectorType::value_type;
488 const MatrixType & A;
489 const PreconditionerType &preconditioner;
509 Number r_dot_preconditioner_dot_r;
512 double residual_norm;
513 Number previous_alpha;
515 IterationWorkerBase(
const MatrixType & A,
516 const PreconditionerType &preconditioner,
521 , preconditioner(preconditioner)
532 , r_dot_preconditioner_dot_r(Number())
536 , previous_alpha(Number())
540 startup(
const VectorType &b)
560 residual_norm = r.l2_norm();
568 template <
typename VectorType,
570 typename PreconditionerType,
572 struct IterationWorker
573 :
public IterationWorkerBase<VectorType, MatrixType, PreconditionerType>
576 IterationWorkerBase<VectorType, MatrixType, PreconditionerType>;
578 IterationWorker(
const MatrixType & A,
579 const PreconditionerType &preconditioner,
583 : BaseClass(A, preconditioner, flexible, memory, x)
587 using BaseClass::alpha;
588 using BaseClass::beta;
590 using BaseClass::preconditioner;
592 using BaseClass::r_dot_preconditioner_dot_r;
593 using BaseClass::residual_norm;
599 do_iteration(
const unsigned int iteration_index)
601 using Number =
typename VectorType::value_type;
603 const Number previous_r_dot_preconditioner_dot_r =
604 r_dot_preconditioner_dot_r;
606 if (std::is_same<PreconditionerType, PreconditionIdentity>::value ==
609 preconditioner.vmult(v, r);
610 r_dot_preconditioner_dot_r = r * v;
613 r_dot_preconditioner_dot_r = residual_norm * residual_norm;
615 const VectorType &direction =
616 std::is_same<PreconditionerType, PreconditionIdentity>::value ? r : v;
618 if (iteration_index > 1)
623 r_dot_preconditioner_dot_r / previous_r_dot_preconditioner_dot_r;
625 beta -= (r * z) / previous_r_dot_preconditioner_dot_r;
626 p.sadd(beta, 1., direction);
629 p.equ(1., direction);
636 const Number p_dot_A_dot_p = p * v;
639 this->previous_alpha = alpha;
640 alpha = r_dot_preconditioner_dot_r / p_dot_A_dot_p;
647 finalize_after_convergence(
const unsigned int)
659 template <
typename MatrixType,
typename VectorType>
660 using vmult_functions_t =
decltype(std::declval<MatrixType const>().vmult(
661 std::declval<VectorType &>(),
662 std::declval<const VectorType &>(),
664 const std::function<
void(
const unsigned int,
const unsigned int)> &>(),
665 std::declval<
const std::function<
void(
const unsigned int,
666 const unsigned int)> &>()));
668 template <
typename MatrixType,
typename VectorType>
669 constexpr bool has_vmult_functions =
670 is_supported_operation<vmult_functions_t, MatrixType, VectorType>;
675 template <
typename PreconditionerType>
676 using apply_to_subrange_t =
677 decltype(std::declval<PreconditionerType const>()
678 .apply_to_subrange(0U, 0U,
nullptr,
nullptr));
680 template <
typename PreconditionerType>
681 constexpr bool has_apply_to_subrange =
682 is_supported_operation<apply_to_subrange_t, PreconditionerType>;
687 template <
typename PreconditionerType>
689 decltype(std::declval<PreconditionerType const>().apply(0U, 0.0));
691 template <
typename PreconditionerType>
692 constexpr bool has_apply =
693 is_supported_operation<apply_t, PreconditionerType>;
699 template <
typename VectorType,
701 typename PreconditionerType>
702 struct IterationWorker<
706 std::enable_if_t<has_vmult_functions<MatrixType, VectorType> &&
707 (has_apply_to_subrange<PreconditionerType> ||
708 has_apply<PreconditionerType>)&&std::
710 LinearAlgebra::distributed::Vector<
711 typename VectorType::value_type,
712 MemorySpace::Host>>::value,
714 :
public IterationWorkerBase<VectorType, MatrixType, PreconditionerType>
716 using Number =
typename VectorType::value_type;
718 Number next_r_dot_preconditioner_dot_r;
719 Number previous_beta;
721 IterationWorker(
const MatrixType & A,
722 const PreconditionerType &preconditioner,
726 : IterationWorkerBase<VectorType, MatrixType, PreconditionerType>(
732 , next_r_dot_preconditioner_dot_r(0.)
739 do_iteration(
const unsigned int iteration_index)
741 if (iteration_index > 1)
743 previous_beta = this->beta;
744 this->beta = next_r_dot_preconditioner_dot_r /
745 this->r_dot_preconditioner_dot_r;
748 std::array<VectorizedArray<Number>, 7> vectorized_sums = {};
753 [&](
const unsigned int begin,
const unsigned int end) {
754 operation_before_loop(iteration_index, begin, end);
756 [&](
const unsigned int begin,
const unsigned int end) {
757 operation_after_loop(begin, end, vectorized_sums);
760 std::array<Number, 7> scalar_sums;
761 for (
unsigned int i = 0; i < 7; ++i)
762 scalar_sums[i] = vectorized_sums[i][0];
763 for (
unsigned int l = 1; l < VectorizedArray<Number>::size(); ++
l)
764 for (
unsigned int i = 0; i < 7; ++i)
765 scalar_sums[i] += vectorized_sums[i][l];
769 this->r.get_mpi_communicator(),
772 this->r_dot_preconditioner_dot_r = scalar_sums[6];
774 const Number p_dot_A_dot_p = scalar_sums[0];
777 this->previous_alpha = this->alpha;
778 this->alpha = this->r_dot_preconditioner_dot_r / p_dot_A_dot_p;
784 this->alpha * (-2. * scalar_sums[2] + this->alpha * scalar_sums[1])));
786 next_r_dot_preconditioner_dot_r =
std::abs(
787 this->r_dot_preconditioner_dot_r +
788 this->alpha * (-2. * scalar_sums[4] + this->alpha * scalar_sums[5]));
793 template <
typename U =
void>
794 std::enable_if_t<has_apply<PreconditionerType>, U>
795 operation_before_loop(
const unsigned int iteration_index,
796 const unsigned int start_range,
797 const unsigned int end_range)
const
799 Number * x = this->x.begin();
800 Number * r = this->r.begin();
801 Number * p = this->p.begin();
802 Number * v = this->v.begin();
803 const Number alpha = this->alpha;
804 const Number beta = this->beta;
806 const unsigned int end_regular =
807 start_range + (end_range - start_range) / n_lanes * n_lanes;
808 if (iteration_index == 1)
813 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
818 for (
unsigned int l = 0;
l < n_lanes; ++
l)
819 pj[l] = this->preconditioner.apply(j + l, rj[l]);
824 for (
unsigned int j = end_regular; j < end_range; ++j)
826 p[j] = this->preconditioner.apply(j, r[j]);
830 else if (iteration_index % 2 == 0 && beta != Number())
832 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
840 for (
unsigned int l = 0;
l < n_lanes; ++
l)
841 prec_rj[l] = this->preconditioner.apply(j + l, rj[l]);
843 pj = beta * pj + prec_rj;
848 for (
unsigned int j = end_regular; j < end_range; ++j)
850 r[j] -= alpha * v[j];
851 p[j] = beta * p[j] + this->preconditioner.apply(j, r[j]);
855 else if (iteration_index % 2 == 0 && beta == Number())
860 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
872 for (
unsigned int l = 0;
l < n_lanes; ++
l)
873 prec_rj[l] = this->preconditioner.apply(j + l, rj[l]);
874 prec_rj.
store(p + j);
878 for (
unsigned int j = end_regular; j < end_range; ++j)
880 r[j] -= alpha * v[j];
881 x[j] += alpha * p[j];
882 p[j] = this->preconditioner.apply(j, r[j]);
888 const Number alpha_plus_previous_alpha_over_beta =
889 alpha + this->previous_alpha / this->previous_beta;
890 const Number previous_alpha_over_beta =
891 this->previous_alpha / this->previous_beta;
892 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
897 xj += alpha_plus_previous_alpha_over_beta * pj;
901 for (
unsigned int l = 0;
l < n_lanes; ++
l)
903 prec_rj[
l] = this->preconditioner.apply(j + l, rj[l]);
904 prec_vj[
l] = this->preconditioner.apply(j + l, vj[l]);
906 xj -= previous_alpha_over_beta * prec_rj;
910 prec_rj -= alpha * prec_vj;
911 pj = beta * pj + prec_rj;
916 for (
unsigned int j = end_regular; j < end_range; ++j)
918 x[j] += alpha_plus_previous_alpha_over_beta * p[j];
919 x[j] -= previous_alpha_over_beta *
920 this->preconditioner.apply(j, r[j]);
921 r[j] -= alpha * v[j];
922 p[j] = beta * p[j] + this->preconditioner.apply(j, r[j]);
930 template <
typename U =
void>
931 std::enable_if_t<has_apply<PreconditionerType>, U>
932 operation_after_loop(
933 const unsigned int start_range,
934 const unsigned int end_range,
937 const Number * r = this->r.begin();
938 const Number * p = this->p.begin();
939 const Number * v = this->v.begin();
940 std::array<VectorizedArray<Number>, 7> my_sums = {};
942 const unsigned int end_regular =
943 start_range + (end_range - start_range) / n_lanes * n_lanes;
944 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
951 for (
unsigned int l = 0;
l < n_lanes; ++
l)
953 prec_vj[
l] = this->preconditioner.apply(j + l, vj[l]);
954 prec_rj[
l] = this->preconditioner.apply(j + l, rj[l]);
956 my_sums[0] += pj * vj;
957 my_sums[1] += vj * vj;
958 my_sums[2] += rj * vj;
959 my_sums[3] += rj * rj;
960 my_sums[4] += rj * prec_vj;
961 my_sums[5] += vj * prec_vj;
962 my_sums[6] += rj * prec_rj;
964 for (
unsigned int j = end_regular; j < end_range; ++j)
966 const Number prec_v = this->preconditioner.apply(j, v[j]);
967 const Number prec_r = this->preconditioner.apply(j, r[j]);
968 my_sums[0][0] += p[j] * v[j];
969 my_sums[1][0] += v[j] * v[j];
970 my_sums[2][0] += r[j] * v[j];
971 my_sums[3][0] += r[j] * r[j];
972 my_sums[4][0] += r[j] * prec_v;
973 my_sums[5][0] += v[j] * prec_v;
974 my_sums[6][0] += r[j] * prec_r;
976 for (
unsigned int i = 0; i < vectorized_sums.size(); ++i)
977 vectorized_sums[i] += my_sums[i];
982 template <
typename U =
void>
983 std::enable_if_t<has_apply<PreconditionerType>, U>
984 finalize_after_convergence(
const unsigned int iteration_index)
986 if (iteration_index % 2 == 1 || this->beta == Number())
987 this->x.add(this->alpha, this->p);
990 using Number =
typename VectorType::value_type;
991 const unsigned int end_range = this->x.locally_owned_size();
993 Number *x = this->x.
begin();
994 Number *r = this->r.begin();
995 Number *p = this->p.begin();
1001 const Number alpha_plus_previous_alpha_over_beta =
1002 this->alpha + this->previous_alpha / this->beta;
1003 const Number previous_alpha_over_beta =
1004 this->previous_alpha / this->beta;
1007 for (
unsigned int j = 0; j < end_range; ++j)
1009 x[j] += alpha_plus_previous_alpha_over_beta * p[j] -
1010 previous_alpha_over_beta *
1011 this->preconditioner.apply(j, r[j]);
1019 template <
typename U =
void>
1020 std::enable_if_t<!has_apply<PreconditionerType>, U>
1021 operation_before_loop(
const unsigned int iteration_index,
1022 const unsigned int start_range,
1023 const unsigned int end_range)
const
1025 Number * x = this->x.begin() + start_range;
1026 Number * r = this->r.begin() + start_range;
1027 Number * p = this->p.begin() + start_range;
1028 Number * v = this->v.begin() + start_range;
1029 const Number alpha = this->alpha;
1030 const Number beta = this->beta;
1031 constexpr unsigned int grain_size = 128;
1032 std::array<Number, grain_size> prec_r;
1033 if (iteration_index == 1)
1035 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1037 const unsigned int length =
std::min(grain_size, end_range - j);
1038 this->preconditioner.apply_to_subrange(j,
1043 for (
unsigned int i = 0; i < length; ++i)
1053 else if (iteration_index % 2 == 0 && beta != Number())
1055 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1057 const unsigned int length =
std::min(grain_size, end_range - j);
1059 for (
unsigned int i = 0; i < length; ++i)
1060 r[i] -= this->alpha * v[i];
1061 this->preconditioner.apply_to_subrange(j,
1066 for (
unsigned int i = 0; i < length; ++i)
1068 p[i] = this->beta * p[i] + prec_r[i];
1076 else if (iteration_index % 2 == 0 && beta == Number())
1081 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1083 const unsigned int length =
std::min(grain_size, end_range - j);
1085 for (
unsigned int i = 0; i < length; ++i)
1086 r[i] -= this->alpha * v[i];
1087 this->preconditioner.apply_to_subrange(j,
1092 for (
unsigned int i = 0; i < length; ++i)
1094 x[i] += this->alpha * p[i];
1105 const Number alpha_plus_previous_alpha_over_beta =
1106 this->alpha + this->previous_alpha / this->previous_beta;
1107 const Number previous_alpha_over_beta =
1108 this->previous_alpha / this->previous_beta;
1109 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1111 const unsigned int length =
std::min(grain_size, end_range - j);
1112 this->preconditioner.apply_to_subrange(j,
1117 for (
unsigned int i = 0; i < length; ++i)
1119 x[i] += alpha_plus_previous_alpha_over_beta * p[i] -
1120 previous_alpha_over_beta * prec_r[i];
1121 r[i] -= this->alpha * v[i];
1123 this->preconditioner.apply_to_subrange(j,
1128 for (
unsigned int i = 0; i < length; ++i)
1130 p[i] = this->beta * p[i] + prec_r[i];
1144 template <
typename U =
void>
1145 std::enable_if_t<!has_apply<PreconditionerType>, U>
1146 operation_after_loop(
1147 const unsigned int start_range,
1148 const unsigned int end_range,
1151 const Number * r = this->r.begin();
1152 const Number * p = this->p.begin();
1153 const Number * v = this->v.begin();
1154 std::array<VectorizedArray<Number>, 7> my_sums = {};
1155 constexpr unsigned int grain_size = 128;
1158 const unsigned int end_regular =
1159 start_range + (end_range - start_range) / grain_size * grain_size;
1160 std::array<Number, grain_size> prec_r;
1161 std::array<Number, grain_size> prec_v;
1162 for (
unsigned int j = start_range; j < end_regular; j += grain_size)
1164 this->preconditioner.apply_to_subrange(j,
1168 this->preconditioner.apply_to_subrange(j,
1173 for (
unsigned int i = 0; i < grain_size;
1179 prec_rj.
load(prec_r.data() + i);
1180 prec_vj.
load(prec_v.data() + i);
1182 my_sums[0] += pj * vj;
1183 my_sums[1] += vj * vj;
1184 my_sums[2] += rj * vj;
1185 my_sums[3] += rj * rj;
1186 my_sums[4] += rj * prec_vj;
1187 my_sums[5] += vj * prec_vj;
1188 my_sums[6] += rj * prec_rj;
1191 const unsigned int length = end_range - end_regular;
1193 this->preconditioner.apply_to_subrange(end_regular,
1194 end_regular + length,
1197 this->preconditioner.apply_to_subrange(end_regular,
1198 end_regular + length,
1201 for (
unsigned int j = end_regular; j < end_range; ++j)
1203 my_sums[0][0] += p[j] * v[j];
1204 my_sums[1][0] += v[j] * v[j];
1205 my_sums[2][0] += r[j] * v[j];
1206 my_sums[3][0] += r[j] * r[j];
1207 my_sums[4][0] += r[j] * prec_v[j - end_regular];
1208 my_sums[5][0] += v[j] * prec_v[j - end_regular];
1209 my_sums[6][0] += r[j] * prec_r[j - end_regular];
1211 for (
unsigned int i = 0; i < vectorized_sums.size(); ++i)
1212 vectorized_sums[i] += my_sums[i];
1218 template <
typename U =
void>
1219 std::enable_if_t<!has_apply<PreconditionerType>, U>
1220 finalize_after_convergence(
const unsigned int iteration_index)
1222 if (iteration_index % 2 == 1 || this->beta == Number())
1223 this->x.add(this->alpha, this->p);
1226 const unsigned int end_range = this->x.locally_owned_size();
1228 Number * x = this->x.begin();
1229 Number * r = this->r.begin();
1230 Number * p = this->p.begin();
1231 const Number alpha_plus_previous_alpha_over_beta =
1232 this->alpha + this->previous_alpha / this->beta;
1233 const Number previous_alpha_over_beta =
1234 this->previous_alpha / this->beta;
1236 constexpr unsigned int grain_size = 128;
1237 std::array<Number, grain_size> prec_r;
1238 for (
unsigned int j = 0; j < end_range; j += grain_size)
1240 const unsigned int length =
std::min(grain_size, end_range - j);
1241 this->preconditioner.apply_to_subrange(j,
1246 for (
unsigned int i = 0; i < length; ++i)
1247 x[i] += alpha_plus_previous_alpha_over_beta * p[i] -
1248 previous_alpha_over_beta * prec_r[i];
1262template <
typename VectorType>
1263template <
typename MatrixType,
typename PreconditionerType>
1267 const VectorType & b,
1268 const PreconditionerType &preconditioner)
1270 using number =
typename VectorType::value_type;
1277 const bool do_eigenvalues =
1278 !condition_number_signal.empty() || !all_condition_numbers_signal.empty() ||
1279 !eigenvalues_signal.empty() || !all_eigenvalues_signal.empty();
1282 std::vector<typename VectorType::value_type>
diagonal;
1283 std::vector<typename VectorType::value_type> offdiagonal;
1285 typename VectorType::value_type eigen_beta_alpha = 0;
1289 internal::SolverCG::
1290 IterationWorker<VectorType, MatrixType, PreconditionerType>
1292 A, preconditioner, determine_beta_by_flexible_formula, this->memory, x);
1296 solver_state = this->iteration_status(0, worker.residual_norm, x);
1304 worker.do_iteration(it);
1306 print_vectors(it, x, worker.r, worker.p);
1310 this->coefficients_signal(worker.previous_alpha, worker.beta);
1315 diagonal.push_back(number(1.) / worker.previous_alpha +
1317 eigen_beta_alpha = worker.beta / worker.previous_alpha;
1318 offdiagonal.push_back(
std::sqrt(worker.beta) /
1319 worker.previous_alpha);
1321 compute_eigs_and_cond(diagonal,
1323 all_eigenvalues_signal,
1324 all_condition_numbers_signal);
1327 solver_state = this->iteration_status(it, worker.residual_norm, x);
1330 worker.finalize_after_convergence(it);
1332 compute_eigs_and_cond(diagonal,
1335 condition_number_signal);
1343template <
typename VectorType>
1344boost::signals2::connection
1346 const std::function<
void(
typename VectorType::value_type,
1347 typename VectorType::value_type)> &slot)
1349 return coefficients_signal.connect(slot);
1354template <
typename VectorType>
1355boost::signals2::connection
1357 const std::function<
void(
double)> &slot,
1358 const bool every_iteration)
1360 if (every_iteration)
1362 return all_condition_numbers_signal.connect(slot);
1366 return condition_number_signal.connect(slot);
1372template <
typename VectorType>
1373boost::signals2::connection
1375 const std::function<
void(
const std::vector<double> &)> &slot,
1376 const bool every_iteration)
1378 if (every_iteration)
1380 return all_eigenvalues_signal.connect(slot);
1384 return eigenvalues_signal.connect(slot);
1390template <
typename VectorType>
1393 const AdditionalData &)
1401template <
typename VectorType>
1403 const AdditionalData &)
SolverCG(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
boost::signals2::signal< void(const std::vector< double > &)> eigenvalues_signal
bool determine_beta_by_flexible_formula
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
boost::signals2::signal< void(typename VectorType::value_type, typename VectorType::value_type)> coefficients_signal
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::signal< void(double)> condition_number_signal
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)
boost::signals2::signal< void(double)> all_condition_numbers_signal
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
boost::signals2::signal< void(const std::vector< double > &)> all_eigenvalues_signal
AdditionalData additional_data
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()
VectorizedArrayIterator< T > begin()
void store(OtherNumber *ptr) const
void load(const OtherNumber *ptr)
#define DEAL_II_OPENMP_SIMD_PRAGMA
#define DEAL_II_NAMESPACE_OPEN
#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)