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 
26 namespace libMesh
27 {
28 class System;
29 }
30 
32 
43 {
44 public:
45  ParameterMesh(const libMesh::FEType & param_type,
46  const std::string & exodus_mesh,
47  const std::vector<std::string> & var_names = {},
48  const bool find_closest = false,
49  const unsigned int kdtree_candidates = 5);
50 
54  dof_id_type size() const { return _param_dofs; }
62  void getIndexAndWeight(const Point & pt,
63  std::vector<dof_id_type> & dof_indices,
64  std::vector<Real> & weights) const;
72  void getIndexAndWeight(const Point & pt,
73  std::vector<dof_id_type> & dof_indices,
74  std::vector<RealGradient> & weights) const;
82  std::vector<Real> getParameterValues(std::string var_name, unsigned int timestep) const;
83 
84 protected:
88  const bool _find_closest;
89  std::unique_ptr<libMesh::EquationSystems> _eq;
91  std::unique_ptr<libMesh::PointLocatorBase> _point_locator;
92  std::unique_ptr<libMesh::ExodusII_IO> _exodusII_io;
93 
95 
97  std::vector<Point> _mesh_nodes;
98  std::unique_ptr<KDTree> _node_kdtree;
99  std::unordered_map<dof_id_type, std::set<const libMesh::Elem *>> _node_to_elements;
100  unsigned int _kdtree_candidates;
101 
102 private:
108  Point projectToMesh(const Point & p) const;
109 
116  Point closestPoint(const Elem & elem, const Point & p) const;
117 };
RealVectorValue RealGradient
libMesh::ReplicatedMesh _mesh
Definition: ParameterMesh.h:86
Utility class to use an Exodus mesh to define controllable parameters for optimization problems This ...
Definition: ParameterMesh.h:42
std::unique_ptr< libMesh::ExodusII_IO > _exodusII_io
Definition: ParameterMesh.h:92
const bool _find_closest
Find closest projection points.
Definition: ParameterMesh.h:88
std::unique_ptr< libMesh::PointLocatorBase > _point_locator
Definition: ParameterMesh.h:91
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
libMesh::System * _sys
Definition: ParameterMesh.h:90
libMesh::Parallel::Communicator _communicator
Definition: ParameterMesh.h:85
std::unordered_map< dof_id_type, std::set< const libMesh::Elem * > > _node_to_elements
Definition: ParameterMesh.h:99
unsigned int _kdtree_candidates
std::unique_ptr< libMesh::EquationSystems > _eq
Definition: ParameterMesh.h:89
dof_id_type size() const
Definition: ParameterMesh.h:54
ParameterMesh(const libMesh::FEType &param_type, const std::string &exodus_mesh, const std::vector< std::string > &var_names={}, const bool find_closest=false, const unsigned int kdtree_candidates=5)
Definition: ParameterMesh.C:37
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
Definition: ParameterMesh.h:98
dof_id_type _param_dofs
Definition: ParameterMesh.h:94
std::vector< Point > _mesh_nodes
Node-based KDTree optimization.
Definition: ParameterMesh.h:97
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...
std::vector< Real > getParameterValues(std::string var_name, unsigned int timestep) const
Initializes parameter data and sets bounds in the main optmiization application getParameterValues is...
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