dune-fem  2.4.1-rc
linesegmentsampler.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_LINESEGMENTSAMPLER_HH
2 #define DUNE_FEM_LINESEGMENTSAMPLER_HH
3 
4 #include <limits>
5 #include <vector>
6 
7 #include <dune/common/fvector.hh>
8 
9 #include <dune/geometry/referenceelements.hh>
10 
11 namespace Dune
12 {
13 
14  namespace Fem
15  {
16 
17  // LineSegmentSampler
18  // ------------------
19 
33  template< class GridPart >
35  {
36  typedef GridPart GridPartType;
37 
38  typedef typename GridPartType::GridType::ctype DomainFieldType;
39 
40  static const int dimDomain = GridPartType::dimensionworld;
41  static const int dimGrid = GridPartType::dimension;
42 
43  typedef FieldVector< DomainFieldType, dimDomain > DomainType;
44  typedef FieldVector< DomainFieldType, dimGrid > LocalDomainType;
45 
46  private:
47  static_assert( dimDomain == dimGrid, "LineSegmentSampler supports only flat grids." );
48 
49  template< class Vector >
50  struct Reduce;
51 
52  typedef Dune::ReferenceElement< DomainFieldType, dimGrid > ReferenceElement;
53  typedef Dune::ReferenceElements< DomainFieldType, dimGrid > ReferenceElements;
54 
55  public:
66  LineSegmentSampler ( const GridPart &gridPart, const DomainType &left, const DomainType &right )
67  : gridPart_( gridPart ), left_( left ), right_( right )
68  {}
69 
80  template< class GridFunction >
81  void operator() ( const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples ) const;
82 
92  void samplePoints( std::vector< DomainType > &points ) const;
93 
95  const GridPart &gridPart () const { return gridPart_; }
96 
97  private:
98  const GridPart &gridPart_;
99  DomainType left_, right_;
100  };
101 
102 
103 
104  // LineSegmentSampler::Reduce
105  // --------------------------
106 
107  template< class GridPart >
108  template< class Vector >
109  struct LineSegmentSampler< GridPart >::Reduce
110  {
111  Vector operator() ( const Vector &a, const Vector &b ) const
112  {
113  Vector c;
114  for( int k = 0; k < Vector::dimension; ++k )
115  c[ k ] = std::min( a[ k ], b[ k ] );
116  return c;
117  }
118  };
119 
120 
121 
122  // Implementation of LineSegmentSampler
123  // ------------------------------------
124 
125  template< class GridPart >
126  template< class GridFunction >
128  ::operator() ( const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples ) const
129  {
130  typedef typename GridPartType::template Codim< 0 >::IteratorType IteratorType;
131  typedef typename GridPartType::template Codim< 0 >::EntityType EntityType;
132 
133  typedef typename GridFunction::LocalFunctionType LocalFunctionType;
134 
135  const int numSamples = samples.size();
136  if( numSamples < 2 )
137  DUNE_THROW( InvalidStateException, "LineSegmentSampler cannot sample less than 2 points." );
138  DomainType ds = right_ - left_;
139  ds /= DomainFieldType( numSamples - 1 );
140 
141  const typename GridFunction::RangeFieldType invalid
142  = std::numeric_limits< typename GridFunction::RangeFieldType >::infinity();
143  for( int i = 0; i < numSamples; ++i )
144  samples[ i ] = typename GridFunction::RangeType( invalid );
145 
146  const IteratorType end = gridPart().template end< 0 >();
147  for( IteratorType it = gridPart().template begin< 0 >(); it != end; ++it )
148  {
149  const EntityType &entity = *it;
150  const typename EntityType::Geometry &geometry = entity.geometry();
151 
154 
155  const ReferenceElement &refElement = ReferenceElements::general( geometry.type() );
156  const int numFaces = refElement.size( 1 );
157  for( int face = 0; face < numFaces; ++face )
158  {
159  const LocalDomainType &center = refElement.position( face, 1 );
160 
161  DomainType normal;
162  const LocalDomainType &refNormal = refElement.integrationOuterNormal( face );
163  geometry.jacobianInverseTransposed( center ).mv( refNormal, normal );
164 
165  const DomainFieldType nds = normal * ds;
166  const DomainFieldType ncl = normal * (geometry.global( center ) - left_);
167  if( std::abs( nds ) > 1e-8 )
168  {
169  // ds is not parallel to the face
170  const DomainFieldType alpha = ncl / nds;
171  if( nds > 0 )
172  upper = std::min( upper, alpha );
173  else
174  lower = std::max( lower, alpha );
175  }
176  else if( ncl < -1e-8 )
177  {
178  // ds is parallel to the face and the line lies on the outside
181  }
182  }
183 
184  if( lower <= upper )
185  {
186  const LocalFunctionType &lf = f.localFunction( entity );
187  const int i_upper = std::min( int( std::floor( upper + 1e-8 ) ), numSamples - 1 );
188  const int i_lower = std::max( int( std::ceil( lower - 1e-8 ) ), 0 );
189  for( int i = i_lower; i <= i_upper; ++i )
190  {
191  DomainType xi = left_;
192  xi.axpy( DomainFieldType( i ), ds );
193  lf.evaluate( geometry.local( xi ), samples[ i ] );
194  }
195  }
196  }
197 
198  typedef Reduce< typename GridFunction::RangeType > Op;
199  gridPart().comm().template allreduce< Op >( &(samples[ 0 ]), numSamples );
200 
201  bool valid = true;
202  for( int i = 0; i < numSamples; ++i )
203  valid &= (samples[ i ] != typename GridFunction::RangeType( invalid ));
204  if( !valid )
205  DUNE_THROW( InvalidStateException, "LineSegmentSampler could not find all samples." );
206  }
207 
208 
209  template< class GridPart >
211  :: samplePoints ( std::vector< DomainType > &points ) const
212  {
213  const int numSamples = points.size();
214  if( numSamples < 2 )
215  DUNE_THROW( InvalidStateException, "LineSegmentSampler cannot sample less than 2 points." );
216  DomainType ds = right_ - left_;
217  ds /= DomainFieldType( numSamples - 1 );
218 
219  for( int i = 0; i < numSamples; ++i )
220  {
221  points[ i ] = left_;
222  points[ i ].axpy( DomainFieldType( i ), ds );
223  }
224  }
225 
226  } // namespace Fem
227 
228 } // namespace Dune
229 
230 #endif // #ifndef DUNE_FEM_LINESEGMENTSAMPLER_HH
static const int dimGrid
Definition: linesegmentsampler.hh:41
GridPartType::GridType::ctype DomainFieldType
Definition: linesegmentsampler.hh:38
static constexpr T min(T a)
Definition: utility.hh:81
samples values of a discrete function along a given line segment
Definition: linesegmentsampler.hh:34
FieldVector< DomainFieldType, dimGrid > LocalDomainType
Definition: linesegmentsampler.hh:44
static constexpr T max(T a)
Definition: utility.hh:65
static const int dimDomain
Definition: linesegmentsampler.hh:40
GridPart GridPartType
Definition: linesegmentsampler.hh:36
Double abs(const Double &a)
Definition: double.hh:860
Definition: coordinate.hh:4
const GridPart & gridPart() const
obtain grid part on which the LineSegmentSampler works
Definition: linesegmentsampler.hh:95
LineSegmentSampler(const GridPart &gridPart, const DomainType &left, const DomainType &right)
constructor
Definition: linesegmentsampler.hh:66
FieldVector< DomainFieldType, dimDomain > DomainType
Definition: linesegmentsampler.hh:43
void samplePoints(std::vector< DomainType > &points) const
returns sampling points
Definition: linesegmentsampler.hh:211
void operator()(const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples) const
sample a given function
Definition: linesegmentsampler.hh:128