This file is indexed.

/usr/include/dune/pdelab/constraints/hangingnode.hh is in libdune-pdelab-dev 2.0.0-1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
// -*- tab-width: 4; indent-tabs-mode: nil -*-
#ifndef DUNE_PDELAB_HANGINGNODECONSTRAINTS_HH
#define DUNE_PDELAB_HANGINGNODECONSTRAINTS_HH

#include <cstddef>

#include<dune/common/exceptions.hh>
#include<dune/geometry/referenceelements.hh>
#include<dune/geometry/type.hh>
#include<dune/pdelab/common/geometrywrapper.hh>
#include"conforming.hh"
#include"hangingnodemanager.hh"

namespace Dune {
  namespace PDELab {

    //! \addtogroup Constraints
    //! \ingroup FiniteElementMap
    //! \{

    class HangingNodesConstraintsAssemblers
    {
    public:
      class CubeGridQ1Assembler
      {
      public:
        template<typename IG, typename LFS, typename T, typename FlagVector>
        static void assembleConstraints(const FlagVector & nodeState_e, const FlagVector & nodeState_f,
                                        const bool & e_has_hangingnodes, const bool & f_has_hangingnodes,
                                        const LFS & lfs_e, const LFS & lfs_f,
                                        T& trafo_e, T& trafo_f,
                                        const IG& ig)
        {
          typedef IG Intersection;
          typedef typename Intersection::EntityPointer CellEntityPointer;
          typedef typename Intersection::Entity Cell;
          typedef typename Intersection::Geometry FaceGeometry;
          typedef typename FaceGeometry::ctype DT;
          typedef typename LFS::Traits::SizeType SizeType;

          typedef typename LFS::Traits::GridFunctionSpace::Traits::GridView::IndexSet IndexSet;
          const CellEntityPointer e = ig.inside();
          const CellEntityPointer f = ! ig.boundary() ? ig.outside() : ig.inside();

          const std::size_t dimension = Intersection::dimension;

          typedef Dune::ReferenceElement<DT,dimension> GRE;
          const GRE& refelement_e = Dune::ReferenceElements<DT,dimension>::general(e->type());
          const GRE& refelement_f = Dune::ReferenceElements<DT,dimension>::general(f->type());

          // If both entities have hangingnodes, then the face is
          // conforming and no constraints have to be applied.
          if(e_has_hangingnodes && f_has_hangingnodes)
            return;

          // Choose local function space etc for element with hanging nodes
          const LFS & lfs = e_has_hangingnodes ? lfs_e : lfs_f;
          const IndexSet& indexSet = lfs.gridFunctionSpace().gridView().indexSet();

          const Cell& cell = *(e_has_hangingnodes ? e : f);
          const int faceindex = e_has_hangingnodes ? ig.indexInInside() : ig.indexInOutside();
          const GRE & refelement = e_has_hangingnodes ? refelement_e : refelement_f;
          const FlagVector & nodeState = e_has_hangingnodes ? nodeState_e : nodeState_f;
          T & trafo = e_has_hangingnodes ? trafo_e : trafo_f;

          // A map mapping the local indices from the face to local
          // indices of the cell
          std::vector<int> m(refelement.size(faceindex,1,dimension));
          for (int j=0; j<refelement.size(faceindex,1,dimension); j++)
            m[j] = refelement.subEntity(faceindex,1,j,dimension);

          // A map mapping the local indices from the face to global gridview indices
          std::vector<std::size_t> global_vertex_idx(refelement.size(faceindex,1,dimension));
          for (int j=0; j<refelement.size(faceindex,1,dimension); ++j)
            global_vertex_idx[j] = indexSet.subIndex(cell,refelement.subEntity(faceindex,1,j,dimension),dimension);

          // Create a DOFIndex that we will use to manually craft the correct dof indices for the constraints trafo
          // We copy one of the indices from the LocalFunctionSpace; that way, we automatically get the correct
          // TreeIndex into the DOFIndex and only have to fiddle with the EntityIndex.
          typename LFS::Traits::DOFIndex dof_index(lfs.dofIndex(0));

          typedef typename LFS::Traits::GridFunctionSpace::Ordering::Traits::DOFIndexAccessor::GeometryIndex GeometryIndexAccessor;
          const GeometryType vertex_gt(0);

          // Find the corresponding entity in the reference element
          for (int j=0; j<refelement.size(faceindex,1,dimension); j++){

            // The contribution factors of the base function bound to entity j
            typename T::RowType contribution;

            if(dimension == 3){

              assert(nodeState.size() == 8);

              const SizeType i = 4*j;

              // Neigbor relations in local indices on a quadrilateral face:
              // {Node, Direct Neighbor, Direct Neighbor, Diagonal Neighbor, Node, ...}
              const unsigned int fi[16] = {0,1,2,3, 1,0,3,2, 2,0,3,1, 3,1,2,0};

              // Only hanging nodes have contribution to other nodes
              if(nodeState[m[j]].isHanging()){

                // If both neighbors are hanging nodes, then this node
                // is diagonal to the target of the contribution
                if(nodeState[m[fi[i+1]]].isHanging() && nodeState[m[fi[i+2]]].isHanging())
                  {
                    GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                 vertex_gt,
                                                 global_vertex_idx[fi[i+3]]);

                    contribution[dof_index] = 0.25;

                    GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                 vertex_gt,
                                                 global_vertex_idx[j]);

                    trafo[dof_index] = contribution;
                  }
                // Direct neigbor
                else if(!nodeState[m[fi[i+1]]].isHanging())
                  {
                    GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                 vertex_gt,
                                                 global_vertex_idx[fi[i+1]]);

                    contribution[dof_index] = 0.5;

                    GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                 vertex_gt,
                                                 global_vertex_idx[j]);

                    trafo[dof_index] = contribution;
                  }
                // Direct neigbor
                else if(!nodeState[m[fi[i+2]]].isHanging())
                  {
                    GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                 vertex_gt,
                                                 global_vertex_idx[fi[i+2]]);

                    contribution[dof_index] = 0.5;

                    GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                 vertex_gt,
                                                 global_vertex_idx[j]);

                    trafo[dof_index] = contribution;
                  }
              }

            } else if(dimension == 2){

              assert(nodeState.size() == 4);


              // Only hanging nodes have contribution to other nodes
              if(nodeState[m[j]].isHanging()){

                const SizeType n_j = 1-j;

                assert( !nodeState[m[n_j]].isHanging() );

                GeometryIndexAccessor::store(dof_index.entityIndex(),
                                             vertex_gt,
                                             global_vertex_idx[n_j]);

                contribution[dof_index] = 0.5;

                GeometryIndexAccessor::store(dof_index.entityIndex(),
                                             vertex_gt,
                                             global_vertex_idx[j]);

                trafo[dof_index] = contribution;
              }

            } // end if(dimension==3)

          } // j

        } // end of static void assembleConstraints

      }; // end of class CubeGridQ1Assembler


      class SimplexGridP1Assembler
      {
      public:
        template<typename IG,
                 typename LFS,
                 typename T,
                 typename FlagVector>
        static void assembleConstraints( const FlagVector & nodeState_e,
                                         const FlagVector & nodeState_f,
                                         const bool & e_has_hangingnodes,
                                         const bool & f_has_hangingnodes,
                                         const LFS & lfs_e, const LFS & lfs_f,
                                         T& trafo_e, T& trafo_f,
                                         const IG& ig)
        {
          typedef IG Intersection;
          typedef typename Intersection::EntityPointer CellEntityPointer;
          typedef typename Intersection::Entity Cell;
          typedef typename Intersection::Geometry FaceGeometry;
          typedef typename FaceGeometry::ctype DT;
          typedef typename LFS::Traits::SizeType SizeType;
          typedef typename LFS::Traits::GridFunctionSpace::Traits::GridView::IndexSet IndexSet;

          const CellEntityPointer e = ig.inside();
          const CellEntityPointer f = ! ig.boundary() ? ig.outside() : ig.inside();

          const std::size_t dimension = Intersection::dimension;

          typedef Dune::ReferenceElement<DT,dimension> GRE;
          const GRE& refelement_e = Dune::ReferenceElements<DT,dimension>::general(e->type());
          const GRE& refelement_f = Dune::ReferenceElements<DT,dimension>::general(f->type());

          // If both entities have hangingnodes, then the face is
          // conforming and no constraints have to be applied.
          if(e_has_hangingnodes && f_has_hangingnodes)
            return;

          // Choose local function space etc for element with hanging nodes
          const LFS & lfs = e_has_hangingnodes ? lfs_e : lfs_f;
          const IndexSet& indexSet = lfs.gridFunctionSpace().gridView().indexSet();

          const Cell& cell = *(e_has_hangingnodes ? e : f);
          const int faceindex = e_has_hangingnodes ? ig.indexInInside() : ig.indexInOutside();
          const GRE & refelement = e_has_hangingnodes ? refelement_e : refelement_f;
          const FlagVector & nodeState = e_has_hangingnodes ? nodeState_e : nodeState_f;
          T & trafo = e_has_hangingnodes ? trafo_e : trafo_f;

          // A map mapping the local indices from the face to local
          // indices of the cell
          std::vector<int> m(refelement.size(faceindex,1,dimension));
          for (int j=0; j<refelement.size(faceindex,1,dimension); j++)
            m[j] = refelement.subEntity(faceindex,1,j,dimension);

          // A map mapping the local indices from the face to global gridview indices
          std::vector<std::size_t> global_vertex_idx(refelement.size(faceindex,1,dimension));
          for (int j=0; j<refelement.size(faceindex,1,dimension); ++j)
            global_vertex_idx[j] = indexSet.subIndex(cell,refelement.subEntity(faceindex,1,j,dimension),dimension);

          // Create a DOFIndex that we will use to manually craft the correct dof indices for the constraints trafo
          // We copy one of the indices from the LocalFunctionSpace; that way, we automatically get the correct
          // TreeIndex into the DOFIndex and only have to fiddle with the EntityIndex.
          typename LFS::Traits::DOFIndex dof_index(lfs.dofIndex(0));

          typedef typename LFS::Traits::GridFunctionSpace::Ordering::Traits::DOFIndexAccessor::GeometryIndex GeometryIndexAccessor;
          const GeometryType vertex_gt(0);

          // Find the corresponding entity in the reference element
          for (int j=0; j<refelement.size(faceindex,1,dimension); j++){

            // The contribution factors of the base function bound to entity j
            typename T::RowType contribution;

            if(dimension == 3){

              assert(nodeState.size() == 4);
              // Only hanging nodes have contribution to other nodes
              if(nodeState[m[j]].isHanging()){
                for( int k=1; k<=2; ++k ){

                  const SizeType n_j = (j+k)%3;

                  if( !nodeState[m[n_j]].isHanging() )
                    {
                      GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                   vertex_gt,
                                                   global_vertex_idx[n_j]);

                      contribution[dof_index] = 0.5;

                      GeometryIndexAccessor::store(dof_index.entityIndex(),
                                                   vertex_gt,
                                                   global_vertex_idx[j]);

                      trafo[dof_index] = contribution;
                    }
                }
              }
            } else if(dimension == 2){

              assert(nodeState.size() == 3);
              // Only hanging nodes have contribution to other nodes
              if(nodeState[m[j]].isHanging()){
                const SizeType n_j = 1-j;
                assert( !nodeState[m[n_j]].isHanging() );
                // If both neighbors are hanging nodes, then this node
                // is diagonal to the target of the contribution
                GeometryIndexAccessor::store(dof_index.entityIndex(),
                                             vertex_gt,
                                             global_vertex_idx[n_j]);

                contribution[dof_index] = 0.5;

                GeometryIndexAccessor::store(dof_index.entityIndex(),
                                             vertex_gt,
                                             global_vertex_idx[j]);

                trafo[dof_index] = contribution;
              }


            } // end if(dimension==3)

          } // j

        } // end of static void assembleConstraints

      }; // end of class SimplexGridP1Assembler

    }; // end of class HangingNodesConstraintsAssemblers


    //! Hanging Node constraints construction
    // works in any dimension and on all element types
    template <class Grid, class HangingNodesConstraintsAssemblerType, class BoundaryFunction>
    class HangingNodesDirichletConstraints : public ConformingDirichletConstraints
    {
    private:
      typedef Dune::PDELab::HangingNodeManager<Grid,BoundaryFunction> HangingNodeManager;
      HangingNodeManager manager;

    public:
      enum { doBoundary = true };
      enum { doSkeleton = true };
      enum { doVolume = false };
      enum { dimension = Grid::dimension };

      HangingNodesDirichletConstraints( Grid & grid,
                                        bool adaptToIsolatedHangingNodes,
                                        const BoundaryFunction & _boundaryFunction )
        : manager(grid, _boundaryFunction)
      {
        if(adaptToIsolatedHangingNodes)
          manager.adaptToIsolatedHangingNodes();
      }

      void update( Grid & grid ){
        manager.analyzeView();
        manager.adaptToIsolatedHangingNodes();
      }


      //! skeleton constraints
      /**
       * \tparam I  intersection geometry
       * \tparam LFS local function space
       * \tparam T   TransformationType
       */
      template<typename IG, typename LFS, typename T>
      void skeleton (const IG& ig,
                     const LFS& lfs_e, const LFS& lfs_f,
                     T& trafo_e, T& trafo_f) const
      {
        typedef IG Intersection;
        typedef typename Intersection::EntityPointer CellEntityPointer;
        typedef typename Intersection::Geometry FaceGeometry;
        typedef typename FaceGeometry::ctype DT;

        const CellEntityPointer e = ig.inside();
        const CellEntityPointer f = ig.outside();

        const Dune::ReferenceElement<DT,dimension>& refelem_e
          = Dune::ReferenceElements<DT,dimension>::general(e->type());
        const Dune::ReferenceElement<DT,dimension>& refelem_f
          = Dune::ReferenceElements<DT,dimension>::general(f->type());

        // the return values of the hanging node manager
        typedef typename std::vector<typename HangingNodeManager::NodeState> FlagVector;
        const FlagVector isHangingNode_e(manager.hangingNodes(*e));
        const FlagVector isHangingNode_f(manager.hangingNodes(*f));

        // just to make sure that the hanging node manager is doing
        // what is expected of him
        assert(std::size_t(refelem_e.size(dimension))
               == isHangingNode_e.size());
        assert(std::size_t(refelem_f.size(dimension))
               == isHangingNode_f.size());

        // the LOCAL indices of the faces in the reference element
        const int faceindex_e = ig.indexInInside();
        const int faceindex_f = ig.indexInOutside();

        bool e_has_hangingnodes = false;
        {
          for (int j=0; j<refelem_e.size(faceindex_e,1,dimension); j++){
            const int index = refelem_e.subEntity(faceindex_e,1,j,dimension);
            if(isHangingNode_e[index].isHanging())
              {
                e_has_hangingnodes = true;
                break;
              }
          }
        }
        bool f_has_hangingnodes = false;
        {
          for (int j=0; j<refelem_f.size(faceindex_f,1,dimension); j++){
            const int index = refelem_f.subEntity(faceindex_f,1,j,dimension);
            if(isHangingNode_f[index].isHanging())
              {
                f_has_hangingnodes = true;
                break;
              }
          }
        }

        if(! e_has_hangingnodes && ! f_has_hangingnodes)
          return;

        HangingNodesConstraintsAssemblerType::
          assembleConstraints(isHangingNode_e, isHangingNode_f,
                              e_has_hangingnodes, f_has_hangingnodes,
                              lfs_e,lfs_f,
                              trafo_e, trafo_f,
                              ig);
      } // skeleton

    }; // end of class HangingNodesDirichletConstraints
    //! \}

  }
}

#endif