libMesh
radial_basis_interpolation.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_RADIAL_BASIS_INTERPOLATION_H
21 #define LIBMESH_RADIAL_BASIS_INTERPOLATION_H
22 
23 // Local includes
24 #include "libmesh/libmesh_config.h"
25 #include "libmesh/libmesh_common.h"
26 #include "libmesh/meshfree_interpolation.h"
27 #include "libmesh/radial_basis_functions.h"
28 #include "libmesh/bounding_box.h"
29 
30 
31 
32 namespace libMesh
33 {
34 
42 template <unsigned int KDDim, class RBF = WendlandRBF<KDDim, 2>>
44 {
51 
52 protected:
53 
58 
62  std::vector<Number> _weights;
63 
68 
73 
74 public:
75 
80  Real radius=-1) :
81  InverseDistanceInterpolation<KDDim> (comm_in,8,2),
82  _r_bbox(0.),
84  { }
85 
90  virtual void clear() override;
91 
95  virtual void prepare_for_use () override;
96 
101  virtual void interpolate_field_data (const std::vector<std::string> & field_names,
102  const std::vector<Point> & tgt_pts,
103  std::vector<Number> & tgt_vals) const override;
104 };
105 
106 } // namespace libMesh
107 
108 
109 #endif // #define LIBMESH_RADIAL_BASIS_INTERPOLATION_H
Inverse distance interpolation.
const Real radius
Radial Basis Function interpolation.
virtual void prepare_for_use() override
Prepares data structures for use.
BoundingBox _src_bbox
Bounding box for our source points.
virtual void interpolate_field_data(const std::vector< std::string > &field_names, const std::vector< Point > &tgt_pts, std::vector< Number > &tgt_vals) const override
Interpolate source data at target points.
Real _r_bbox
Diagonal of the bounding box.
The libMesh namespace provides an interface to certain functionality in the library.
std::vector< Number > _weights
basis coefficients.
RadialBasisInterpolation(const libMesh::Parallel::Communicator &comm_in, Real radius=-1)
Constructor.
virtual void clear() override
Clears all internal data structures and restores to a pristine state.
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real