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 typename
std::enable_if<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 typename std::enable_if<has_apply<PreconditionerType>, U>::type
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)
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]);
857 const Number alpha_plus_previous_alpha_over_beta =
858 alpha + this->previous_alpha / this->previous_beta;
859 const Number previous_alpha_over_beta =
860 this->previous_alpha / this->previous_beta;
861 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
866 xj += alpha_plus_previous_alpha_over_beta * pj;
870 for (
unsigned int l = 0;
l < n_lanes; ++
l)
872 prec_rj[
l] = this->preconditioner.apply(j + l, rj[l]);
873 prec_vj[
l] = this->preconditioner.apply(j + l, vj[l]);
875 xj -= previous_alpha_over_beta * prec_rj;
879 prec_rj -= alpha * prec_vj;
880 pj = beta * pj + prec_rj;
885 for (
unsigned int j = end_regular; j < end_range; ++j)
887 x[j] += alpha_plus_previous_alpha_over_beta * p[j];
888 x[j] -= previous_alpha_over_beta *
889 this->preconditioner.apply(j, r[j]);
890 r[j] -= alpha * v[j];
891 p[j] = beta * p[j] + this->preconditioner.apply(j, r[j]);
899 template <
typename U =
void>
900 typename std::enable_if<has_apply<PreconditionerType>, U>::type
901 operation_after_loop(
902 const unsigned int start_range,
903 const unsigned int end_range,
906 const Number * r = this->r.begin();
907 const Number * p = this->p.begin();
908 const Number * v = this->v.begin();
909 std::array<VectorizedArray<Number>, 7> my_sums = {};
911 const unsigned int end_regular =
912 start_range + (end_range - start_range) / n_lanes * n_lanes;
913 for (
unsigned int j = start_range; j < end_regular; j += n_lanes)
920 for (
unsigned int l = 0;
l < n_lanes; ++
l)
922 prec_vj[
l] = this->preconditioner.apply(j + l, vj[l]);
923 prec_rj[
l] = this->preconditioner.apply(j + l, rj[l]);
925 my_sums[0] += pj * vj;
926 my_sums[1] += vj * vj;
927 my_sums[2] += rj * vj;
928 my_sums[3] += rj * rj;
929 my_sums[4] += rj * prec_vj;
930 my_sums[5] += vj * prec_vj;
931 my_sums[6] += rj * prec_rj;
933 for (
unsigned int j = end_regular; j < end_range; ++j)
935 const Number prec_v = this->preconditioner.apply(j, v[j]);
936 const Number prec_r = this->preconditioner.apply(j, r[j]);
937 my_sums[0][0] += p[j] * v[j];
938 my_sums[1][0] += v[j] * v[j];
939 my_sums[2][0] += r[j] * v[j];
940 my_sums[3][0] += r[j] * r[j];
941 my_sums[4][0] += r[j] * prec_v;
942 my_sums[5][0] += v[j] * prec_v;
943 my_sums[6][0] += r[j] * prec_r;
945 for (
unsigned int i = 0; i < vectorized_sums.size(); ++i)
946 vectorized_sums[i] += my_sums[i];
951 template <
typename U =
void>
952 typename std::enable_if<has_apply<PreconditionerType>, U>::type
953 finalize_after_convergence(
const unsigned int iteration_index)
955 if (iteration_index % 2 == 1 || iteration_index == 2)
956 this->x.add(this->alpha, this->p);
959 using Number =
typename VectorType::value_type;
960 const unsigned int end_range = this->x.locally_owned_size();
962 Number *x = this->x.
begin();
963 Number *r = this->r.begin();
964 Number *p = this->p.begin();
970 const Number alpha_plus_previous_alpha_over_beta =
971 this->alpha + this->previous_alpha / this->previous_beta;
972 const Number previous_alpha_over_beta =
973 this->previous_alpha / this->previous_beta;
976 for (
unsigned int j = 0; j < end_range; ++j)
978 x[j] += alpha_plus_previous_alpha_over_beta * p[j] -
979 previous_alpha_over_beta *
980 this->preconditioner.apply(j, r[j]);
988 template <
typename U =
void>
989 typename std::enable_if<!has_apply<PreconditionerType>, U>::type
990 operation_before_loop(
const unsigned int iteration_index,
991 const unsigned int start_range,
992 const unsigned int end_range)
const
994 Number * x = this->x.begin() + start_range;
995 Number * r = this->r.begin() + start_range;
996 Number * p = this->p.begin() + start_range;
997 Number * v = this->v.begin() + start_range;
998 const Number alpha = this->alpha;
999 const Number beta = this->beta;
1000 constexpr unsigned int grain_size = 128;
1001 std::array<Number, grain_size> prec_r;
1002 if (iteration_index == 1)
1004 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1006 const unsigned int length =
std::min(grain_size, end_range - j);
1007 this->preconditioner.apply_to_subrange(j,
1012 for (
unsigned int i = 0; i < length; ++i)
1022 else if (iteration_index % 2 == 0)
1024 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1026 const unsigned int length =
std::min(grain_size, end_range - j);
1028 for (
unsigned int i = 0; i < length; ++i)
1029 r[i] -= this->alpha * v[i];
1030 this->preconditioner.apply_to_subrange(j,
1035 for (
unsigned int i = 0; i < length; ++i)
1037 p[i] = this->beta * p[i] + prec_r[i];
1047 const Number alpha_plus_previous_alpha_over_beta =
1048 this->alpha + this->previous_alpha / this->previous_beta;
1049 const Number previous_alpha_over_beta =
1050 this->previous_alpha / this->previous_beta;
1051 for (
unsigned int j = start_range; j < end_range; j += grain_size)
1053 const unsigned int length =
std::min(grain_size, end_range - j);
1054 this->preconditioner.apply_to_subrange(j,
1059 for (
unsigned int i = 0; i < length; ++i)
1061 x[i] += alpha_plus_previous_alpha_over_beta * p[i] -
1062 previous_alpha_over_beta * prec_r[i];
1063 r[i] -= this->alpha * v[i];
1065 this->preconditioner.apply_to_subrange(j,
1070 for (
unsigned int i = 0; i < length; ++i)
1072 p[i] = this->beta * p[i] + prec_r[i];
1086 template <
typename U =
void>
1087 typename std::enable_if<!has_apply<PreconditionerType>, U>::type
1088 operation_after_loop(
1089 const unsigned int start_range,
1090 const unsigned int end_range,
1093 const Number * r = this->r.begin();
1094 const Number * p = this->p.begin();
1095 const Number * v = this->v.begin();
1096 std::array<VectorizedArray<Number>, 7> my_sums = {};
1097 constexpr unsigned int grain_size = 128;
1100 const unsigned int end_regular =
1101 start_range + (end_range - start_range) / grain_size * grain_size;
1102 std::array<Number, grain_size> prec_r;
1103 std::array<Number, grain_size> prec_v;
1104 for (
unsigned int j = start_range; j < end_regular; j += grain_size)
1106 this->preconditioner.apply_to_subrange(j,
1110 this->preconditioner.apply_to_subrange(j,
1115 for (
unsigned int i = 0; i < grain_size;
1121 prec_rj.
load(prec_r.data() + i);
1122 prec_vj.
load(prec_v.data() + i);
1124 my_sums[0] += pj * vj;
1125 my_sums[1] += vj * vj;
1126 my_sums[2] += rj * vj;
1127 my_sums[3] += rj * rj;
1128 my_sums[4] += rj * prec_vj;
1129 my_sums[5] += vj * prec_vj;
1130 my_sums[6] += rj * prec_rj;
1133 const unsigned int length = end_range - end_regular;
1135 this->preconditioner.apply_to_subrange(end_regular,
1136 end_regular + length,
1139 this->preconditioner.apply_to_subrange(end_regular,
1140 end_regular + length,
1143 for (
unsigned int j = end_regular; j < end_range; ++j)
1145 my_sums[0][0] += p[j] * v[j];
1146 my_sums[1][0] += v[j] * v[j];
1147 my_sums[2][0] += r[j] * v[j];
1148 my_sums[3][0] += r[j] * r[j];
1149 my_sums[4][0] += r[j] * prec_v[j - end_regular];
1150 my_sums[5][0] += v[j] * prec_v[j - end_regular];
1151 my_sums[6][0] += r[j] * prec_r[j - end_regular];
1153 for (
unsigned int i = 0; i < vectorized_sums.size(); ++i)
1154 vectorized_sums[i] += my_sums[i];
1160 template <
typename U =
void>
1161 typename std::enable_if<!has_apply<PreconditionerType>, U>::type
1162 finalize_after_convergence(
const unsigned int iteration_index)
1164 if (iteration_index % 2 == 1 || iteration_index == 2)
1165 this->x.add(this->alpha, this->p);
1168 const unsigned int end_range = this->x.locally_owned_size();
1170 Number * x = this->x.begin();
1171 Number * r = this->r.begin();
1172 Number * p = this->p.begin();
1173 const Number alpha_plus_previous_alpha_over_beta =
1174 this->alpha + this->previous_alpha / this->previous_beta;
1175 const Number previous_alpha_over_beta =
1176 this->previous_alpha / this->previous_beta;
1178 constexpr unsigned int grain_size = 128;
1179 std::array<Number, grain_size> prec_r;
1180 for (
unsigned int j = 0; j < end_range; j += grain_size)
1182 const unsigned int length =
std::min(grain_size, end_range - j);
1183 this->preconditioner.apply_to_subrange(j,
1188 for (
unsigned int i = 0; i < length; ++i)
1189 x[i] += alpha_plus_previous_alpha_over_beta * p[i] -
1190 previous_alpha_over_beta * prec_r[i];
1204template <
typename VectorType>
1205template <
typename MatrixType,
typename PreconditionerType>
1209 const VectorType & b,
1210 const PreconditionerType &preconditioner)
1212 using number =
typename VectorType::value_type;
1219 const bool do_eigenvalues =
1220 !condition_number_signal.empty() || !all_condition_numbers_signal.empty() ||
1221 !eigenvalues_signal.empty() || !all_eigenvalues_signal.empty();
1224 std::vector<typename VectorType::value_type>
diagonal;
1225 std::vector<typename VectorType::value_type> offdiagonal;
1227 typename VectorType::value_type eigen_beta_alpha = 0;
1231 internal::SolverCG::
1232 IterationWorker<VectorType, MatrixType, PreconditionerType>
1234 A, preconditioner, determine_beta_by_flexible_formula, this->memory, x);
1238 solver_state = this->iteration_status(0, worker.residual_norm, x);
1246 worker.do_iteration(it);
1248 print_vectors(it, x, worker.r, worker.p);
1252 this->coefficients_signal(worker.previous_alpha, worker.beta);
1257 diagonal.push_back(number(1.) / worker.previous_alpha +
1259 eigen_beta_alpha = worker.beta / worker.previous_alpha;
1260 offdiagonal.push_back(
std::sqrt(worker.beta) /
1261 worker.previous_alpha);
1263 compute_eigs_and_cond(diagonal,
1265 all_eigenvalues_signal,
1266 all_condition_numbers_signal);
1269 solver_state = this->iteration_status(it, worker.residual_norm, x);
1272 worker.finalize_after_convergence(it);
1274 compute_eigs_and_cond(diagonal,
1277 condition_number_signal);
1285template <
typename VectorType>
1286boost::signals2::connection
1288 const std::function<
void(
typename VectorType::value_type,
1289 typename VectorType::value_type)> &slot)
1291 return coefficients_signal.connect(slot);
1296template <
typename VectorType>
1297boost::signals2::connection
1299 const std::function<
void(
double)> &slot,
1300 const bool every_iteration)
1302 if (every_iteration)
1304 return all_condition_numbers_signal.connect(slot);
1308 return condition_number_signal.connect(slot);
1314template <
typename VectorType>
1315boost::signals2::connection
1317 const std::function<
void(
const std::vector<double> &)> &slot,
1318 const bool every_iteration)
1320 if (every_iteration)
1322 return all_eigenvalues_signal.connect(slot);
1326 return eigenvalues_signal.connect(slot);
1332template <
typename VectorType>
1335 const AdditionalData &)
1343template <
typename VectorType>
1345 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(Number *ptr) const
void load(const Number *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)