This file is indexed.

/usr/include/dune/grid/test/checkindexset.hh is in libdune-grid-dev 2.5.1-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
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_GRID_TEST_CHECKINDEXSET_HH
#define DUNE_GRID_TEST_CHECKINDEXSET_HH

#include <algorithm>
#include <iostream>
#include <limits>
#include <map>
#include <set>
#include <vector>

#include <dune/common/fvector.hh>
#include <dune/common/stdstreams.hh>
#include <dune/geometry/referenceelements.hh>
#include <dune/grid/common/capabilities.hh>
#include <dune/grid/common/exceptions.hh>

/** @file
   @author Robert Kloefkorn
   @brief Provides a check of the grids index set.
 */

template< class Grid >
struct EnableLevelIntersectionIteratorCheck
{
  static const bool v = true;
};

template< class Grid >
struct EnableLevelIntersectionIteratorCheck< const Grid >
{
  static const bool v = EnableLevelIntersectionIteratorCheck< Grid >::v;
};


namespace Dune
{

  // Check whether two FieldVectors are equal in the max norm
  template <typename ctype, int dim>
  bool compareVec(const FieldVector<ctype,dim> & vx1 , const FieldVector<ctype,dim> & vx2 )
  {
    const ctype eps = 1e5 * std::numeric_limits<ctype>::epsilon();
    return (vx1-vx2).infinity_norm() < eps;
  }

  /** \brief Check various features of codim-0 entities (elements)
      \param grid The grid we are testing
      \param en The grid element we are testing
      \param lset The corresponding level set
      \param sout The output stream that status (non-error) messages will go to
      \param setOfVerticesPerSubEntity A map that associates a vector of vertex indices to a subEntity
      \param subEntityPerSetOfVertices The reverse: a map that associates a subEntity to a vector of vertex indices
      \param vertexCoordsMap Map that associates vertex positions to indices (integers)
   */
  template <int codim, class GridType,
      class IndexSetType, class OutputStreamImp,
      class MapType1 , class MapType2 , class MapType3 >
  void checkSubEntity ( const GridType &,
                        const typename GridType::template Codim<0>::Entity &en , const IndexSetType & lset,
                        OutputStreamImp & sout , MapType1 & setOfVerticesPerSubEntity ,
                        MapType2 & subEntityPerSetOfVertices,
                        const MapType3 & vertexCoordsMap )
  {
    enum { dim = GridType::dimension };
    const int dimworld = GridType::dimensionworld;
    typedef typename GridType::ctype coordType;

    const GeometryType type = en.type();
    assert( type == en.geometry().type() );

    if( type.isNone() )
      return;

    const ReferenceElement< coordType, dim > &refElem
      = ReferenceElements< coordType, dim >::general( type );

    // check whether the element has the number of codim-subentities mandated by the reference element
    if(int(en.subEntities(codim)) != refElem.size(0,0,codim))
    {
      std::cerr << "entity index = " << lset.index(en)
                << ", type = " << type
                << std::endl
                << "codim = " << codim
                << std::endl
                << "subEntities(codim) = " << en.subEntities(codim)
                << std::endl
                << "refElem.size(codim) = " << refElem.size(0,0,codim)
                << std::endl;
      DUNE_THROW(GridError,
                 "wrong number of subEntities of codim " << codim);
    }

    for( int subEntity = 0; subEntity < refElem.size( 0, 0, codim ); ++subEntity )
    {

      // Number of vertices of the current subentity
      int numVertices = refElem.size( subEntity, codim, dim );

      // every entity must have at least one vertex
      assert( numVertices > 0 );

      // The local vertex numbers
      std::vector<int> local (numVertices);

      for( int j = 0; j < numVertices; ++j )
        local[ j ] = refElem.subEntity ( subEntity, codim, j, dim );

      sout << numVertices << " vertices on subEntity< codim = " << codim << " >" << std::endl;
      sout << "check subentity [" << local[ 0 ];
      for( int j = 1; j < numVertices; ++j )
        sout << ", " << local[ j ];
      sout << "]" << std::endl;

      // The global vertex numbers
      std::vector<int> global(numVertices);
      for( int j = 0; j < numVertices; ++j )
        global[ j ] = lset.subIndex( en, local[ j ], dim );

      sout << "Found global numbers of entity [ ";
      for( int j = 0; j < numVertices; ++j )
        sout << global[ j ] << " ";
      sout << "]" << std::endl;

      // Check whether we get the same index if we
      // -- ask for the subEntity itself and then its index
      // -- ask for the subIndex directly
      typedef typename GridType::template Codim< codim >::Entity SubE;
      SubE subE = en.template subEntity< codim >( subEntity );

      if( lset.subIndex( en, subEntity, codim ) != lset.index( subE) )
      {
        std::cerr << "Index for subEntity does not match." << std::endl;
        assert( lset.subIndex( en, subEntity, codim ) == lset.index( subE) );
      }

      // Make a unique identifier for the subEntity.  Since indices are unique only per GeometryType,
      // we need a (index,GeometryType)-pair.
      std::pair<int, GeometryType> globalSubEntity( lset.index( subE ), subE.type() );
      assert( globalSubEntity.first >= 0 );
      sout << "local subentity " << subEntity << " consider subentity with global key (" << globalSubEntity.first << ","
           << globalSubEntity.second << ") on en = " << lset.index(en) << std::endl;

      if( subE.type() != subE.geometry().type() )
      {
        std::cerr << "Geometry types for subEntity don't match." << std::endl;
        assert( subE.type() == subE.geometry().type() );
      }

      // assert that all sub entities have the same level
      assert( subE.level() == en.level() );

      // Loop over all vertices
      for( int j = 0; j < numVertices; ++j )
      {

        // get entity pointer to subEntity vertex
        typedef typename GridType::template Codim< dim >::Entity VertexE;
        VertexE vxE = en.template subEntity< dim >( local[ j ] );

        // Find the global coordinate of the vertex by its index
        if(vertexCoordsMap.find(global[j]) != vertexCoordsMap.end())
        {
          // Check whether index and coordinate match
          FieldVector<coordType,dimworld> vxcheck ( vertexCoordsMap.find(global[j])->second );
          FieldVector< coordType, dimworld > vx1 = vxE.geometry().corner( 0 );
          if( ! compareVec( vxcheck, vx1 ) )
          {
            std::cerr << "ERROR map global vertex [" << global[j] << "] vx " << vxcheck << " is not " << vx1 << "\n";
            assert( compareVec( vxcheck, vx1 ) );
          }
        }
        sout << "vx[" << global[j] << "] = "  <<  subE.geometry().corner( j ) << "\n";
      }
      sout << "sort vector of global vertex\n";

      // sort vector of global vertex number for storage in map
      // the smallest entry is the first entry
      std::sort( global.begin(), global.end() );

      // check whether vertex key is already stored in map
      if(subEntityPerSetOfVertices.find(global) == subEntityPerSetOfVertices.end())
      {
        subEntityPerSetOfVertices[global] = globalSubEntity;
      }
      else
      {
        assert( globalSubEntity == subEntityPerSetOfVertices[global] );
      }

      // check whether subEntity is already stored in map
      if(setOfVerticesPerSubEntity.find(globalSubEntity) == setOfVerticesPerSubEntity.end() )
      {
        setOfVerticesPerSubEntity[globalSubEntity] = global;
      }
      else
      {
        std::vector<int> globalcheck = setOfVerticesPerSubEntity[globalSubEntity];
        if(! (global == globalcheck ))
        {
          std::cerr << "For subEntity key (" << globalSubEntity.first << "," << globalSubEntity.second << ") \n";
          std::cerr << "Got ";
          for(int j=0 ; j<numVertices; j++ )
          {
            std::cerr << global[j] << " ";
          }
          std::cerr << "\n";
          std::cerr << "Found ";
          for(int j=0 ; j<numVertices; j++ )
          {
            std::cerr << globalcheck [j] << " ";
          }
          std::cerr << "\n";
          DUNE_THROW(Dune::GridError, "global != globalcheck");
        }
      }

    }   // end check sub entities
    sout << "end check sub entities\n";
  }


  // check some functionality of grid
  template< int codim, class Grid, class GridView, class OutputStream >
  void checkIndexSetForCodim ( const Grid &grid, const GridView &view,
                               OutputStream &sout, bool levelIndex )
  {
    enum { dim = Grid :: dimension };
    enum { dimworld = Grid :: dimensionworld };

    typedef typename Grid :: ctype coordType;

    //typedef typename GridView :: template Codim< 0 > :: Entity EntityCodim0Type;
    typedef typename GridView :: IndexSet IndexSetType;
    typedef typename GridView :: template Codim< codim > :: Iterator IteratorType;

    const IndexSetType &lset = view.indexSet();

    sout <<"\n\nStart consistency check of index set \n\n";

    // ////////////////////////////////////////////////////////////////
    //   Check whether geomTypes() returns correct result
    // ////////////////////////////////////////////////////////////////

    std :: set< GeometryType > geometryTypes;

    const IteratorType endit = view.template end< codim >();
    IteratorType it = view.template begin< codim >();

    if (it == endit) return;

    for (; it!=endit; ++it) {
      // while we're here: check whether the GeometryTypes returned by the entity
      // and the Geometry match
      assert(it->type()==it->type());
      geometryTypes.insert(it->type());
    }

    const typename IndexSetType::Types types = lset.types( codim );
    bool geomTypesError = false;

    // Check whether all entries in the official geometry types list are contained in our self-computed one
    for( auto it = types.begin(); it != types.end(); ++it )
      if( geometryTypes.find( *it ) == geometryTypes.end() )
        geomTypesError = true;


    // And vice versa
    for( std::set<GeometryType>::iterator it = geometryTypes.begin(); it != geometryTypes.end(); ++it )
    {
      if( std::find( types.begin(), types.end(), *it ) == types.end() )
        geomTypesError = true;
    }



    if (geomTypesError) {

      std::cerr << "There is a mismatch in the list of geometry types of codim " << codim << "." << std::endl;
      std::cerr << "Geometry types present in the grid are:" << std::endl;
      for (std::set<GeometryType>::iterator it = geometryTypes.begin(); it!=geometryTypes.end(); ++it)
        std::cerr << "  " << *it << std::endl;

      std::cerr << std::endl << "but the method types() returned:" << std::endl;
      for( auto it = types.begin(); it != types.end(); ++it )
        std::cerr << "  " << *it << std::endl;

      DUNE_THROW(GridError, "!");
    }

    //*****************************************************************
    // check size of index set
    int gridsize = 0;
    {
      typedef typename GridView :: template Codim< codim > :: Iterator IteratorType;

      int count = 0;
      const IteratorType endit = view.template end< codim >();
      for( IteratorType it = view.template begin< codim >(); it != endit; ++it )
        ++count;

      int lsetsize = lset.size(codim);
      if( count != lsetsize)
      {
        derr << "WARNING: walk = "<< count << " entities | set = "
             << lsetsize << " for codim " << codim << std::endl;
      }
      gridsize = count;
      // lsetsize should be at least the size of iterated entities
      assert( count <= gridsize );
    }

    {
      typedef typename GridView :: template Codim< 0 > :: Iterator Iterator;
      typedef typename Grid :: Traits :: LocalIdSet LocalIdSetType;
      typedef typename LocalIdSetType :: IdType IdType;

      std::set< IdType > entityfound;

      const Iterator endit = view.template end< 0 >();
      Iterator it = view.template begin< 0 >();
      if( it == endit )
        return;

      const LocalIdSetType &localIdSet = grid.localIdSet();
      for( ; it != endit; ++it )
      {
        const typename std::iterator_traits<Iterator>::value_type &entity = *it;
        if( !lset.contains( entity ) )
        {
          std::cerr << "Error: IndexSet does not contain all entities." << std::endl;
          assert( false );
        }
        const int subcount = entity.subEntities(codim);
        for( int i = 0; i < subcount; ++i )
        {
          const IdType id = localIdSet.id( entity.template subEntity< codim >( i ) );
          entityfound.insert( id );
        }
      }

      if( (size_t)gridsize != entityfound.size() )
      {
        std::cerr << "Warning: gridsize = " << gridsize << " entities"
                  << ", set of entities = " << entityfound.size()
                  << " [codim " << codim << "]" << std::endl;
      }

      // gridsize should be at least the size of found entities
      //assert( gridsize <= (int) entityfound.size() );
    }

    //******************************************************************

    typedef std::pair < int , GeometryType > SubEntityKeyType;
    //typedef std::map < int , std::pair<int,int> > subEntitymapType;
    std::map < SubEntityKeyType , std::vector<int> > setOfVerticesPerSubEntity;
    std::map < std::vector<int> , SubEntityKeyType > subEntityPerSetOfVertices;

    std::map < int , FieldVector<coordType,dimworld> > vertexCoordsMap;
    // setup vertex map , store vertex coords for vertex number
    {
      unsigned int count = 0;
      typedef typename GridView :: template Codim< dim > :: Iterator VxIterator;
      const VxIterator end = view.template end< dim >();
      for( VxIterator it = view.template begin< dim >(); it != end; ++it )
      {
        ++count;
        // get coordinates of vertex
        FieldVector< coordType, dimworld > vx ( it->geometry().corner( 0 ) );

        // get index of vertex
        sout << "Vertex " << vx << "\n";
        assert( lset.contains ( *it ) );
        int idx = lset.index( *it );

        sout << "Vertex " << idx << " = [" << vx << "]\n";

        // if vertex not in map insert it
        if( vertexCoordsMap.find(idx) == vertexCoordsMap.end())
          vertexCoordsMap[idx] = vx;
      }
      sout << "Found " << vertexCoordsMap.size() << " vertices for that index set!\n\n";

      // check whether size of map equals all found vertices
      assert( vertexCoordsMap.size() == count );

      // check whether size of vertices of set equals all found vertices
      sout << "Checking size of vertices "
           << count
           << " equals all found vertices "
           << (unsigned int)lset.size(Dune::GeometryType(0))
           << "\n";
      // assertion goes wrong:
      // - for parallel grid since no iteration over ghost subentities
      // - if there are vertices with geometry type 'none'
      //assert( count == (unsigned int)lset.size(Dune::GeometryType(0)) );
    }

    {
      typedef typename GridView :: template Codim< 0 > :: Iterator Iterator;
      // choose the right reference element
    #ifndef NDEBUG
      const Iterator refend = view.template end< 0 >();
    #endif
      Iterator refit = view.template begin< 0 >();
      assert( refit != refend );

      GeometryType type = refit->type();

      const ReferenceElement< coordType, dim > &refElem
        = ReferenceElements< coordType, dim >::general( type );

      // print dune reference element
      sout << "Dune reference element provides: \n";
      for(int i = 0; i < refElem.size( codim ); ++i )
      {
        sout << i << " subEntity [";
        int s = refElem.size(i,codim,dim);
        for(int j=0; j<s; j++)
        {
          sout << refElem.subEntity(i , codim , j , dim );
          if(j != s-1) sout << ",";
        }
        sout << "]\n";
      }
    }

    {
      typedef typename GridView :: template Codim< 0 > :: Iterator Iterator;
      const Iterator endit = view.template end< 0 >();
      for( Iterator it = view.template begin< 0 >(); it != endit; ++it )
      {
        // if (it->partitionType()==4) continue;
        sout << "****************************************\n";
        sout << "Element = " << lset.index(*it) << " on level " << it->level () << "\n";
        sout << "Vertices      = [";
        int svx = it->subEntities(dim);

        // print all vertex numbers
        for( int i = 0; i < svx; ++i )
        {
          const typename IndexSetType::IndexType idx = lset.subIndex( *it, i, dim );
          sout << idx << (i < svx-1 ? ", " : "]\n");
        }

        // print all vertex coordinates
        sout << "Vertex Coords = [";
        for( int i = 0; i < svx; ++i )
        {
          // get entity pointer of sub entity codim=dim (Vertex)
          typedef typename GridView::template Codim< dim >::Entity VertexE;
          VertexE vxE = it->template subEntity< dim >( i );

          // get coordinates of entity pointer
          FieldVector< coordType, dimworld > vx( vxE.geometry().corner( 0 ) );

          // output vertex coordinates
          sout << vx << (i < svx-1 ? ", " : "]\n");

          const typename IndexSetType::IndexType vxidx = lset.subIndex( *it, i, dim );

          // the subIndex and the index for subEntity must be the same
          if( vxidx != lset.index( vxE ) )
          {
            std::cerr << "Error: index( *subEntity< dim >( i ) ) != subIndex( entity, i, dim )" << std::endl;
            assert( vxidx == lset.index( vxE ) );
          }

          // check whether the coordinates are the same
          assert(vertexCoordsMap.find(vxidx)!=vertexCoordsMap.end());
          FieldVector<coordType,dimworld> vxcheck ( vertexCoordsMap[vxidx] );
          if( !compareVec( vxcheck, vx ) )
          {
            std::cerr << "Error: Inconsistent map of global vertex " << vxidx
                      << ": " << vxcheck << " != " << vx
                      << "  (type = " << it->partitionType() << ")" << std::endl;
            assert( compareVec( vxcheck, vx ) );
          }
        }

        ////////////////////////////////////////////////////////////
        // check sub entities
        ////////////////////////////////////////////////////////////
        checkSubEntity< codim >( grid, *it, lset, sout,
                                 setOfVerticesPerSubEntity, subEntityPerSetOfVertices, vertexCoordsMap );

        // check neighbors
        if( codim == 1 )
        {
          typedef typename GridView :: IntersectionIterator IntersectionIterator;

          if( !levelIndex || EnableLevelIntersectionIteratorCheck< typename GridView::Grid >::v )
          {
            const IntersectionIterator endnit = view.iend( *it );
            for( IntersectionIterator nit = view.ibegin( *it ); nit != endnit; ++nit )
            {
              if( !nit->neighbor() )
                continue;

              checkSubEntity< codim >( grid, nit->outside(), lset, sout,
                                       setOfVerticesPerSubEntity, subEntityPerSetOfVertices, vertexCoordsMap );
            }
          }
          else
          {
            static int called = 0;
            if( called++ == 0 )
              std::cerr << "Warning: Skipping index test using LevelIntersectionIterator." << std::endl;
          }
        }
      }
    }
  }


  template< class Grid, class GridView, class OutputStream, int codim, bool hasCodim >
  struct CheckIndexSet
  {
    static void checkIndexSet ( const Grid &grid, const GridView &view,
                                OutputStream &sout, bool levelIndex )
    {
      checkIndexSetForCodim< codim >( grid, view, sout, levelIndex );
      typedef Dune :: Capabilities :: hasEntity< Grid, codim-1 > hasNextCodim;
      CheckIndexSet< Grid, GridView, OutputStream, codim-1, hasNextCodim :: v >
      :: checkIndexSet( grid, view, sout, levelIndex );
    }
  };

  template< class Grid, class GridView, class OutputStream, int codim >
  struct CheckIndexSet< Grid, GridView, OutputStream, codim, false >
  {
    static void checkIndexSet ( const Grid &grid, const GridView &view,
                                OutputStream &sout, bool levelIndex )
    {
      derr << "WARNING: Entities for codim " << codim
           << " are not being tested!" << std::endl;
      typedef Dune :: Capabilities :: hasEntity< Grid, codim-1 > hasNextCodim;
      CheckIndexSet< Grid, GridView, OutputStream, codim-1, hasNextCodim :: v >
      :: checkIndexSet( grid, view, sout, levelIndex );
    }
  };

  template< class Grid, class GridView, class OutputStream >
  struct CheckIndexSet< Grid, GridView, OutputStream, 0, true >
  {
    static void checkIndexSet ( const Grid &grid, const GridView &view,
                                OutputStream &sout, bool levelIndex )
    {
      checkIndexSetForCodim< 0 >( grid, view, sout, levelIndex );
    }
  };

  template< class Grid, class GridView, class OutputStream >
  void checkIndexSet ( const Grid &grid, const GridView &view,
                       OutputStream &sout,  bool levelIndex = false )
  {
    CheckIndexSet< Grid, GridView, OutputStream, Grid :: dimension, true >
    :: checkIndexSet ( grid, view, sout, levelIndex );
  }

} // end namespace Dune

#endif // DUNE_GRID_TEST_CHECKINDEXSET_HH