Reference documentation for deal.II version 9.3.3
|
#include <deal.II/base/tensor.h>
#include <deal.II/lac/exceptions.h>
#include <deal.II/lac/lapack_templates.h>
Go to the source code of this file.
Functions | |
template<int dim, typename Number > | |
Tensor< 2, dim, Number > | project_onto_orthogonal_tensors (const Tensor< 2, dim, Number > &A) |
template Tensor< 2, 1, float > | project_onto_orthogonal_tensors (const Tensor< 2, 1, float > &) |
template Tensor< 2, 2, float > | project_onto_orthogonal_tensors (const Tensor< 2, 2, float > &) |
template Tensor< 2, 3, float > | project_onto_orthogonal_tensors (const Tensor< 2, 3, float > &) |
template Tensor< 2, 1, double > | project_onto_orthogonal_tensors (const Tensor< 2, 1, double > &) |
template Tensor< 2, 2, double > | project_onto_orthogonal_tensors (const Tensor< 2, 2, double > &) |
template Tensor< 2, 3, double > | project_onto_orthogonal_tensors (const Tensor< 2, 3, double > &) |
Tensor< 2, dim, Number > project_onto_orthogonal_tensors | ( | const Tensor< 2, dim, Number > & | A | ) |
Return the nearest orthogonal matrix \(\hat {\mathbf A}=\mathbf U \mathbf{V}^T\) by combining the products of the singular value decomposition (SVD) \({\mathbf A}=\mathbf U \mathbf S \mathbf V^T\) for a given input \({\mathbf A}\), effectively replacing \(\mathbf S\) with the identity matrix.
This is a (nonlinear) projection operation since when applied twice, we have \(\hat{\hat{\mathbf A}}=\hat{\mathbf A}\) as is easy to see. (That is because the SVD of \(\hat {\mathbf A}\) is simply \(\mathbf U \mathbf I \mathbf{V}^T\).) Furthermore, \(\hat {\mathbf A}\) is really an orthogonal matrix because orthogonal matrices have to satisfy \({\hat {\mathbf A}}^T \hat {\mathbf A}={\mathbf I}\), which here implies that
\begin{align*} {\hat {\mathbf A}}^T \hat {\mathbf A} &= \left(\mathbf U \mathbf{V}^T\right)^T\left(\mathbf U \mathbf{V}^T\right) \\ &= \mathbf V \mathbf{U}^T \mathbf U \mathbf{V}^T \\ &= \mathbf V \left(\mathbf{U}^T \mathbf U\right) \mathbf{V}^T \\ &= \mathbf V \mathbf I \mathbf{V}^T \\ &= \mathbf V \mathbf{V}^T \\ &= \mathbf I \end{align*}
due to the fact that the \(\mathbf U\) and \(\mathbf V\) factors that come out of the SVD are themselves orthogonal matrices.
A | The tensor for which to find the closest orthogonal tensor. |
Number | The type used to store the entries of the tensor. Must be either float or double . |
A
must not be singular. This is because, conceptually, the problem to be solved here is trying to find a matrix \(\hat{\mathbf A}\) that minimizes some kind of distance from \(\mathbf A\) while satisfying the quadratic constraint \({\hat {\mathbf A}}^T \hat {\mathbf A}={\mathbf I}\). This is not so dissimilar to the kind of problem where one wants to find a vector \(\hat{\mathbf x}\in{\mathbb R}^n\) that minimizes the quadratic objective function \(\|\hat {\mathbf x} - \mathbf x\|^2\) for a given \(\mathbf x\) subject to the constraint \(\|\mathbf x\|^2=1\) – in other words, we are seeking the point \(\hat{\mathbf x}\) on the unit sphere that is closest to \(\mathbf x\). This problem has a solution for all \(\mathbf x\) except if \(\mathbf x=0\). The corresponding condition for the problem we are considering here is that \(\mathbf A\) must not have a zero eigenvalue.