File: dof-permutations.cpp

package info (click to toggle)
basix 0.0.1~git20210122.4f10ef2-2
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 696 kB
  • sloc: cpp: 3,987; python: 1,918; makefile: 33
file content (114 lines) | stat: -rw-r--r-- 3,029 bytes parent folder | download
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
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
// Copyright (c) 2020 Chris Richardson & Matthew Scroggs
// FEniCS Project
// SPDX-License-Identifier:    MIT

#include "dof-permutations.h"

using namespace basix;

//-----------------------------------------------------------------------------
Eigen::ArrayXi dofperms::interval_reflection(int degree)
{
  Eigen::ArrayXi perm(degree);
  for (int i = 0; i < degree; ++i)
    perm(i) = degree - 1 - i;
  return perm;
}
//-----------------------------------------------------------------------------
Eigen::ArrayXi dofperms::triangle_reflection(int degree)
{
  const int n = degree * (degree + 1) / 2;
  Eigen::ArrayXi perm(n);
  int p = 0;

  for (int st = 0; st < degree; ++st)
  {
    int dof = st;
    for (int add = degree; add > st; --add)
    {
      perm(p++) = dof;
      dof += add;
    }
  }

  return perm;
}
//-----------------------------------------------------------------------------
Eigen::ArrayXi dofperms::triangle_rotation(int degree)
{
  const int n = degree * (degree + 1) / 2;
  Eigen::ArrayXi perm(n);
  int p = 0;
  int st = n - 1;
  for (int i = 1; i <= degree; ++i)
  {
    int dof = st;
    for (int sub = i; sub <= degree; ++sub)
    {
      perm(p++) = dof;
      dof -= sub + 1;
    }
    st -= i;
  }

  return perm;
}
//-----------------------------------------------------------------------------
Eigen::ArrayXi dofperms::quadrilateral_reflection(int degree)
{
  const int n = degree * degree;
  Eigen::ArrayXi perm(n);
  int p = 0;

  for (int st = 0; st < degree; ++st)
    for (int i = 0; i < degree; ++i)
      perm(p++) = st + i * degree;

  return perm;
}
//-----------------------------------------------------------------------------
Eigen::ArrayXi dofperms::quadrilateral_rotation(int degree)
{
  const int n = degree * degree;
  Eigen::ArrayXi perm(n);
  int p = 0;

  for (int st = degree - 1; st >= 0; --st)
    for (int i = 0; i < degree; ++i)
      perm(st + degree * i) = p++;

  return perm;
}
//-----------------------------------------------------------------------------
Eigen::ArrayXXd dofperms::interval_reflection_tangent_directions(int degree)
{
  return -Eigen::MatrixXd::Identity(degree, degree);
}
//-----------------------------------------------------------------------------
Eigen::ArrayXXd dofperms::triangle_reflection_tangent_directions(int degree)
{
  Eigen::ArrayXXd dirs
      = Eigen::ArrayXXd::Zero(degree * (degree + 1), degree * (degree + 1));
  for (int i = 0; i < degree * (degree + 1); i += 2)
  {
    dirs(i, i + 1) = 1;
    dirs(i + 1, i) = 1;
  }

  return dirs;
}
//-----------------------------------------------------------------------------
Eigen::ArrayXXd dofperms::triangle_rotation_tangent_directions(int degree)
{
  Eigen::ArrayXXd dirs
      = Eigen::ArrayXXd::Zero(degree * (degree + 1), degree * (degree + 1));
  for (int i = 0; i < degree * (degree + 1); i += 2)
  {
    dirs(i, i + 1) = -1;
    dirs(i + 1, i) = 1;
    dirs(i, i) = -1;
  }

  return dirs;
}
//-----------------------------------------------------------------------------