DUNE-FEM (unstable)

linesegmentsampler.hh
1#ifndef DUNE_FEM_LINESEGMENTSAMPLER_HH
2#define DUNE_FEM_LINESEGMENTSAMPLER_HH
3
4#include <limits>
5#include <vector>
6#include <cmath>
7#include <type_traits>
8#include <new>
9
10#if __GNUC__ >= 13
11#pragma GCC diagnostic push
12#pragma GCC diagnostic ignored "-Wdangling-reference"
13#endif
14
16
17#include <dune/geometry/referenceelements.hh>
18
19#include <dune/fem/function/localfunction/const.hh>
20
21namespace Dune
22{
23
24 namespace Fem
25 {
26
27 // LineSegmentSampler
28 // ------------------
29
43 template< class GridPart >
45 {
46 typedef GridPart GridPartType;
47
48 typedef typename GridPartType::GridType::ctype DomainFieldType;
49
50 static const int dimDomain = GridPartType::dimensionworld;
51 static const int dimGrid = GridPartType::dimension;
52
55
56 private:
57 static_assert( dimDomain == dimGrid, "LineSegmentSampler supports only flat grids." );
58
59 template< class Vector >
60 struct Reduce;
61
63
64 public:
75 LineSegmentSampler ( const GridPart &gridPart, const DomainType &left, const DomainType &right )
76 : gridPart_( gridPart ), left_( left ), right_( right )
77 {}
78
89 template< class GridFunction >
90 void operator() ( const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples, std::nothrow_t ) const;
91 template< class GridFunction >
92 void operator() ( const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples ) const;
93
103 void samplePoints( std::vector< DomainType > &points ) const;
104
106 const GridPart &gridPart () const { return gridPart_; }
107
108 private:
109 const GridPart &gridPart_;
110 DomainType left_, right_;
111 };
112
113
114
115 // LineSegmentSampler::Reduce
116 // --------------------------
117
118 template< class GridPart >
119 template< class Vector >
120 struct LineSegmentSampler< GridPart >::Reduce
121 {
122 Vector operator() ( const Vector &a, const Vector &b ) const
123 {
124 Vector c;
125 for( int k = 0; k < Vector::dimension; ++k )
126 c[ k ] = std::min( a[ k ], b[ k ] );
127 return c;
128 }
129 };
130
131
132
133 // Implementation of LineSegmentSampler
134 // ------------------------------------
135
136 template< class GridPart >
137 template< class GridFunction >
139 ::operator() ( const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples,
140 std::nothrow_t) const
141 {
142 typedef typename GridPartType::template Codim< 0 >::IteratorType IteratorType;
143 typedef typename GridPartType::template Codim< 0 >::EntityType EntityType;
144 typedef typename GridFunction::RangeFieldType RangeFieldType;
145
146 const int numSamples = samples.size();
147 DomainType ds = right_ - left_;
148 ds /= DomainFieldType( numSamples - 1 );
149
150 // use max value as invalid to fill values from other cores
151 const RangeFieldType invalid
153 for( int i = 0; i < numSamples; ++i )
154 samples[ i ] = typename GridFunction::RangeType( invalid );
155
156 ConstLocalFunction< GridFunction > lf( f );
157 const IteratorType end = gridPart().template end< 0 >();
158 for( IteratorType it = gridPart().template begin< 0 >(); it != end; ++it )
159 {
160 const EntityType &entity = *it;
161 const typename EntityType::Geometry &geometry = entity.geometry();
162
163 DomainFieldType lower = std::numeric_limits< DomainFieldType >::min();
164 DomainFieldType upper = std::numeric_limits< DomainFieldType >::max();
165
166 auto &refElement = ReferenceElements::general( geometry.type() );
167 const int numFaces = refElement.size( 1 );
168 for( int face = 0; face < numFaces; ++face )
169 {
170 const LocalDomainType &center = refElement.position( face, 1 );
171
172 DomainType normal;
173 const LocalDomainType &refNormal = refElement.integrationOuterNormal( face );
174 geometry.jacobianInverseTransposed( center ).mv( refNormal, normal );
175
176 const DomainFieldType nds = normal * ds;
177 const DomainFieldType ncl = normal * (geometry.global( center ) - left_);
178 if( std::abs( nds ) > 1e-8 )
179 {
180 // ds is not parallel to the face
181 const DomainFieldType alpha = ncl / nds;
182 if( nds > 0 )
183 upper = std::min( upper, alpha );
184 else
185 lower = std::max( lower, alpha );
186 }
187 else if( ncl < -1e-8 )
188 {
189 // ds is parallel to the face and the line lies on the outside
192 }
193 }
194
195 if( lower <= upper )
196 {
197 lf.bind( entity );
198 const int i_upper = std::min( int( std::floor( upper + 1e-8 ) ), numSamples - 1 );
199 const int i_lower = std::max( int( std::ceil( lower - 1e-8 ) ), 0 );
200 for( int i = i_lower; i <= i_upper; ++i )
201 {
202 DomainType xi = left_;
203 xi.axpy( DomainFieldType( i ), ds );
204 lf.evaluate( geometry.local( xi ), samples[ i ] );
205 }
206 lf.unbind();
207 }
208 }
209
210 typedef Reduce< typename GridFunction::RangeType > Op;
211 gridPart().comm().template allreduce< Op >( &(samples[ 0 ]), numSamples );
212 }
213
214 template< class GridPart >
215 template< class GridFunction >
217 ::operator() ( const GridFunction &f, std::vector< typename GridFunction::RangeType > &samples ) const
218 {
219 typedef typename GridFunction::RangeFieldType RangeFieldType;
220
221 const int numSamples = samples.size();
222 if( numSamples < 2 )
223 DUNE_THROW( InvalidStateException, "LineSegmentSampler cannot sample less than 2 points." );
224
225 (*this)(f,samples,std::nothrow);
226
227 bool valid = true;
228
229 // only use isnan when field type is a floating point
230 if constexpr ( std::is_floating_point< RangeFieldType >::value )
231 {
232 for( int i = 0; i < numSamples; ++i )
233 {
234 for( int d=0; d<GridFunction::RangeType::dimension; ++d )
235 {
236 valid &= ! std::isnan( samples[ i ][ d ] );
237 }
238 }
239 }
240
241 if( !valid )
242 DUNE_THROW( InvalidStateException, "LineSegmentSampler could not find all samples." );
243 }
244
245
246 template< class GridPart >
248 :: samplePoints ( std::vector< DomainType > &points ) const
249 {
250 const int numSamples = points.size();
251 if( numSamples < 2 )
252 DUNE_THROW( InvalidStateException, "LineSegmentSampler cannot sample less than 2 points." );
253 DomainType ds = right_ - left_;
254 ds /= DomainFieldType( numSamples - 1 );
255
256 for( int i = 0; i < numSamples; ++i )
257 {
258 points[ i ] = left_;
259 points[ i ].axpy( DomainFieldType( i ), ds );
260 }
261 }
262
263 } // namespace Fem
264
265} // namespace Dune
266
267#if __GNUC__ >= 13
268#pragma GCC diagnostic pop
269#endif
270
271#endif // #ifndef DUNE_FEM_LINESEGMENTSAMPLER_HH
constexpr derived_type & axpy(const field_type &a, const DenseVector< Other > &x)
vector space axpy operation ( *this += a x )
Definition: densevector.hh:590
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:489
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:511
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:248
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Jan 9, 23:34, 2026)