15 #include "libmesh/id_types.h" 16 #include <unordered_map> 17 #include "libmesh/parallel.h" 18 #include "libmesh/fe_type.h" 19 #include "libmesh/point.h" 20 #include "libmesh/replicated_mesh.h" 21 #include "libmesh/equation_systems.h" 22 #include "libmesh/point_locator_base.h" 23 #include "libmesh/exodusII_io.h" 47 const std::string & exodus_mesh,
48 const bool find_closest =
false,
49 const unsigned int kdtree_candidates = 5);
73 std::vector<dof_id_type> & dof_indices,
74 std::vector<Real> & weights)
const;
83 std::vector<dof_id_type> & dof_indices,
84 std::vector<RealGradient> & weights)
const;
110 template <
typename T>
119 std::unique_ptr<libMesh::EquationSystems>
_eq;
146 Point
closestPoint(
const Elem & elem,
const Point & p)
const;
161 const std::vector<std::vector<Real>> & phi,
162 const std::vector<std::vector<RealGradient>> & dphi,
163 const unsigned int qp,
164 const std::vector<dof_id_type> & dof_indices,
165 const std::vector<Real> & JxW,
181 const std::vector<std::vector<Real>> & phi,
182 const std::vector<std::vector<RealGradient>> & dphi,
183 const unsigned int qp,
184 const std::vector<dof_id_type> & dof_indices,
185 const std::vector<Real> & JxW,
187 std::vector<Real> & gradient)
const;
const unsigned short int _param_var_id
RealVectorValue RealGradient
libMesh::ReplicatedMesh _mesh
RegularizationType
Enumerations for regularization computations.
Utility class to use an Exodus mesh to define controllable parameters for optimization problems This ...
std::unique_ptr< libMesh::ExodusII_IO > _exodusII_io
ParameterMesh(const libMesh::FEType ¶m_type, const std::string &exodus_mesh, const bool find_closest=false, const unsigned int kdtree_candidates=5)
const bool _find_closest
Find closest projection points.
std::unique_ptr< libMesh::PointLocatorBase > _point_locator
T computeRegularizationLoop(const std::vector< Real > ¶meter_values, RegularizationType reg_type) const
Template method containing the element loop for regularization computations.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
Real computeRegularizationQp(const std::vector< Real > ¶meter_values, const std::vector< std::vector< Real >> &phi, const std::vector< std::vector< RealGradient >> &dphi, const unsigned int qp, const std::vector< dof_id_type > &dof_indices, const std::vector< Real > &JxW, RegularizationType reg_type) const
Compute regularization objective for a single quadrature point This is the main function users should...
libMesh::Parallel::Communicator _communicator
std::unordered_map< dof_id_type, std::set< const libMesh::Elem * > > _node_to_elements
std::vector< Real > computeRegularizationGradient(const std::vector< Real > ¶meter_values, RegularizationType reg_type) const
Computes regularization gradient for a given regularization type.
unsigned int _kdtree_candidates
Real computeRegularizationObjective(const std::vector< Real > ¶meter_values, RegularizationType reg_type) const
Computes regularization objective value for a given regularization type.
std::unique_ptr< libMesh::EquationSystems > _eq
Point closestPoint(const Elem &elem, const Point &p) const
Find closest point on the element to the given point.
std::unique_ptr< KDTree > _node_kdtree
void computeRegularizationGradientQp(const std::vector< Real > ¶meter_values, const std::vector< std::vector< Real >> &phi, const std::vector< std::vector< RealGradient >> &dphi, const unsigned int qp, const std::vector< dof_id_type > &dof_indices, const std::vector< Real > &JxW, RegularizationType reg_type, std::vector< Real > &gradient) const
Compute regularization gradient for a single quadrature point This is the main function users should ...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const libMesh::DofMap * _dof_map
std::vector< Point > _mesh_nodes
Node-based KDTree optimization.
void getIndexAndWeight(const Point &pt, std::vector< dof_id_type > &dof_indices, std::vector< Real > &weights) const
Interpolate parameters onto the computational mesh getIndexAndWeight is only used by ParameterMeshFun...
Point projectToMesh(const Point &p) const
Returns the point on the parameter mesh that is projected from the test point.
const libMesh::FEType _fe_type