void MeshfreeSolutionTransfer::transfer(const Variable & from_var, const Variable & to_var) { libmesh_experimental(); System * from_sys = from_var.system(); System * to_sys = to_var.system(); EquationSystems & from_es = from_sys->get_equation_systems(); MeshBase & from_mesh = from_es.get_mesh(); InverseDistanceInterpolation<LIBMESH_DIM> idi (from_mesh.comm(), 4, 2); std::vector<Point> & src_pts (idi.get_source_points()); std::vector<Number> & src_vals (idi.get_source_vals()); std::vector<std::string> field_vars; field_vars.push_back(from_var.name()); idi.set_field_variables(field_vars); // We now will loop over every node in the source mesh // and add it to a source point list, along with the solution { MeshBase::const_node_iterator nd = from_mesh.local_nodes_begin(); MeshBase::const_node_iterator end = from_mesh.local_nodes_end(); for (; nd!=end; ++nd) { const Node * node = *nd; src_pts.push_back(*node); src_vals.push_back((*from_sys->solution)(node->dof_number(from_sys->number(),from_var.number(),0))); } } // We have only set local values - prepare for use by gathering remote gata idi.prepare_for_use(); // Create a MeshlessInterpolationFunction that uses our // InverseDistanceInterpolation object. Since each // MeshlessInterpolationFunction shares the same // InverseDistanceInterpolation object in a threaded environment we // must also provide a locking mechanism. Threads::spin_mutex mutex; MeshlessInterpolationFunction mif(idi, mutex); // project the solution to_sys->project_solution(&mif); }
int main(int argc, char** argv) { // Skip this example if we do not meet certain requirements libmesh_example_requires(3 <= LIBMESH_DIM, "3D support"); #ifndef LIBMESH_HAVE_EIGEN libmesh_example_requires(false, "--enable-eigen"); #endif #ifndef LIBMESH_HAVE_ZLIB_H libmesh_example_requires(false, "--enable-zlib"); #endif // Initialize libMesh. LibMeshInit init (argc, argv); { // Demonstration case 1 { std::vector<Point> tgt_pts; std::vector<Number> tgt_data_idi, tgt_data_rbi; std::vector<std::string> field_vars; field_vars.push_back("u"); field_vars.push_back("v"); InverseDistanceInterpolation<3> idi (init.comm(), /* n_interp_pts = */ 8, /* power = */ 2); RadialBasisInterpolation<3> rbi (init.comm()); idi.set_field_variables (field_vars); rbi.set_field_variables (field_vars); create_random_point_cloud (100, idi.get_source_points()); // Explicitly set the data values we will interpolate from { const std::vector<Point> &src_pts (idi.get_source_points()); std::vector<Number> &src_vals (idi.get_source_vals()); src_vals.clear(); src_vals.reserve(2*src_pts.size()); for (std::vector<Point>::const_iterator pt_it=src_pts.begin(); pt_it != src_pts.end(); ++pt_it) { src_vals.push_back (exact_solution_u (*pt_it)); src_vals.push_back (exact_solution_v (*pt_it)); } } // give rbi the same info as idi rbi.get_source_points() = idi.get_source_points(); rbi.get_source_vals() = idi.get_source_vals(); idi.prepare_for_use(); rbi.prepare_for_use(); std::cout << idi; // Interpolate to some other random points, and evaluate the result { create_random_point_cloud (10, tgt_pts); //tgt_pts = rbi.get_source_points(); idi.interpolate_field_data (field_vars, tgt_pts, tgt_data_idi); rbi.interpolate_field_data (field_vars, tgt_pts, tgt_data_rbi); std::vector<Number>::const_iterator v_idi=tgt_data_idi.begin(), v_rbi=tgt_data_rbi.begin(); for (std::vector<Point>::const_iterator p_it=tgt_pts.begin(); p_it!=tgt_pts.end(); ++p_it) { std::cout << "\nAt target point " << *p_it << "\n u_interp_idi=" << *v_idi << ", u_interp_rbi=" << *v_rbi << ", u_exact=" << exact_solution_u(*p_it); ++v_idi; ++v_rbi; std::cout << "\n v_interp_idi=" << *v_idi << ", v_interp_rbi=" << *v_rbi << ", v_exact=" << exact_solution_v(*p_it) << std::endl; ++v_idi; ++v_rbi; } } } // Demonstration case 2 { Mesh mesh_a(init.comm()), mesh_b(init.comm()); mesh_a.read("struct.ucd.gz"); mesh_b.read("unstruct.ucd.gz"); // Create equation systems objects. EquationSystems es_a(mesh_a), es_b(mesh_b); System &sys_a = es_a.add_system<System>("src_system"), &sys_b = es_b.add_system<System>("dest_system"); sys_a.add_variable ("Cp", FIRST); sys_b.add_variable ("Cp", FIRST); sys_a.attach_init_function (init_sys); es_a.init(); // Write out the initial conditions. TecplotIO(mesh_a).write_equation_systems ("src.dat", es_a); InverseDistanceInterpolation<3> idi (init.comm(), /* n_interp_pts = */ 4, /* power = */ 2); RadialBasisInterpolation<3> rbi (init.comm()); std::vector<Point> &src_pts (idi.get_source_points()); std::vector<Number> &src_vals (idi.get_source_vals()); std::vector<std::string> field_vars; field_vars.push_back("Cp"); idi.set_field_variables(field_vars); // We now will loop over every node in the source mesh // and add it to a source point list, along with the solution { MeshBase::const_node_iterator nd = mesh_a.local_nodes_begin(); MeshBase::const_node_iterator end = mesh_a.local_nodes_end(); for (; nd!=end; ++nd) { const Node *node(*nd); src_pts.push_back(*node); src_vals.push_back(sys_a.current_solution(node->dof_number(0,0,0))); } rbi.set_field_variables(field_vars); rbi.get_source_points() = idi.get_source_points(); rbi.get_source_vals() = idi.get_source_vals(); } // We have only set local values - prepare for use by gathering remote gata idi.prepare_for_use(); rbi.prepare_for_use(); // Create a MeshlessInterpolationFunction that uses our InverseDistanceInterpolation // object. Since each MeshlessInterpolationFunction shares the same InverseDistanceInterpolation // object in a threaded environment we must also provide a locking mechanism. { Threads::spin_mutex mutex; MeshlessInterpolationFunction mif(idi, mutex); // project the solution onto system b es_b.init(); sys_b.project_solution (&mif); // Write the result TecplotIO(mesh_b).write_equation_systems ("dest_idi.dat", es_b); } // Create a MeshlessInterpolationFunction that uses our RadialBasisInterpolation // object. Since each MeshlessInterpolationFunction shares the same RadialBasisInterpolation // object in a threaded environment we must also provide a locking mechanism. { Threads::spin_mutex mutex; MeshlessInterpolationFunction mif(rbi, mutex); // project the solution onto system b sys_b.project_solution (&mif); // Write the result TecplotIO(mesh_b).write_equation_systems ("dest_rbi.dat", es_b); } } } return 0; }