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 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154
|
// Copyright (c) 2020 Chris Richardson
// FEniCS Project
// SPDX-License-Identifier: MIT
#include "regge.h"
#include "lattice.h"
#include "polyset.h"
#include <iostream>
using namespace basix;
namespace
{
//-----------------------------------------------------------------------------
Eigen::MatrixXd create_regge_space(cell::type celltype, int degree)
{
if (celltype != cell::type::triangle and celltype != cell::type::tetrahedron)
throw std::runtime_error("Unsupported celltype");
const int tdim = cell::topological_dimension(celltype);
const int nc = tdim * (tdim + 1) / 2;
const int basis_size = polyset::dim(celltype, degree);
const int ndofs = basis_size * nc;
const int psize = basis_size * tdim * tdim;
Eigen::ArrayXXd wcoeffs = Eigen::ArrayXXd::Zero(ndofs, psize);
int s = basis_size;
for (int i = 0; i < tdim; ++i)
{
for (int j = 0; j < tdim; ++j)
{
int xoff = i + tdim * j;
int yoff = i + j;
if (tdim == 3 and i > 0 and j > 0)
++yoff;
wcoeffs.block(yoff * s, xoff * s, s, s) = Eigen::MatrixXd::Identity(s, s);
}
}
return wcoeffs;
}
//-----------------------------------------------------------------------------
Eigen::MatrixXd create_regge_dual(cell::type celltype, int degree)
{
const int tdim = cell::topological_dimension(celltype);
const int basis_size = polyset::dim(celltype, degree);
const int ndofs = basis_size * (tdim + 1) * tdim / 2;
const int space_size = basis_size * tdim * tdim;
Eigen::ArrayXXd dualmat(ndofs, space_size);
std::vector<std::vector<std::vector<int>>> topology
= cell::topology(celltype);
const Eigen::ArrayXXd geometry = cell::geometry(celltype);
// dof counter
int dof = 0;
for (std::size_t dim = 1; dim < topology.size(); ++dim)
{
for (std::size_t i = 0; i < topology[dim].size(); ++i)
{
const Eigen::ArrayXXd entity_geom
= cell::sub_entity_geometry(celltype, dim, i);
Eigen::ArrayXd point = entity_geom.row(0);
cell::type ct = cell::sub_entity_type(celltype, dim, i);
Eigen::ArrayXXd lattice
= lattice::create(ct, degree + 2, lattice::type::equispaced, false);
Eigen::ArrayXXd pts(lattice.rows(), entity_geom.cols());
for (int j = 0; j < lattice.rows(); ++j)
{
pts.row(j) = entity_geom.row(0);
for (int k = 0; k < entity_geom.rows() - 1; ++k)
{
pts.row(j)
+= (entity_geom.row(k + 1) - entity_geom.row(0)) * lattice(j, k);
}
}
Eigen::MatrixXd basis = polyset::tabulate(celltype, degree, 0, pts)[0];
// Store up outer(t, t) for all tangents
std::vector<int>& vert_ids = topology[dim][i];
int ntangents = dim * (dim + 1) / 2;
std::vector<Eigen::MatrixXd> vvt(ntangents);
int c = 0;
for (std::size_t s = 0; s < dim; ++s)
{
for (std::size_t d = s + 1; d < dim + 1; ++d)
{
const Eigen::VectorXd edge_t
= geometry.row(vert_ids[d]) - geometry.row(vert_ids[s]);
// outer product v.v^T
vvt[c++] = edge_t * edge_t.transpose();
}
}
for (int k = 0; k < pts.rows(); ++k)
{
for (int j = 0; j < ntangents; ++j)
{
Eigen::Map<Eigen::VectorXd> vvt_flat(vvt[j].data(),
vvt[j].rows() * vvt[j].cols());
// outer product: outer(outer(t, t), basis)
const Eigen::MatrixXd vvt_b = vvt_flat * basis.row(k);
// Copy tensor values row by row into dualmat
for (int r = 0; r < vvt_b.rows(); ++r)
dualmat.block(dof, r * vvt_b.cols(), 1, vvt_b.cols())
= vvt_b.row(r);
++dof;
}
}
}
}
return dualmat;
}
//-----------------------------------------------------------------------------
} // namespace
//-----------------------------------------------------------------------------
FiniteElement basix::create_regge(cell::type celltype, int degree,
const std::string& name)
{
const int tdim = cell::topological_dimension(celltype);
Eigen::MatrixXd wcoeffs = create_regge_space(celltype, degree);
Eigen::MatrixXd dual = create_regge_dual(celltype, degree);
// TODO
const int ndofs = dual.rows();
int perm_count = tdim == 2 ? 3 : 14;
std::vector<Eigen::MatrixXd> base_permutations(
perm_count, Eigen::MatrixXd::Identity(ndofs, ndofs));
Eigen::MatrixXd coeffs = compute_expansion_coefficients(wcoeffs, dual);
// Regge has (d+1) dofs on each edge, 3d(d+1)/2 on each face
// and d(d-1)(d+1) on the interior in 3D
const std::vector<std::vector<std::vector<int>>> topology
= cell::topology(celltype);
std::vector<std::vector<int>> entity_dofs(topology.size());
entity_dofs[0].resize(topology[0].size(), 0);
entity_dofs[1].resize(topology[1].size(), degree + 1);
entity_dofs[2].resize(topology[2].size(), 3 * (degree + 1) * degree / 2);
if (tdim > 2)
entity_dofs[3] = {(degree + 1) * degree * (degree - 1)};
return FiniteElement(name, celltype, degree, {tdim, tdim}, coeffs,
entity_dofs, base_permutations, {}, {}, "double covariant piola");
}
//-----------------------------------------------------------------------------
|