File: e-hhj.cpp

package info (click to toggle)
fenics-basix 0.10.0.post0-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 3,156 kB
  • sloc: cpp: 23,435; python: 10,829; makefile: 43; sh: 26
file content (309 lines) | stat: -rw-r--r-- 12,530 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
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
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
// Copyright (c) 2022 Matthew Scroggs
// FEniCS Project
// SPDX-License-Identifier:    MIT

#include "e-hhj.h"
#include "e-lagrange.h"
#include "element-families.h"
#include "maps.h"
#include "math.h"
#include "polyset.h"
#include "quadrature.h"
#include "sobolev-spaces.h"
#include <cmath>

using namespace basix;

//-----------------------------------------------------------------------------
template <std::floating_point T>
FiniteElement<T> basix::element::create_hhj(cell::type celltype, int degree,
                                            bool discontinuous)
{
  if (celltype != cell::type::triangle and celltype != cell::type::tetrahedron)
    throw std::runtime_error("Unsupported celltype");

  const std::size_t tdim = cell::topological_dimension(celltype);

  const int nc = tdim * (tdim + 1) / 2;
  const int basis_size
      = polyset::dim(celltype, polyset::type::standard, degree);

  impl::mdarray_t<T, 2> wcoeffs(basis_size * nc, basis_size * tdim * tdim);
  for (std::size_t i = 0; i < tdim; ++i)
  {
    for (std::size_t j = 0; j < tdim; ++j)
    {
      int xoff = i + tdim * j;
      int yoff = i + j;
      if (tdim == 3 and i > 0 and j > 0)
        ++yoff;

      const std::size_t s = basis_size;
      for (std::size_t k = 0; k < s; ++k)
        wcoeffs(yoff * s + k, xoff * s + k) = i == j ? 1.0 : std::sqrt(0.5);
    }
  }

  const std::vector<std::vector<std::vector<int>>> topology
      = cell::topology(celltype);
  const auto [gbuffer, gshape] = cell::geometry<T>(celltype);
  impl::mdspan_t<const T, 2> geometry(gbuffer.data(), gshape);

  std::array<std::vector<impl::mdarray_t<T, 2>>, 4> x;
  std::array<std::vector<impl::mdarray_t<T, 4>>, 4> M;

  const std::size_t facet_dim = tdim - 1;

  // No DOFs on entities of lower dimension than facets
  for (std::size_t d = 0; d < facet_dim; ++d)
  {
    for (std::size_t e = 0; e < topology[d].size(); ++e)
    {
      x[d].emplace_back(0, tdim);
      M[d].emplace_back(0, tdim * tdim, 0, 1);
    }
  }
  // DOFs on facets
  // These dofs are defined by l(F) = integral(n^t @ F @ n * p),
  // where n is normal to the facet, p is polynomial on the facet
  // of degree `degree`
  auto [_data, _shape] = cell::scaled_facet_normals<T>(celltype);
  impl::mdspan_t<const T, 2> normals(_data.data(), _shape);

  for (std::size_t e = 0; e < topology[facet_dim].size(); ++e)
  {
    cell::type ct = cell::sub_entity_type(celltype, facet_dim, e);

    const std::size_t ndofs = polyset::dim(ct, polyset::type::standard, degree);
    const auto [ptsbuffer, wts] = quadrature::make_quadrature<T>(
        quadrature::type::Default, ct, polyset::type::standard, 2 * degree);
    impl::mdspan_t<const T, 2> pts(ptsbuffer.data(), wts.size(), facet_dim);

    FiniteElement<T> moment_space = create_lagrange<T>(
        ct, degree, element::lagrange_variant::legendre, true);
    const auto [phib, phishape] = moment_space.tabulate(0, pts);
    impl::mdspan_t<const T, 4> moment_values(phib.data(), phishape);

    // Entity coordinates
    const auto [entity_x_buffer, eshape]
        = cell::sub_entity_geometry<T>(celltype, facet_dim, e);
    std::span<const T> x0(entity_x_buffer.data(), eshape[1]);
    impl::mdspan_t<const T, 2> entity_x(entity_x_buffer.data(), eshape);

    // Copy points
    auto& _x = x[facet_dim].emplace_back(pts.extent(0), tdim);
    for (std::size_t p = 0; p < pts.extent(0); ++p)
    {
      for (std::size_t k = 0; k < _x.extent(1); ++k)
        _x(p, k) = x0[k];
      for (std::size_t k0 = 0; k0 + 1 < entity_x.extent(0); ++k0)
        for (std::size_t k1 = 0; k1 < _x.extent(1); ++k1)
          _x(p, k1) += (entity_x(k0 + 1, k1) - x0[k1]) * pts(p, k0);
    }

    auto& _M = M[facet_dim].emplace_back(ndofs, tdim * tdim, pts.extent(0), 1);
    for (int n = 0; n < moment_space.dim(); ++n)
    {
      for (std::size_t q = 0; q < pts.extent(0); ++q)
      {
        for (std::size_t k0 = 0; k0 < tdim; ++k0)
        {
          for (std::size_t k1 = 0; k1 < tdim; ++k1)
          {
            _M(n, tdim * k0 + k1, q, 0) = normals(e, k0) * normals(e, k1)
                                          * wts[q] * moment_values(0, q, n, 0);
          }
        }
      }
    }
  }

  // Interior
  // These DOFs are defined as in section 5.2 of
  // https://doi.org/10.1007/s00211-017-0933-3 (Pechstein & Schoberl, 2017)
  if (tdim == 2)
  {
    if (degree == 0)
    {
      x[tdim].emplace_back(0, tdim);
      M[tdim].emplace_back(0, tdim * tdim, 0, 1);
    }
    else
    {
      const std::size_t ndofs
          = polyset::dim(celltype, polyset::type::standard, degree - 1);

      const auto [ptsbuffer, wts] = quadrature::make_quadrature<T>(
          quadrature::type::Default, celltype, polyset::type::standard,
          2 * degree - 1);
      impl::mdspan_t<const T, 2> pts(ptsbuffer.data(), wts.size(), tdim);

      // Copy points
      auto& _x = x[tdim].emplace_back(pts.extent(0), tdim);
      for (std::size_t p = 0; p < pts.extent(0); ++p)
      {
        for (std::size_t k = 0; k < pts.extent(1); ++k)
          _x(p, k) += pts(p, k);
      }

      auto& _M = M[tdim].emplace_back(ndofs * 3, tdim * tdim, pts.extent(0), 1);

      FiniteElement<T> moment_space = create_lagrange<T>(
          celltype, degree - 1, element::lagrange_variant::legendre, true);
      const auto [phib, phishape] = moment_space.tabulate(0, pts);
      impl::mdspan_t<const T, 4> moment_values(phib.data(), phishape);

      for (int n = 0; n < moment_space.dim(); ++n)
      {
        for (std::size_t q = 0; q < pts.extent(0); ++q)
        {
          // [0, 1], [1, 0]
          _M(n * 3, 1, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 3, 2, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          // [-2, 1], [1, 0]
          _M(n * 3 + 1, 0, q, 0) = -2.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 3 + 1, 1, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 3 + 1, 2, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          // [0, 1], [1, -2]
          _M(n * 3 + 2, 1, q, 0) = -1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 3 + 2, 2, q, 0) = -1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 3 + 2, 3, q, 0) = 2.0 * wts[q] * moment_values(0, q, n, 0);
        }
      }
    }
  }
  else
  {
    const std::size_t ndofs0
        = degree == 0
              ? 0
              : polyset::dim(celltype, polyset::type::standard, degree - 1);
    const std::size_t ndofs1
        = polyset::dim(celltype, polyset::type::standard, degree);

    const auto [ptsbuffer, wts]
        = quadrature::make_quadrature<T>(quadrature::type::Default, celltype,
                                         polyset::type::standard, 2 * degree);
    impl::mdspan_t<const T, 2> pts(ptsbuffer.data(), wts.size(), tdim);

    // Copy points
    auto& _x = x[tdim].emplace_back(pts.extent(0), tdim);
    for (std::size_t p = 0; p < pts.extent(0); ++p)
    {
      for (std::size_t k = 0; k < pts.extent(1); ++k)
        _x(p, k) += pts(p, k);
    }

    auto& _M = M[tdim].emplace_back(ndofs0 * 4 + ndofs1 * 2, tdim * tdim,
                                    pts.extent(0), 1);

    if (degree > 0)
    {
      FiniteElement<T> moment_space = create_lagrange<T>(
          celltype, degree - 1, element::lagrange_variant::legendre, true);
      const auto [phib, phishape] = moment_space.tabulate(0, pts);
      impl::mdspan_t<const T, 4> moment_values(phib.data(), phishape);

      for (int n = 0; n < moment_space.dim(); ++n)
      {
        for (std::size_t q = 0; q < pts.extent(0); ++q)
        {
          // [0, 1, 1], [1, 0, 1], [1, 1, 0]
          _M(n * 4, 1, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4, 2, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4, 3, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4, 5, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4, 6, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4, 7, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);

          // [-6, 1, 1], [1, 0, 1], [1, 1, 0]
          _M(n * 4 + 1, 0, q, 0) = -6.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 1, 1, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 1, 2, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 1, 3, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 1, 5, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 1, 6, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 1, 7, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);

          // [0, 1, 1], [1, -6, 1], [1, 1, 0]
          _M(n * 4 + 2, 1, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 2, 2, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 2, 3, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 2, 4, q, 0) = -6.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 2, 5, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 2, 6, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 2, 7, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);

          // [0, 1, 1], [1, 0, 1], [1, 1, -6]
          _M(n * 4 + 3, 1, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 3, 2, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 3, 3, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 3, 5, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 3, 6, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 3, 7, q, 0) = 1.0 * wts[q] * moment_values(0, q, n, 0);
          _M(n * 4 + 3, 8, q, 0) = -6.0 * wts[q] * moment_values(0, q, n, 0);
        }
      }
    }
    FiniteElement<T> moment_space = create_lagrange<T>(
        celltype, degree, element::lagrange_variant::legendre, true);
    const auto [phib, phishape] = moment_space.tabulate(0, pts);
    impl::mdspan_t<const T, 4> moment_values(phib.data(), phishape);

    for (int n = 0; n < moment_space.dim(); ++n)
    {
      for (std::size_t q = 0; q < pts.extent(0); ++q)
      {
        // [0, 0, -1], [0, 0, 1], [-1, 1, 0]
        _M(ndofs0 * 4 + n * 2, 2, q, 0)
            = -1.0 * wts[q] * moment_values(0, q, n, 0);
        _M(ndofs0 * 4 + n * 2, 5, q, 0)
            = 1.0 * wts[q] * moment_values(0, q, n, 0);
        _M(ndofs0 * 4 + n * 2, 6, q, 0)
            = -1.0 * wts[q] * moment_values(0, q, n, 0);
        _M(ndofs0 * 4 + n * 2, 7, q, 0)
            = 1.0 * wts[q] * moment_values(0, q, n, 0);
        // [0, -1, 0], [-1, 0, 1], [0, 1, 0]
        _M(ndofs0 * 4 + n * 2 + 1, 1, q, 0)
            = -1.0 * wts[q] * moment_values(0, q, n, 0);
        _M(ndofs0 * 4 + n * 2 + 1, 3, q, 0)
            = -1.0 * wts[q] * moment_values(0, q, n, 0);
        _M(ndofs0 * 4 + n * 2 + 1, 5, q, 0)
            = 1.0 * wts[q] * moment_values(0, q, n, 0);
        _M(ndofs0 * 4 + n * 2 + 1, 7, q, 0)
            = 1.0 * wts[q] * moment_values(0, q, n, 0);
      }
    }
  }

  std::array<std::vector<impl::mdspan_t<const T, 2>>, 4> xview
      = impl::to_mdspan(x);
  std::array<std::vector<impl::mdspan_t<const T, 4>>, 4> Mview
      = impl::to_mdspan(M);

  std::array<std::vector<std::vector<T>>, 4> xbuffer;
  std::array<std::vector<std::vector<T>>, 4> Mbuffer;
  if (discontinuous)
  {
    std::array<std::vector<std::array<std::size_t, 2>>, 4> xshape;
    std::array<std::vector<std::array<std::size_t, 4>>, 4> Mshape;
    std::tie(xbuffer, xshape, Mbuffer, Mshape)
        = element::make_discontinuous(xview, Mview, tdim, tdim * tdim);
    xview = impl::to_mdspan(xbuffer, xshape);
    Mview = impl::to_mdspan(Mbuffer, Mshape);
  }

  sobolev::space space
      = discontinuous ? sobolev::space::L2 : sobolev::space::HDivDiv;

  return FiniteElement<T>(
      element::family::HHJ, celltype, polyset::type::standard, degree,
      {tdim, tdim}, impl::mdspan_t<T, 2>(wcoeffs.data(), wcoeffs.extents()),
      xview, Mview, 0, maps::type::doubleContravariantPiola, space,
      discontinuous, -1, degree, element::lagrange_variant::unset,
      element::dpc_variant::unset);
}
//-----------------------------------------------------------------------------
template FiniteElement<float> element::create_hhj(cell::type, int, bool);
template FiniteElement<double> element::create_hhj(cell::type, int, bool);
//-----------------------------------------------------------------------------