dune-composites (2.5.1)

zembasis.hh
1
2#ifndef DUNE_GENEO_ZEMBASIS_HH
3#define DUNE_GENEO_ZEMBASIS_HH
4
5#include "subdomainbasis.hh"
6
7namespace Dune{
8 namespace Geneo{
9
10template<class GFS, class LFS, class X, int dim, const int d>
11class ZEMBasis : public SubdomainBasis<X>
12{
13
14public:
15 ZEMBasis(const GFS& gfs, LFS& lfs, X& part_unity) {
16 using Dune::PDELab::Backend::native;
17
18 this->local_basis.resize(6);
19
20 //First add displacements in each direction
21 std::vector<Dune::FieldVector<double,dim>> a(3);
22 a[0] = {1., 0. ,0.};
23 a[1] = {0., 1., 0.};
24 a[2] = {0., 0., 1.};
25 for (int k = 0; k < 3; k++){
26 this->local_basis[k] = std::make_shared<X>(part_unity);
27 for (auto iter = native(*this->local_basis[k]).begin(); iter != native(*this->local_basis[k]).end(); iter++){
28 for(int j = 0; j< dim; j++)
29 (*iter)[j] *= a[k][j];
30 }
31 }
32
33 //Then add rotations
34 for(int k = 0; k < 3; k++){
35 auto zv = X(gfs,0.0);
36 this->local_basis[k+3] = std::make_shared<X>(part_unity);
37
38 for (auto it = gfs.gridView().template begin<0>(); it != gfs.gridView().template end<0>(); ++it) {
39 lfs.bind(*it);
40
41 auto geo = it->geometry();
42 const auto gt = geo.type();
43 const auto& ref_el = Dune::ReferenceElements<double, d>::general(gt);
44
45 auto& coeffs = lfs.finiteElement().localCoefficients();
46 for (std::size_t i = 0; i < coeffs.size(); ++i) {
47
48 auto local_pos = ref_el.position (coeffs.localKey(i).subEntity(), coeffs.localKey(i).codim());
49 auto global_pos = geo.global(local_pos);
50
51 auto subindex = gfs.entitySet().indexSet().subIndex(*it, coeffs.localKey(i).subEntity(), coeffs.localKey(i).codim());
52
53 //ensure each subindex is only written once
54 if(!(native(zv)[subindex][0] > 0.0)){
55 auto val = global_pos;
56 double val_h = val[k];
57 val[k] = val[(k+1)%3];
58 val[(k+1)%2] = val_h;
59
60 for(int j = 0; j < dim; j++){
61 native(*this->local_basis[k+3])[subindex][j] *= (global_pos[j] - val[j]);
62 }
63 native(zv)[subindex][0] = 1;
64 }
65 }
66 }
67 }
68 }
69};
70
71}
72}
73
74#endif //DUNE_GENEO_ZEMBASIS_HH
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Sep 5, 22:35, 2025)