28 #include "libmesh/libmesh.h" 29 #include "libmesh/meshfree_interpolation.h" 30 #include "libmesh/radial_basis_interpolation.h" 31 #include "libmesh/mesh.h" 32 #include "libmesh/mesh_refinement.h" 33 #include "libmesh/equation_systems.h" 34 #include "libmesh/numeric_vector.h" 35 #include "libmesh/tecplot_io.h" 36 #include "libmesh/threads.h" 37 #include "libmesh/node.h" 38 #include "libmesh/getpot.h" 39 #include "libmesh/meshfree_interpolation_function.h" 50 std::vector<Point> & pts,
51 const Real max_range = 10)
53 libMesh::out <<
"Generating "<< Npts <<
" local point cloud...";
56 for (
size_t i=0;i<Npts;i++)
58 pts[i](0) = max_range * (std::rand() % 1000) /
Real(1000);
59 pts[i](1) = max_range * (std::rand() % 1000) /
Real(1000);
60 pts[i](2) = max_range * (std::rand() % 1000) /
Real(1000);
106 const std::string & system_name)
118 int main(
int argc,
char ** argv)
121 libmesh_example_requires(3 <= LIBMESH_DIM,
"3D support");
122 #ifndef LIBMESH_HAVE_EIGEN 123 libmesh_example_requires(
false,
"--enable-eigen");
125 #ifndef LIBMESH_HAVE_ZLIB_H 126 libmesh_example_requires(
false,
"--enable-zlib");
129 #ifdef LIBMESH_DEFAULT_TRIPLE_PRECISION 130 libmesh_example_requires(
false,
"--disable-triple-precision");
132 #ifdef LIBMESH_DEFAULT_QUADRUPLE_PRECISION 133 libmesh_example_requires(
false,
"--disable-quadruple-precision");
139 GetPot input(argc, argv);
143 std::vector<Point> tgt_pts;
144 std::vector<Number> tgt_data_idi, tgt_data_rbi;
145 std::vector<std::string> field_vars {
"u",
"v"};
154 rbi.set_field_variables (field_vars);
156 const int n_source_points = input(
"n_source_points", 100);
159 const int my_n_source_points =
160 n_source_points /
init.comm().size() +
161 (
init.comm().rank() < n_source_points %
init.comm().size());
164 idi.get_source_points());
169 const std::vector<Point> & src_pts (idi.get_source_points());
170 std::vector<Number> & src_vals (idi.get_source_vals());
172 src_vals.clear(); src_vals.reserve(2*src_pts.size());
174 for (std::vector<Point>::const_iterator pt_it=src_pts.begin();
175 pt_it != src_pts.end(); ++pt_it)
183 rbi.get_source_points() = idi.get_source_points();
184 rbi.get_source_vals() = idi.get_source_vals();
186 idi.prepare_for_use();
187 rbi.prepare_for_use();
193 const int n_target_points = input(
"n_target_points", 100);
200 idi.interpolate_field_data (field_vars,
204 rbi.interpolate_field_data (field_vars,
208 std::vector<Number>::const_iterator
209 v_idi = tgt_data_idi.begin(),
210 v_rbi = tgt_data_rbi.begin();
212 for (std::vector<Point>::const_iterator p_it=tgt_pts.begin();
213 p_it!=tgt_pts.end(); ++p_it)
216 <<
"\n u_interp_idi=" << *v_idi
217 <<
", u_interp_rbi=" << *v_rbi
222 <<
", v_interp_rbi=" << *v_rbi
236 mesh_a.
read(
"struct.ucd.gz");
237 mesh_b.read(
"unstruct.ucd.gz");
240 const int n_refinements = input(
"n_refinements", 0);
243 #ifndef LIBMESH_ENABLE_AMR 244 libmesh_example_requires(n_refinements==0,
"--enable-amr");
254 es_a(mesh_a), es_b(mesh_b);
260 sys_b.add_variable (
"Cp",
FIRST);
274 std::vector<Point> & src_pts (idi.get_source_points());
275 std::vector<Number> & src_vals (idi.get_source_vals());
280 for (
const auto & node : mesh_a.local_node_ptr_range())
282 src_pts.push_back(*node);
286 rbi.set_field_variables({
"Cp"});
287 rbi.get_source_points() = idi.get_source_points();
288 rbi.get_source_vals() = idi.get_source_vals();
291 idi.prepare_for_use();
292 rbi.prepare_for_use();
303 sys_b.project_solution (&mif);
318 sys_b.project_solution (&mif);
This is the EquationSystems class.
Inverse distance interpolation.
Radial Basis Function interpolation.
This class provides the ability to map between arbitrary, user-defined strings and several data types...
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Real exact_solution_v(const Point &p)
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
The libMesh namespace provides an interface to certain functionality in the library.
const T_sys & get_system(std::string_view name) const
Number current_solution(const dof_id_type global_dof_number) const
void project_solution(FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr) const
Projects arbitrary functions onto the current solution.
void init_sys(EquationSystems &es, const std::string &system_name)
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Manages consistently variables, degrees of freedom, and coefficient vectors.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Real exact_solution_u(const Point &p)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void set_field_variables(std::vector< std::string > names)
Defines the field variable(s) we are responsible for, and importantly their assumed ordering...
Parameters parameters
Data structure holding arbitrary parameters.
This class implements writing meshes in the Tecplot format.
virtual void init()
Initialize all the systems.
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false) override
Reads the file specified by name.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
A Point defines a location in LIBMESH_DIM dimensional Real space.
void attach_init_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in initializing the system.
Number exact_value(const Point &p, const Parameters &, const std::string &, const std::string &)
int main(int argc, char **argv)
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
void create_random_point_cloud(const unsigned int Npts, std::vector< Point > &pts, const Real max_range=10)