https://mooseframework.inl.gov
ParameterMesh.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "MooseTypes.h"
13 #include "KDTree.h"
14 
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"
24 #include <vector>
25 #include <functional>
26 
27 namespace libMesh
28 {
29 class System;
30 class DofMap;
31 }
32 
34 
44 {
45 public:
46  ParameterMesh(const libMesh::FEType & param_type,
47  const std::string & exodus_mesh,
48  const bool find_closest = false,
49  const unsigned int kdtree_candidates = 5);
50 
52  enum class RegularizationType
53  {
55  // Future regularization types can be added here:
56  // L1,
57  // H1,
58  // TV (Total Variation)
59  };
60 
64  dof_id_type size() const { return _param_dofs; }
72  void getIndexAndWeight(const Point & pt,
73  std::vector<dof_id_type> & dof_indices,
74  std::vector<Real> & weights) const;
82  void getIndexAndWeight(const Point & pt,
83  std::vector<dof_id_type> & dof_indices,
84  std::vector<RealGradient> & weights) const;
91  Real computeRegularizationObjective(const std::vector<Real> & parameter_values,
92  RegularizationType reg_type) const;
93 
100  std::vector<Real> computeRegularizationGradient(const std::vector<Real> & parameter_values,
101  RegularizationType reg_type) const;
102 
103 private:
110  template <typename T>
111  T computeRegularizationLoop(const std::vector<Real> & parameter_values,
112  RegularizationType reg_type) const;
113 
114 protected:
118  const bool _find_closest;
119  std::unique_ptr<libMesh::EquationSystems> _eq;
121  std::unique_ptr<libMesh::PointLocatorBase> _point_locator;
122  std::unique_ptr<libMesh::ExodusII_IO> _exodusII_io;
123 
125 
127  std::vector<Point> _mesh_nodes;
128  std::unique_ptr<KDTree> _node_kdtree;
129  std::unordered_map<dof_id_type, std::set<const libMesh::Elem *>> _node_to_elements;
130  unsigned int _kdtree_candidates;
131 
132 private:
138  Point projectToMesh(const Point & p) const;
139 
146  Point closestPoint(const Elem & elem, const Point & p) const;
147 
160  Real computeRegularizationQp(const std::vector<Real> & parameter_values,
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,
166  RegularizationType reg_type) const;
167 
180  void computeRegularizationGradientQp(const std::vector<Real> & parameter_values,
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,
186  RegularizationType reg_type,
187  std::vector<Real> & gradient) const;
188  // Cached values for gradient computations
189  const unsigned short int _param_var_id;
192 };
const unsigned short int _param_var_id
RealVectorValue RealGradient
libMesh::ReplicatedMesh _mesh
RegularizationType
Enumerations for regularization computations.
Definition: ParameterMesh.h:52
Utility class to use an Exodus mesh to define controllable parameters for optimization problems This ...
Definition: ParameterMesh.h:43
std::unique_ptr< libMesh::ExodusII_IO > _exodusII_io
ParameterMesh(const libMesh::FEType &param_type, const std::string &exodus_mesh, const bool find_closest=false, const unsigned int kdtree_candidates=5)
Definition: ParameterMesh.C:40
const bool _find_closest
Find closest projection points.
std::unique_ptr< libMesh::PointLocatorBase > _point_locator
T computeRegularizationLoop(const std::vector< Real > &parameter_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...
libMesh::System * _sys
Real computeRegularizationQp(const std::vector< Real > &parameter_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 > &parameter_values, RegularizationType reg_type) const
Computes regularization gradient for a given regularization type.
unsigned int _kdtree_candidates
Real computeRegularizationObjective(const std::vector< Real > &parameter_values, RegularizationType reg_type) const
Computes regularization objective value for a given regularization type.
std::unique_ptr< libMesh::EquationSystems > _eq
dof_id_type size() const
Definition: ParameterMesh.h:64
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 > &parameter_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
dof_id_type _param_dofs
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.
uint8_t dof_id_type
const libMesh::FEType _fe_type