dune-composites (unstable)

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