1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73
|
// Copyright (c) 2020 Chris Richardson & Matthew Scroggs
// FEniCS Project
// SPDX-License-Identifier: MIT
#pragma once
#include <Eigen/Dense>
namespace basix
{
/// Functions to help with the creation of DOF permutation and direction
/// correction.
namespace dofperms
{
/// Reflect the DOFs on an interval
/// @param degree The number of DOFs on the interval
/// @return A reordering of the numbers 0 to degree-1 representing the
/// permutation
Eigen::ArrayXi interval_reflection(int degree);
/// Reflect the DOFs on a triangle
/// @param degree The number of DOFs along one side of the triangle
/// @return A reordering of the numbers 0 to (degree)*(degree+1)/2-1
/// representing the permutation
Eigen::ArrayXi triangle_reflection(int degree);
/// Rotate the DOFs on a triangle
/// @param degree The number of DOFs along one side of the triangle
/// @return A reordering of the numbers 0 to (degree)*(degree+1)/2-1
/// representing the permutation
Eigen::ArrayXi triangle_rotation(int degree);
/// Reflect the DOFs on a quadrilateral
/// @param degree The number of DOFs along one side of the quadrilateral
/// @return A reordering of the numbers 0 to degree*degree-1 representing the
/// permutation
Eigen::ArrayXi quadrilateral_reflection(int degree);
/// Rotate the DOFs on a quadrilateral
/// @param degree The number of DOFs along one side of the quadrilateral
/// @return A reordering of the numbers 0 to degree*degree-1 representing the
/// permutation
Eigen::ArrayXi quadrilateral_rotation(int degree);
//-----------------------------------------------------------------------------
/// Generate a matrix to correct the direction of tangent vector-values DOFs on
/// an interval when that interval is reflected
/// @param degree The number of DOFs on the interval
/// @return A matrix representing the effect of reversing the edge on the DOF
/// values
Eigen::ArrayXXd interval_reflection_tangent_directions(int degree);
/// Generate a matrix to correct the direction of tangent vector-values DOFs on
/// a triangle when that triangle is reflected
/// @param degree The number of DOFs along one side of the triangle
/// @return A matrix representing the effect of reflecting the triangle edge on
/// the DOF values
Eigen::ArrayXXd triangle_reflection_tangent_directions(int degree);
/// Generate a matrix to correct the direction of tangent vector-values DOFs on
/// a triangle when that triangle is rotated
/// @param degree The number of DOFs along one side of the triangle
/// @return A matrix representing the effect of rotating the triangle edge on
/// the DOF values
Eigen::ArrayXXd triangle_rotation_tangent_directions(int degree);
// TODO: quad tangent directions
}; // namespace dofperms
} // namespace basix
|