1#ifndef DUNE_FEM_LINESEGMENTSAMPLER_HH
2#define DUNE_FEM_LINESEGMENTSAMPLER_HH
11#pragma GCC diagnostic push
12#pragma GCC diagnostic ignored "-Wdangling-reference"
17#include <dune/geometry/referenceelements.hh>
19#include <dune/fem/function/localfunction/const.hh>
43 template<
class Gr
idPart >
46 typedef GridPart GridPartType;
48 typedef typename GridPartType::GridType::ctype DomainFieldType;
50 static const int dimDomain = GridPartType::dimensionworld;
51 static const int dimGrid = GridPartType::dimension;
57 static_assert( dimDomain == dimGrid,
"LineSegmentSampler supports only flat grids." );
59 template<
class Vector >
76 : gridPart_(
gridPart ), left_( left ), right_( right )
89 template<
class Gr
idFunction >
90 void operator() (
const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples, std::nothrow_t )
const;
91 template<
class Gr
idFunction >
92 void operator() (
const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples )
const;
103 void samplePoints( std::vector< DomainType > &points )
const;
106 const GridPart &
gridPart ()
const {
return gridPart_; }
109 const GridPart &gridPart_;
110 DomainType left_, right_;
118 template<
class Gr
idPart >
119 template<
class Vector >
120 struct LineSegmentSampler< GridPart >::Reduce
122 Vector
operator() (
const Vector &a,
const Vector &b )
const
125 for(
int k = 0; k < Vector::dimension; ++k )
126 c[ k ] = std::min( a[ k ], b[ k ] );
136 template<
class Gr
idPart >
137 template<
class Gr
idFunction >
139 ::operator() (
const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples,
140 std::nothrow_t)
const
144 typedef typename GridFunction::RangeFieldType RangeFieldType;
146 const int numSamples = samples.size();
148 ds /= DomainFieldType( numSamples - 1 );
150 const RangeFieldType invalid
151 = std::numeric_limits< RangeFieldType >::quiet_NaN();
152 for(
int i = 0; i < numSamples; ++i )
153 samples[ i ] =
typename GridFunction::RangeType( invalid );
155 ConstLocalFunction< GridFunction > lf( f );
156 const IteratorType end = gridPart().template end< 0 >();
157 for( IteratorType it = gridPart().
template begin< 0 >(); it != end; ++it )
159 const EntityType &entity = *it;
160 const typename EntityType::Geometry &geometry = entity.geometry();
166 const int numFaces = refElement.size( 1 );
167 for(
int face = 0; face < numFaces; ++face )
172 const LocalDomainType &refNormal = refElement.integrationOuterNormal( face );
173 geometry.jacobianInverseTransposed( center ).mv( refNormal, normal );
175 const DomainFieldType nds = normal * ds;
176 const DomainFieldType ncl = normal * (geometry.global( center ) - left_);
177 if( std::abs( nds ) > 1e-8 )
180 const DomainFieldType alpha = ncl / nds;
182 upper = std::min( upper, alpha );
184 lower = std::max( lower, alpha );
186 else if( ncl < -1e-8 )
197 const int i_upper = std::min(
int( std::floor( upper + 1e-8 ) ), numSamples - 1 );
198 const int i_lower = std::max(
int( std::ceil( lower - 1e-8 ) ), 0 );
199 for(
int i = i_lower; i <= i_upper; ++i )
202 xi.
axpy( DomainFieldType( i ), ds );
203 lf.evaluate( geometry.local( xi ), samples[ i ] );
209 typedef Reduce< typename GridFunction::RangeType > Op;
210 gridPart().comm().template allreduce< Op >( &(samples[ 0 ]), numSamples );
213 template<
class Gr
idPart >
214 template<
class Gr
idFunction >
216 ::operator() (
const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples )
const
218 typedef typename GridFunction::RangeFieldType RangeFieldType;
220 const int numSamples = samples.size();
224 (*this)(f,samples,std::nothrow);
229 if constexpr ( std::is_floating_point< RangeFieldType >::value )
231 for(
int i = 0; i < numSamples; ++i )
233 for(
int d=0; d<GridFunction::RangeType::dimension; ++d )
235 valid &= ! std::isnan( samples[ i ][ d ] );
241 DUNE_THROW( InvalidStateException,
"LineSegmentSampler could not find all samples." );
245 template<
class Gr
idPart >
249 const int numSamples = points.size();
253 ds /= DomainFieldType( numSamples - 1 );
255 for(
int i = 0; i < numSamples; ++i )
258 points[ i ].axpy( DomainFieldType( i ), ds );
267#pragma GCC diagnostic pop
constexpr derived_type & axpy(const field_type &a, const DenseVector< Other > &x)
vector space axpy operation ( *this += a x )
Definition: densevector.hh:576
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:375
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:485
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:507
Dune namespace.
Definition: alignedallocator.hh:13
Static tag representing a codimension.
Definition: dimension.hh:24
samples values of a discrete function along a given line segment
Definition: linesegmentsampler.hh:45
void samplePoints(std::vector< DomainType > &points) const
returns sampling points
Definition: linesegmentsampler.hh:247
const GridPart & gridPart() const
obtain grid part on which the LineSegmentSampler works
Definition: linesegmentsampler.hh:106
void operator()(const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples, std::nothrow_t) const
sample a given function
Definition: linesegmentsampler.hh:139
LineSegmentSampler(const GridPart &gridPart, const DomainType &left, const DomainType &right)
constructor
Definition: linesegmentsampler.hh:75
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:128
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:156