libMesh
|
This class provides function-like objects for data distributed over a mesh. More...
#include <mesh_function.h>
Public Member Functions | |
MeshFunction (const EquationSystems &eqn_systems, const NumericVector< Number > &vec, const DofMap &dof_map, std::vector< unsigned int > vars, const FunctionBase< Number > *master=nullptr) | |
Constructor for mesh based functions with vectors as return value. More... | |
MeshFunction (const EquationSystems &eqn_systems, const NumericVector< Number > &vec, const DofMap &dof_map, const unsigned int var, const FunctionBase< Number > *master=nullptr) | |
Constructor for mesh based functions with a number as return value. More... | |
MeshFunction (const MeshFunction &mf) | |
A regular copy constructor. More... | |
MeshFunction (MeshFunction &&)=default | |
Special functions. More... | |
MeshFunction & | operator= (const MeshFunction &)=delete |
MeshFunction & | operator= (MeshFunction &&)=delete |
~MeshFunction () | |
Destructor. More... | |
virtual void | init () override |
Override the FunctionBase::init() member function. More... | |
void | init (const Trees::BuildType point_locator_build_type) |
The actual initialization process. More... | |
virtual void | clear () override |
Clears the function. More... | |
virtual std::unique_ptr< FunctionBase< Number > > | clone () const override |
virtual Number | operator() (const Point &p, const Real time=0.) override |
std::map< const Elem *, Number > | discontinuous_value (const Point &p, const Real time=0.) |
Gradient | gradient (const Point &p, const Real time=0.) |
std::map< const Elem *, Gradient > | discontinuous_gradient (const Point &p, const Real time=0.) |
Tensor | hessian (const Point &p, const Real time=0.) |
virtual void | operator() (const Point &p, const Real time, DenseVector< Number > &output) override |
Computes values at coordinate p and for time time , which defaults to zero, optionally restricting the point to the MeshFunction subdomain_ids. More... | |
void | operator() (const Point &p, const Real time, DenseVector< Number > &output, const std::set< subdomain_id_type > *subdomain_ids) |
Computes values at coordinate p and for time time , restricting the point to the passed subdomain_ids, which parameter overrides the internal subdomain_ids. More... | |
void | discontinuous_value (const Point &p, const Real time, std::map< const Elem *, DenseVector< Number >> &output) |
Similar to operator() with the same parameter list, but with the difference that multiple values on faces are explicitly permitted. More... | |
void | discontinuous_value (const Point &p, const Real time, std::map< const Elem *, DenseVector< Number >> &output, const std::set< subdomain_id_type > *subdomain_ids) |
Similar to operator() with the same parameter list, but with the difference that multiple values on faces are explicitly permitted. More... | |
void | gradient (const Point &p, const Real time, std::vector< Gradient > &output) |
Computes gradients at coordinate p and for time time , which defaults to zero, optionally restricting the point to the MeshFunction subdomain_ids. More... | |
void | gradient (const Point &p, const Real time, std::vector< Gradient > &output, const std::set< subdomain_id_type > *subdomain_ids) |
Computes gradients at coordinate p and for time time , which defaults to zero, optionally restricting the point to the passed subdomain_ids, which parameter overrides the internal subdomain_ids. More... | |
void | discontinuous_gradient (const Point &p, const Real time, std::map< const Elem *, std::vector< Gradient >> &output) |
Similar to gradient, but with the difference that multiple values on faces are explicitly permitted. More... | |
void | discontinuous_gradient (const Point &p, const Real time, std::map< const Elem *, std::vector< Gradient >> &output, const std::set< subdomain_id_type > *subdomain_ids) |
Similar to gradient, but with the difference that multiple values on faces are explicitly permitted. More... | |
void | hessian (const Point &p, const Real time, std::vector< Tensor > &output) |
Computes gradients at coordinate p and for time time , which defaults to zero, optionally restricting the point to the MeshFunction subdomain_ids. More... | |
void | hessian (const Point &p, const Real time, std::vector< Tensor > &output, const std::set< subdomain_id_type > *subdomain_ids) |
Computes gradients at coordinate p and for time time , which defaults to zero, optionally restricting the point to the passed subdomain_ids. More... | |
const PointLocatorBase & | get_point_locator () const |
PointLocatorBase & | get_point_locator () |
void | enable_out_of_mesh_mode (const DenseVector< Number > &value) |
Enables out-of-mesh mode. More... | |
void | enable_out_of_mesh_mode (const Number &value) |
Enables out-of-mesh mode. More... | |
void | disable_out_of_mesh_mode () |
Disables out-of-mesh mode. More... | |
void | set_point_locator_tolerance (Real tol) |
We may want to specify a tolerance for the PointLocator to use, since in some cases the point we want to evaluate at might be slightly outside the mesh (due to numerical rounding issues, for example). More... | |
void | unset_point_locator_tolerance () |
Turn off the user-specified PointLocator tolerance. More... | |
void | set_subdomain_ids (const std::set< subdomain_id_type > *subdomain_ids) |
Choose a default list of subdomain ids to be searched for points. More... | |
void | operator() (const Point &p, DenseVector< Number > &output) |
Evaluation function for time-independent vector-valued functions. More... | |
virtual Number | component (unsigned int i, const Point &p, Real time=0.) |
bool | initialized () const |
void | set_is_time_dependent (bool is_time_dependent) |
Function to set whether this is a time-dependent function or not. More... | |
bool | is_time_dependent () const |
const Parallel::Communicator & | comm () const |
processor_id_type | n_processors () const |
processor_id_type | processor_id () const |
Protected Member Functions | |
const Elem * | find_element (const Point &p, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const |
Helper function to reduce code duplication. More... | |
std::set< const Elem * > | find_elements (const Point &p, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const |
const Elem * | check_found_elem (const Elem *element, const Point &p) const |
Helper function that is called by MeshFunction::find_element() and MeshFunction::find_elements() to ensure that Elems found by the PointLocator are actually "evaluable" on the processor where they are found. More... | |
void | _gradient_on_elem (const Point &p, const Elem *element, std::vector< Gradient > &output) |
Helper function for finding a gradient as evaluated from a specific element. More... | |
Protected Attributes | |
const EquationSystems & | _eqn_systems |
The equation systems handler, from which the data are gathered. More... | |
const NumericVector< Number > & | _vector |
A reference to the vector that holds the data that is to be interpolated. More... | |
const DofMap & | _dof_map |
Need access to the DofMap of the other system. More... | |
const std::vector< unsigned int > | _system_vars |
The indices of the variables within the other system for which data are to be gathered. More... | |
std::unique_ptr< PointLocatorBase > | _point_locator |
A point locator is needed to locate the points in the mesh. More... | |
std::unique_ptr< std::set< subdomain_id_type > > | _subdomain_ids |
A default set of subdomain ids in which to search for points. More... | |
bool | _out_of_mesh_mode |
true if out-of-mesh mode is enabled. More... | |
DenseVector< Number > | _out_of_mesh_value |
Value to return outside the mesh if out-of-mesh mode is enabled. More... | |
const FunctionBase * | _master |
Const pointer to our master, initialized to nullptr . More... | |
bool | _initialized |
When init() was called so that everything is ready for calls to operator() (...), then this bool is true. More... | |
bool | _is_time_dependent |
Cache whether or not this function is actually time-dependent. More... | |
const Parallel::Communicator & | _communicator |
This class provides function-like objects for data distributed over a mesh.
Definition at line 54 of file mesh_function.h.
libMesh::MeshFunction::MeshFunction | ( | const EquationSystems & | eqn_systems, |
const NumericVector< Number > & | vec, | ||
const DofMap & | dof_map, | ||
std::vector< unsigned int > | vars, | ||
const FunctionBase< Number > * | master = nullptr |
||
) |
Constructor for mesh based functions with vectors as return value.
Optionally takes a master function. If the MeshFunction is to be evaluated outside of the local partition of the mesh, then both the mesh in eqn_systems
and the coefficient vector vec
should be serialized.
Definition at line 42 of file mesh_function.C.
libMesh::MeshFunction::MeshFunction | ( | const EquationSystems & | eqn_systems, |
const NumericVector< Number > & | vec, | ||
const DofMap & | dof_map, | ||
const unsigned int | var, | ||
const FunctionBase< Number > * | master = nullptr |
||
) |
Constructor for mesh based functions with a number as return value.
Optionally takes a master function. If the MeshFunction is to be evaluated outside of the local partition of the mesh, then both the mesh in eqn_systems
and the coefficient vector vec
should be serialized.
Definition at line 59 of file mesh_function.C.
libMesh::MeshFunction::MeshFunction | ( | const MeshFunction & | mf | ) |
A regular copy constructor.
Definition at line 74 of file mesh_function.C.
References _subdomain_ids, libMesh::PointLocatorBase::get_close_to_point_tol(), get_point_locator(), init(), libMesh::PointLocatorBase::initialized(), libMesh::FunctionBase< Output >::initialized(), and set_point_locator_tolerance().
|
default |
Special functions.
|
default |
Destructor.
|
protected |
Helper function for finding a gradient as evaluated from a specific element.
Definition at line 457 of file mesh_function.C.
References _dof_map, _eqn_systems, _out_of_mesh_mode, _out_of_mesh_value, _system_vars, _vector, libMesh::TypeVector< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), libMesh::FEInterface::compute_data(), dim, libMesh::Elem::dim(), libMesh::DofMap::dof_indices(), libMesh::FEComputeData::dshape, libMesh::FEComputeData::enable_derivative(), libMesh::index_range(), libMesh::Elem::infinite(), libMesh::invalid_uint, libMesh::FEMap::inverse_map(), libMesh::libmesh_assert(), libMesh::FEComputeData::local_transform, libMesh::DenseVector< T >::size(), and libMesh::DofMap::variable_type().
Referenced by discontinuous_gradient(), and gradient().
|
protected |
Helper function that is called by MeshFunction::find_element() and MeshFunction::find_elements() to ensure that Elems found by the PointLocator are actually "evaluable" on the processor where they are found.
Definition at line 689 of file mesh_function.C.
References _dof_map, _vector, libMesh::Elem::find_point_neighbors(), libMesh::GHOSTED, libMesh::DofMap::is_evaluable(), libMesh::ParallelObject::processor_id(), libMesh::DofObject::processor_id(), libMesh::SERIAL, and libMesh::NumericVector< T >::type().
Referenced by find_element(), and find_elements().
|
overridevirtual |
Clears the function.
Reimplemented from libMesh::FunctionBase< Number >.
Definition at line 142 of file mesh_function.C.
References libMesh::FunctionBase< Number >::_initialized, libMesh::FunctionBase< Number >::_master, and _point_locator.
|
overridevirtual |
The new copy uses the original as a master function to enable simultaneous evaluations of the copies in different threads.
Implements libMesh::FunctionBase< Number >.
Definition at line 153 of file mesh_function.C.
|
inlineinherited |
Parallel::Communicator
object used by this mesh. Definition at line 97 of file parallel_object.h.
References libMesh::ParallelObject::_communicator.
Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::Partitioner::_find_global_index_by_pid_map(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< libMesh::Number >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< libMesh::Number >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_mult_add(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::DofMap::add_constraints_to_send_list(), add_cube_convex_hull_to_mesh(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::TransientRBConstruction::add_IC_to_RB_space(), libMesh::RBEIMEvaluation::add_interpolation_data(), libMesh::CondensedEigenSystem::add_matrices(), libMesh::EigenSystem::add_matrices(), libMesh::System::add_matrix(), libMesh::RBConstruction::add_scaled_matrix_and_vector(), libMesh::System::add_variable(), libMesh::System::add_variables(), libMesh::System::add_vector(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::TransientRBConstruction::allocate_data_structures(), libMesh::RBConstruction::allocate_data_structures(), libMesh::TransientRBConstruction::assemble_affine_expansion(), libMesh::AdvectionSystem::assemble_claw_rhs(), libMesh::FEMSystem::assemble_qoi(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::MeshCommunication::assign_global_indices(), libMesh::Partitioner::assign_partitioning(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::Partitioner::build_graph(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_data(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::RBConstruction::compute_Fq_representor_innerprods(), libMesh::RBConstruction::compute_max_error_bound(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::RBConstruction::compute_output_dual_innerprods(), libMesh::RBConstruction::compute_residual_dual_norm_slow(), libMesh::RBSCMConstruction::compute_SCM_bounds_on_training_set(), libMesh::DofMap::computed_sparsity_already(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::ContinuationSystem::ContinuationSystem(), libMesh::MeshBase::copy_constraint_rows(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::CondensedEigenSystem::copy_super_to_sub(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::DofMap::create_dof_constraints(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::PetscMatrix< T >::create_submatrix_nosort(), create_wrapped_function(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::RBEIMEvaluation::distribute_bfs(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::DTKSolutionTransfer::DTKSolutionTransfer(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_interiors(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_nodes(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_sides(), libMesh::TransientRBConstruction::enrich_RB_space(), libMesh::EpetraVector< T >::EpetraVector(), AssembleOptimization::equality_constraints(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_error_tolerance(), libMesh::MeshRefinement::flag_elements_by_mean_stddev(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::DofMap::gather_constraints(), libMesh::MeshfreeInterpolation::gather_remote_data(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::RBEIMEvaluation::get_eim_basis_function_node_value(), libMesh::RBEIMEvaluation::get_eim_basis_function_side_value(), libMesh::RBEIMEvaluation::get_eim_basis_function_value(), libMesh::MeshBase::get_info(), libMesh::System::get_info(), libMesh::DofMap::get_info(), libMesh::RBEIMEvaluation::get_interior_basis_functions_as_vecs(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::RBEIMConstruction::get_max_abs_value(), libMesh::RBEIMConstruction::get_node_max_abs_value(), libMesh::RBEIMEvaluation::get_parametrized_function_node_value(), libMesh::RBEIMEvaluation::get_parametrized_function_side_value(), libMesh::RBEIMEvaluation::get_parametrized_function_value(), libMesh::RBEIMConstruction::get_random_point(), AssembleOptimization::inequality_constraints(), AssembleOptimization::inequality_constraints_jacobian(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::StaticCondensation::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::AdvectionSystem::init_data(), libMesh::ClawSystem::init_data(), libMesh::PetscDMWrapper::init_petscdm(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::RBEIMConstruction::initialize_parametrized_functions_in_training_set(), libMesh::RBEIMConstruction::inner_product(), integrate_function(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_equal_connectivity(), libMesh::MeshTools::libmesh_assert_equal_points(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_constraint_rows(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_linesearch_shellfunc(), libMesh::libmesh_petsc_preconditioner_apply(), libMesh::libmesh_petsc_recalculate_monitor(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_interface(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_precheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::LinearImplicitSystem::LinearImplicitSystem(), main(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_bcids_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::TransientRBConstruction::mass_matrix_scaled_matvec(), libMesh::FEMSystem::mesh_position_set(), libMesh::TriangulatorInterface::MeshedHole::MeshedHole(), LinearElasticityWithContact::move_mesh(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::MeshTools::n_connected_components(), libMesh::DofMap::n_constrained_dofs(), libMesh::MeshBase::n_constraint_rows(), libMesh::DofMap::n_dofs(), libMesh::DofMap::n_dofs_per_processor(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), MixedOrderTest::n_neighbor_links(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::SparsityPattern::Build::n_nonzeros(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::RBEIMEvaluation::node_distribute_bfs(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::RBEIMConstruction::node_inner_product(), libMesh::MeshBase::operator==(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::BoundaryInfo::parallel_sync_node_ids(), libMesh::BoundaryInfo::parallel_sync_side_ids(), libMesh::MeshTools::paranoid_n_levels(), libMesh::Partitioner::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::MeshBase::print_constraint_rows(), libMesh::DofMap::print_dof_constraints(), libMesh::DofMap::process_mesh_constraint_rows(), libMesh::Partitioner::processor_pairs_to_interface_nodes(), libMesh::InterMeshProjection::project_system_vectors(), FEMParameters::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::EquationSystems::read(), libMesh::ExodusII_IO::read_header(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::RBEIMEvaluation::read_in_interior_basis_functions(), libMesh::RBEIMEvaluation::read_in_node_basis_functions(), libMesh::RBEIMEvaluation::read_in_side_basis_functions(), libMesh::RBEvaluation::read_in_vectors_from_multiple_files(), libMesh::System::read_legacy_data(), libMesh::TransientRBConstruction::read_riesz_representors_from_files(), libMesh::RBConstruction::read_riesz_representors_from_files(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::Nemesis_IO_Helper::read_var_names_impl(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::SimplexRefiner::refine_via_edges(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DistributedMesh::renumber_nodes_and_elements(), LinearElasticityWithContact::residual_and_jacobian(), OverlappingAlgebraicGhostingTest::run_ghosting_test(), OverlappingCouplingGhostingTest::run_sparsity_pattern_test(), scale_mesh_and_plot(), libMesh::DofMap::scatter_constraints(), libMesh::CheckpointIO::select_split_config(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::send_and_insert_dof_values(), libMesh::TransientRBConstruction::set_error_temporal_data(), libMesh::Partitioner::set_interface_node_processor_ids_BFS(), libMesh::Partitioner::set_interface_node_processor_ids_linear(), libMesh::Partitioner::set_interface_node_processor_ids_petscpartitioner(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Partitioner::set_parent_processor_ids(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::RBEIMEvaluation::side_distribute_bfs(), libMesh::RBEIMEvaluation::side_gather_bfs(), libMesh::RBEIMConstruction::side_inner_product(), libMesh::Partitioner::single_partition(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::ClawSystem::solve_conservation_law(), libMesh::split_mesh(), libMesh::RBEIMConstruction::store_eim_solutions_for_training_set(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), ConstraintOperatorTest::test1DCoarseningNewNodes(), ConstraintOperatorTest::test1DCoarseningOperator(), libMesh::MeshRefinement::test_level_one(), MeshfunctionDFEM::test_mesh_function_dfem(), MeshfunctionDFEM::test_mesh_function_dfem_grad(), MeshFunctionTest::test_p_level(), libMesh::MeshRefinement::test_unflagged(), DofMapTest::testBadElemFECombo(), SystemsTest::testBlockRestrictedVarNDofs(), BoundaryInfoTest::testBoundaryOnChildrenErrors(), ConstraintOperatorTest::testCoreform(), ConnectedComponentsTest::testEdge(), MeshInputTest::testExodusIGASidesets(), MeshTriangulationTest::testFoundCenters(), PointLocatorTest::testLocator(), BoundaryInfoTest::testMesh(), PointLocatorTest::testPlanar(), MeshTriangulationTest::testPoly2TriRefinementBase(), SystemsTest::testProjectCubeWithMeshFunction(), BoundaryInfoTest::testRenumber(), CheckpointIOTest::testSplitter(), MeshInputTest::testTetgenIO(), MeshTriangulationTest::testTriangulatorInterp(), MeshTriangulationTest::testTriangulatorMeshedHoles(), MeshTriangulationTest::testTriangulatorRoundHole(), libMesh::MeshTools::total_weight(), libMesh::RBConstruction::train_reduced_basis_with_POD(), libMesh::MeshFunctionSolutionTransfer::transfer(), libMesh::MeshfreeSolutionTransfer::transfer(), libMesh::Poly2TriTriangulator::triangulate(), libMesh::TransientRBConstruction::truth_assembly(), libMesh::RBConstruction::truth_assembly(), libMesh::MeshRefinement::uniformly_coarsen(), update_current_local_solution(), libMesh::TransientRBConstruction::update_RB_initial_condition_all_N(), libMesh::TransientRBConstruction::update_RB_system_matrices(), libMesh::RBConstruction::update_RB_system_matrices(), libMesh::TransientRBConstruction::update_residual_terms(), libMesh::RBConstruction::update_residual_terms(), libMesh::MeshTools::volume(), libMesh::STLIO::write(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::VTKIO::write_nodal_data(), libMesh::RBEIMEvaluation::write_out_interior_basis_functions(), libMesh::RBEIMEvaluation::write_out_node_basis_functions(), libMesh::RBEIMEvaluation::write_out_side_basis_functions(), libMesh::RBEvaluation::write_out_vectors(), libMesh::TransientRBConstruction::write_riesz_representors_to_files(), libMesh::RBConstruction::write_riesz_representors_to_files(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::RBDataSerialization::RBEvaluationSerialization::write_to_file(), libMesh::RBDataSerialization::TransientRBEvaluationSerialization::write_to_file(), libMesh::RBDataSerialization::RBEIMEvaluationSerialization::write_to_file(), and libMesh::RBDataSerialization::RBSCMEvaluationSerialization::write_to_file().
|
inlinevirtualinherited |
i
at coordinate p
and time time
.i+1
which will result in unexpected behaviour if operator() makes any access beyond that limit. Reimplemented in TripleFunction, SolutionFunction< dim >, PeriodicQuadFunction, SolutionFunction< dim >, SolutionFunction< dim >, SolutionFunction< dim >, SolutionFunction< dim >, SolutionFunction< dim >, and SolutionFunction< dim >.
Definition at line 232 of file function_base.h.
void libMesh::MeshFunction::disable_out_of_mesh_mode | ( | ) |
Disables out-of-mesh mode.
This is also the default.
Definition at line 756 of file mesh_function.C.
References _out_of_mesh_mode, _point_locator, libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
std::map< const Elem *, Gradient > libMesh::MeshFunction::discontinuous_gradient | ( | const Point & | p, |
const Real | time = 0. |
||
) |
p
and for time
. map is from element to Gradient and accounts for double defined values on faces if the gradient is discontinuous Definition at line 195 of file mesh_function.C.
References libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
Referenced by discontinuous_gradient().
void libMesh::MeshFunction::discontinuous_gradient | ( | const Point & | p, |
const Real | time, | ||
std::map< const Elem *, std::vector< Gradient >> & | output | ||
) |
Similar to gradient, but with the difference that multiple values on faces are explicitly permitted.
This is useful for evaluating gradients on faces where the values to the left and right are different.
Definition at line 419 of file mesh_function.C.
References _subdomain_ids, and discontinuous_gradient().
void libMesh::MeshFunction::discontinuous_gradient | ( | const Point & | p, |
const Real | time, | ||
std::map< const Elem *, std::vector< Gradient >> & | output, | ||
const std::set< subdomain_id_type > * | subdomain_ids | ||
) |
Similar to gradient, but with the difference that multiple values on faces are explicitly permitted.
This is useful for evaluating gradients on faces where the values to the left and right are different.
Definition at line 428 of file mesh_function.C.
References _gradient_on_elem(), _system_vars, find_elements(), libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
std::map< const Elem *, Number > libMesh::MeshFunction::discontinuous_value | ( | const Point & | p, |
const Real | time = 0. |
||
) |
p
and for time
.The std::map is from element to Number and accounts for doubly-defined values on faces if discontinuous variables are used.
Definition at line 170 of file mesh_function.C.
References libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
Referenced by discontinuous_value().
void libMesh::MeshFunction::discontinuous_value | ( | const Point & | p, |
const Real | time, | ||
std::map< const Elem *, DenseVector< Number >> & | output | ||
) |
Similar to operator() with the same parameter list, but with the difference that multiple values on faces are explicitly permitted.
This is useful for discontinuous shape functions that are evaluated on faces.
Definition at line 309 of file mesh_function.C.
References _subdomain_ids, and discontinuous_value().
void libMesh::MeshFunction::discontinuous_value | ( | const Point & | p, |
const Real | time, | ||
std::map< const Elem *, DenseVector< Number >> & | output, | ||
const std::set< subdomain_id_type > * | subdomain_ids | ||
) |
Similar to operator() with the same parameter list, but with the difference that multiple values on faces are explicitly permitted.
This is useful for discontinuous shape functions that are evaluated on faces.
Definition at line 318 of file mesh_function.C.
References _dof_map, _eqn_systems, _out_of_mesh_mode, _out_of_mesh_value, _system_vars, _vector, libMesh::FEInterface::compute_data(), dim, libMesh::DofMap::dof_indices(), find_elements(), libMesh::index_range(), libMesh::FunctionBase< Number >::initialized(), libMesh::invalid_uint, libMesh::FEMap::inverse_map(), libMesh::libmesh_assert(), libMesh::FEComputeData::shape, libMesh::DenseVector< T >::size(), value, and libMesh::DofMap::variable_type().
void libMesh::MeshFunction::enable_out_of_mesh_mode | ( | const DenseVector< Number > & | value | ) |
Enables out-of-mesh mode.
In this mode, if asked for a point that is not contained in any element, the MeshFunction
will return the given value
instead of crashing. This mode is off per default. If you use a master mesh function and you want to enable this mode, you will have to enable it for the master mesh function as well and for all mesh functions that have the same master mesh function. You may, however, specify different values.
Definition at line 741 of file mesh_function.C.
References _out_of_mesh_mode, _out_of_mesh_value, _point_locator, libMesh::FunctionBase< Number >::initialized(), libMesh::libmesh_assert(), and value.
Referenced by enable_out_of_mesh_mode().
void libMesh::MeshFunction::enable_out_of_mesh_mode | ( | const Number & | value | ) |
Enables out-of-mesh mode.
In this mode, if asked for a point that is not contained in any element, the MeshFunction
will return the given value
instead of crashing. This mode is off per default. If you use a master mesh function and you want to enable this mode, you will have to enable it for the master mesh function as well and for all mesh functions that have the same master mesh function. You may, however, specify different values.
Definition at line 749 of file mesh_function.C.
References enable_out_of_mesh_mode(), and value.
|
protected |
Helper function to reduce code duplication.
Definition at line 630 of file mesh_function.C.
References libMesh::FunctionBase< Number >::_master, _out_of_mesh_mode, and check_found_elem().
Referenced by gradient(), hessian(), and operator()().
|
protected |
p
.This is similar to find_element()
but covers cases where p
is on the boundary.
Definition at line 656 of file mesh_function.C.
References libMesh::FunctionBase< Number >::_master, _out_of_mesh_mode, and check_found_elem().
Referenced by discontinuous_gradient(), and discontinuous_value().
const PointLocatorBase & libMesh::MeshFunction::get_point_locator | ( | ) | const |
PointLocator
object, for use elsewhere.MeshFunction
object must be initialized before this is called. Definition at line 729 of file mesh_function.C.
References _point_locator, libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
Referenced by MeshFunction().
PointLocatorBase & libMesh::MeshFunction::get_point_locator | ( | ) |
Definition at line 735 of file mesh_function.C.
References _point_locator, libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
p
and for time
, which defaults to zero. Definition at line 185 of file mesh_function.C.
References libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
Referenced by libMesh::ExactErrorEstimator::find_squared_element_error(), gptr(), and gradient().
void libMesh::MeshFunction::gradient | ( | const Point & | p, |
const Real | time, | ||
std::vector< Gradient > & | output | ||
) |
Computes gradients at coordinate p
and for time time
, which defaults to zero, optionally restricting the point to the MeshFunction subdomain_ids.
Definition at line 394 of file mesh_function.C.
References _subdomain_ids, and gradient().
void libMesh::MeshFunction::gradient | ( | const Point & | p, |
const Real | time, | ||
std::vector< Gradient > & | output, | ||
const std::set< subdomain_id_type > * | subdomain_ids | ||
) |
Computes gradients at coordinate p
and for time time
, which defaults to zero, optionally restricting the point to the passed subdomain_ids, which parameter overrides the internal subdomain_ids.
Definition at line 403 of file mesh_function.C.
References _gradient_on_elem(), find_element(), libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
p
and for time
, which defaults to zero. Definition at line 211 of file mesh_function.C.
References libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
Referenced by libMesh::ExactErrorEstimator::find_squared_element_error(), and hessian().
void libMesh::MeshFunction::hessian | ( | const Point & | p, |
const Real | time, | ||
std::vector< Tensor > & | output | ||
) |
Computes gradients at coordinate p
and for time time
, which defaults to zero, optionally restricting the point to the MeshFunction subdomain_ids.
Definition at line 546 of file mesh_function.C.
References _subdomain_ids, and hessian().
void libMesh::MeshFunction::hessian | ( | const Point & | p, |
const Real | time, | ||
std::vector< Tensor > & | output, | ||
const std::set< subdomain_id_type > * | subdomain_ids | ||
) |
Computes gradients at coordinate p
and for time time
, which defaults to zero, optionally restricting the point to the passed subdomain_ids.
This is useful in cases where there are multiple dimensioned elements, for example.
Definition at line 555 of file mesh_function.C.
References _dof_map, _out_of_mesh_mode, _out_of_mesh_value, _system_vars, _vector, libMesh::TypeTensor< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), dim, libMesh::Elem::dim(), libMesh::DofMap::dof_indices(), find_element(), libMesh::index_range(), libMesh::Elem::infinite(), libMesh::FunctionBase< Number >::initialized(), libMesh::invalid_uint, libMesh::FEMap::inverse_map(), libMesh::libmesh_assert(), libMesh::DenseVector< T >::size(), and libMesh::DofMap::variable_type().
|
overridevirtual |
Override the FunctionBase::init() member function.
Reimplemented from libMesh::FunctionBase< Number >.
Definition at line 104 of file mesh_function.C.
References _eqn_systems, libMesh::FunctionBase< Number >::_initialized, _point_locator, _system_vars, libMesh::EquationSystems::get_mesh(), libMesh::libmesh_assert(), mesh, and libMesh::MeshBase::sub_point_locator().
Referenced by init(), main(), MeshFunction(), libMesh::InterMeshProjection::project_system_vectors(), MeshfunctionDFEM::test_mesh_function_dfem(), MeshfunctionDFEM::test_mesh_function_dfem_grad(), MeshFunctionTest::test_subdomain_id_sets(), SystemsTest::testProjectCubeWithMeshFunction(), and libMesh::MeshFunctionSolutionTransfer::transfer().
void libMesh::MeshFunction::init | ( | const Trees::BuildType | point_locator_build_type | ) |
The actual initialization process.
Takes an optional argument which specifies the method to use when building a PointLocator
Definition at line 128 of file mesh_function.C.
References init().
|
inlineinherited |
true
when this object is properly initialized and ready for use, false
otherwise. Definition at line 210 of file function_base.h.
Referenced by disable_out_of_mesh_mode(), discontinuous_gradient(), discontinuous_value(), enable_out_of_mesh_mode(), get_point_locator(), gradient(), hessian(), and operator()().
|
inlineinherited |
true
when the function this object represents is actually time-dependent, false
otherwise. Definition at line 224 of file function_base.h.
|
inlineinherited |
Definition at line 103 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, libMesh::libmesh_assert(), and TIMPI::Communicator::size().
Referenced by libMesh::Partitioner::_find_global_index_by_pid_map(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::DofMap::add_constraints_to_send_list(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::System::add_vector(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::FEMSystem::assembly(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::Partitioner::assign_partitioning(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::Partitioner::build_graph(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::DistributedMesh::clear_elems(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::copy_nodes_and_elements(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::Nemesis_IO::copy_scalar_solution(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_scalar_dofs(), libMesh::DistributedMesh::DistributedMesh(), libMesh::EnsightIO::EnsightIO(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::MeshBase::get_info(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::DofMap::n_dofs_per_processor(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::DofMap::prepare_send_list(), libMesh::MeshBase::print_constraint_rows(), libMesh::DofMap::print_dof_constraints(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::Partitioner::repartition(), OverlappingFunctorTest::run_partitioner_test(), libMesh::DofMap::scatter_constraints(), libMesh::DistributedMesh::set_next_unique_id(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), WriteVecAndScalar::setupTests(), libMesh::RBEIMEvaluation::side_gather_bfs(), DistributedMeshTest::testRemoteElemError(), CheckpointIOTest::testSplitter(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::VTKIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().
|
inlineinherited |
Evaluation function for time-independent vector-valued functions.
Sets output values in the passed-in output
DenseVector.
Definition at line 245 of file function_base.h.
p
and for time
, which defaults to zero. Implements libMesh::FunctionBase< Number >.
Definition at line 160 of file mesh_function.C.
References libMesh::FunctionBase< Number >::initialized(), and libMesh::libmesh_assert().
Referenced by operator()().
|
overridevirtual |
Computes values at coordinate p
and for time time
, which defaults to zero, optionally restricting the point to the MeshFunction subdomain_ids.
This is useful in cases where there are multiple dimensioned elements, for example.
Implements libMesh::FunctionBase< Number >.
Definition at line 222 of file mesh_function.C.
References _subdomain_ids, and operator()().
void libMesh::MeshFunction::operator() | ( | const Point & | p, |
const Real | time, | ||
DenseVector< Number > & | output, | ||
const std::set< subdomain_id_type > * | subdomain_ids | ||
) |
Computes values at coordinate p
and for time time
, restricting the point to the passed subdomain_ids, which parameter overrides the internal subdomain_ids.
Definition at line 229 of file mesh_function.C.
References _dof_map, _eqn_systems, _out_of_mesh_mode, _out_of_mesh_value, _system_vars, _vector, libMesh::FEInterface::compute_data(), dim, libMesh::Elem::dim(), libMesh::DofMap::dof_indices(), find_element(), libMesh::index_range(), libMesh::FunctionBase< Number >::initialized(), libMesh::invalid_uint, libMesh::FEMap::inverse_map(), libMesh::libmesh_assert(), libMesh::DenseVector< T >::resize(), libMesh::FEComputeData::shape, libMesh::DenseVector< T >::size(), value, and libMesh::DofMap::variable_type().
|
delete |
|
delete |
|
inlineinherited |
Definition at line 114 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and TIMPI::Communicator::rank().
Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::MeshTools::Modification::all_tri(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::FEMSystem::assembly(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::Partitioner::assign_partitioning(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Partitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), check_found_elem(), libMesh::DistributedMesh::clear(), libMesh::DistributedMesh::clear_elems(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::RBConstruction::compute_max_error_bound(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::Nemesis_IO::copy_scalar_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_scalar_dofs(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::find_dofs_to_send(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::RBEIMEvaluation::get_interior_basis_functions_as_vecs(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::DofMap::get_local_constraints(), libMesh::MeshBase::get_local_constraints(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::LaplaceMeshSmoother::init(), libMesh::StaticCondensation::init(), libMesh::SystemSubsetBySubdomain::init(), HeatSystem::init_data(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::TransientRBEvaluation::legacy_write_offline_data_to_files(), libMesh::RBSCMEvaluation::legacy_write_offline_data_to_files(), libMesh::RBEvaluation::legacy_write_offline_data_to_files(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), main(), libMesh::MeshRefinement::make_coarsening_compatible(), AugmentSparsityOnInterface::mesh_reinit(), libMesh::TriangulatorInterface::MeshedHole::MeshedHole(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::MeshTools::n_connected_components(), libMesh::MeshBase::n_constraint_rows(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::DistributedMesh::own_node(), libMesh::BoundaryInfo::parallel_sync_node_ids(), libMesh::BoundaryInfo::parallel_sync_side_ids(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::print_constraint_rows(), libMesh::DofMap::print_dof_constraints(), libMesh::DofMap::process_mesh_constraint_rows(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::EquationSystems::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::ExodusII_IO::read_header(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::DynaIO::read_mesh(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::TransientRBConstruction::read_riesz_representors_from_files(), libMesh::RBConstruction::read_riesz_representors_from_files(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::Nemesis_IO_Helper::read_var_names_impl(), libMesh::SimplexRefiner::refine_via_edges(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DistributedMesh::renumber_nodes_and_elements(), libMesh::DofMap::scatter_constraints(), libMesh::CheckpointIO::select_split_config(), libMesh::DistributedMesh::set_next_unique_id(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::RBEIMEvaluation::side_gather_bfs(), ExodusTest< elem_type >::test_read_gold(), ExodusTest< elem_type >::test_write(), MeshInputTest::testAbaqusRead(), MeshInputTest::testBadGmsh(), MeshInputTest::testCopyElementSolutionImpl(), MeshInputTest::testCopyElementVectorImpl(), MeshInputTest::testCopyNodalSolutionImpl(), DefaultCouplingTest::testCoupling(), PointNeighborCouplingTest::testCoupling(), MeshInputTest::testDynaFileMappings(), MeshInputTest::testDynaNoSplines(), MeshInputTest::testDynaReadElem(), MeshInputTest::testDynaReadPatch(), MeshInputTest::testExodusFileMappings(), MeshInputTest::testExodusIGASidesets(), MeshInputTest::testExodusWriteElementDataFromDiscontinuousNodalData(), MeshInputTest::testGoodGmsh(), MeshInputTest::testGoodSTL(), MeshInputTest::testGoodSTLBinary(), MeshInputTest::testLowOrderEdgeBlocks(), SystemsTest::testProjectMatrix3D(), BoundaryInfoTest::testShellFaceConstraints(), MeshInputTest::testSingleElementImpl(), WriteVecAndScalar::testSolution(), CheckpointIOTest::testSplitter(), MeshInputTest::testTetgenIO(), libMesh::MeshTools::total_weight(), libMesh::NetGenMeshInterface::triangulate(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::DTKAdapter::update_variable_values(), libMesh::MeshTools::volume(), libMesh::STLIO::write(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_element_values_element_major(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO_Helper::write_elemset_data(), libMesh::ExodusII_IO_Helper::write_elemsets(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::VTKIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_common(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::ExodusII_IO_Helper::write_nodeset_data(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::RBEIMEvaluation::write_out_interior_basis_functions(), libMesh::RBEIMEvaluation::write_out_node_basis_functions(), libMesh::RBEIMEvaluation::write_out_side_basis_functions(), write_output_solvedata(), libMesh::System::write_parallel_data(), libMesh::RBConstruction::write_riesz_representors_to_files(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::ExodusII_IO_Helper::write_sideset_data(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().
|
inlineinherited |
Function to set whether this is a time-dependent function or not.
This is intended to be only used by subclasses who cannot natively determine time-dependence. In such a case, this function should be used immediately following construction.
Definition at line 217 of file function_base.h.
void libMesh::MeshFunction::set_point_locator_tolerance | ( | Real | tol | ) |
We may want to specify a tolerance for the PointLocator to use, since in some cases the point we want to evaluate at might be slightly outside the mesh (due to numerical rounding issues, for example).
Definition at line 763 of file mesh_function.C.
References _point_locator.
Referenced by MeshFunction().
void libMesh::MeshFunction::set_subdomain_ids | ( | const std::set< subdomain_id_type > * | subdomain_ids | ) |
Choose a default list of subdomain ids to be searched for points.
If the provided list pointer is null or if no list has been provided, then all subdomain ids are searched. This list can be overridden on a per-evaluation basis by using the method overrides with a similar argument.
Definition at line 774 of file mesh_function.C.
References _subdomain_ids.
void libMesh::MeshFunction::unset_point_locator_tolerance | ( | ) |
Turn off the user-specified PointLocator tolerance.
Definition at line 769 of file mesh_function.C.
References _point_locator.
|
protectedinherited |
Definition at line 120 of file parallel_object.h.
Referenced by libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::StaticCondensation::close(), libMesh::ParallelObject::comm(), libMesh::CondensedEigenSystem::initialize_condensed_matrices(), libMesh::ParallelObject::n_processors(), libMesh::ParallelObject::operator=(), libMesh::ParallelObject::processor_id(), and libMesh::BoundaryInfo::regenerate_id_sets().
|
protected |
Need access to the DofMap
of the other system.
Definition at line 389 of file mesh_function.h.
Referenced by _gradient_on_elem(), check_found_elem(), discontinuous_value(), hessian(), and operator()().
|
protected |
The equation systems handler, from which the data are gathered.
Definition at line 378 of file mesh_function.h.
Referenced by _gradient_on_elem(), discontinuous_value(), init(), and operator()().
|
protectedinherited |
When init()
was called so that everything is ready for calls to operator()
(...), then this bool
is true.
Definition at line 184 of file function_base.h.
|
protectedinherited |
Cache whether or not this function is actually time-dependent.
Definition at line 189 of file function_base.h.
|
protectedinherited |
Const pointer to our master, initialized to nullptr
.
There may be cases where multiple functions are required, but to save memory, one master handles some centralized data.
Definition at line 178 of file function_base.h.
Referenced by clear(), find_element(), and find_elements().
|
protected |
true
if out-of-mesh mode is enabled.
See enable_out_of_mesh_mode()
for more details. Default is false
.
Definition at line 412 of file mesh_function.h.
Referenced by _gradient_on_elem(), disable_out_of_mesh_mode(), discontinuous_value(), enable_out_of_mesh_mode(), find_element(), find_elements(), hessian(), and operator()().
|
protected |
Value to return outside the mesh if out-of-mesh mode is enabled.
See enable_out_of_mesh_mode()
for more details.
Definition at line 418 of file mesh_function.h.
Referenced by _gradient_on_elem(), discontinuous_value(), enable_out_of_mesh_mode(), hessian(), and operator()().
|
protected |
A point locator is needed to locate the points in the mesh.
Definition at line 401 of file mesh_function.h.
Referenced by clear(), disable_out_of_mesh_mode(), enable_out_of_mesh_mode(), get_point_locator(), init(), set_point_locator_tolerance(), and unset_point_locator_tolerance().
|
protected |
A default set of subdomain ids in which to search for points.
Definition at line 406 of file mesh_function.h.
Referenced by discontinuous_gradient(), discontinuous_value(), gradient(), hessian(), MeshFunction(), operator()(), and set_subdomain_ids().
|
protected |
The indices of the variables within the other system for which data are to be gathered.
Definition at line 395 of file mesh_function.h.
Referenced by _gradient_on_elem(), discontinuous_gradient(), discontinuous_value(), hessian(), init(), and operator()().
|
protected |
A reference to the vector that holds the data that is to be interpolated.
Definition at line 384 of file mesh_function.h.
Referenced by _gradient_on_elem(), check_found_elem(), discontinuous_value(), hessian(), and operator()().