www.mooseframework.org
LineMaterialSamplerBase.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 // MOOSE includes
14 #include "SamplerBase.h"
15 #include "BlockRestrictable.h"
16 #include "LineSegment.h"
17 #include "RayTracing.h" // Moose::elementsIntersectedByLine()
18 #include "Assembly.h" // Assembly::qRule()
19 #include "MooseMesh.h" // MooseMesh::getMesh()
20 #include "SwapBackSentinel.h"
21 #include "FEProblem.h"
22 
23 #include "libmesh/quadrature.h" // _qrule->n_points()
24 
25 // Forward Declarations
26 class MooseMesh;
27 template <typename T>
29 
30 template <>
31 InputParameters validParams<LineMaterialSamplerBase<Real>>();
32 
41 template <typename T>
43  public SamplerBase,
44  public BlockRestrictable
45 {
46 public:
53 
58  virtual void initialize() override;
59 
64  virtual void execute() override;
65 
70  virtual void finalize() override;
71 
78  virtual Real getScalarFromProperty(const T & property, const Point & curr_point) = 0;
79 
80 protected:
82  Point _start;
83 
85  Point _end;
86 
88  std::vector<const MaterialProperty<T> *> _material_properties;
89 
92 
94  const QBase * const & _qrule;
95 
98 };
99 
100 template <typename T>
102  : GeneralVectorPostprocessor(parameters),
103  SamplerBase(parameters, this, _communicator),
104  BlockRestrictable(this),
105  _start(getParam<Point>("start")),
106  _end(getParam<Point>("end")),
107  _mesh(_subproblem.mesh()),
108  _qrule(_subproblem.assembly(_tid).qRule()),
109  _q_point(_subproblem.assembly(_tid).qPoints())
110 {
111  std::vector<std::string> material_property_names = getParam<std::vector<std::string>>("property");
112  for (unsigned int i = 0; i < material_property_names.size(); ++i)
113  {
114  if (!hasMaterialProperty<T>(material_property_names[i]))
115  mooseError("In LineMaterialSamplerBase material property: " + material_property_names[i] +
116  " does not exist.");
117  _material_properties.push_back(&getMaterialProperty<T>(material_property_names[i]));
118  }
119 
120  SamplerBase::setupVariables(material_property_names);
121 }
122 
123 template <typename T>
124 void
126 {
128 }
129 
130 template <typename T>
131 void
133 {
134  std::vector<Elem *> intersected_elems;
135  std::vector<LineSegment> segments;
136 
137  std::unique_ptr<PointLocatorBase> pl = _mesh.getPointLocator();
138  Moose::elementsIntersectedByLine(_start, _end, _mesh, *pl, intersected_elems, segments);
139 
140  const RealVectorValue line_vec = _end - _start;
141  const Real line_length(line_vec.norm());
142  const RealVectorValue line_unit_vec = line_vec / line_length;
143  std::vector<Real> values(_material_properties.size());
144 
145  std::set<unsigned int> needed_mat_props;
146  const std::set<unsigned int> & mp_deps = getMatPropDependencies();
147  needed_mat_props.insert(mp_deps.begin(), mp_deps.end());
148  _fe_problem.setActiveMaterialProperties(needed_mat_props, _tid);
149 
150  for (const auto & elem : intersected_elems)
151  {
152  if (elem->processor_id() != processor_id())
153  continue;
154 
155  if (!hasBlocks(elem->subdomain_id()))
156  continue;
157 
158  _subproblem.prepare(elem, _tid);
159  _subproblem.reinitElem(elem, _tid);
160 
161  // Set up Sentinel class so that, even if reinitMaterials() throws, we
162  // still remember to swap back during stack unwinding.
163  SwapBackSentinel sentinel(_fe_problem, &FEProblem::swapBackMaterials, _tid);
164  _fe_problem.reinitMaterials(elem->subdomain_id(), _tid);
165 
166  for (unsigned int qp = 0; qp < _qrule->n_points(); ++qp)
167  {
168  const RealVectorValue qp_pos(_q_point[qp]);
169 
170  const RealVectorValue start_to_qp(qp_pos - _start);
171  const Real qp_proj_dist_along_line = start_to_qp * line_unit_vec;
172 
173  if (qp_proj_dist_along_line < 0 || qp_proj_dist_along_line > line_length)
174  continue;
175 
176  for (unsigned int j = 0; j < _material_properties.size(); ++j)
177  values[j] = getScalarFromProperty((*_material_properties[j])[qp], _q_point[qp]);
178 
179  addSample(_q_point[qp], qp_proj_dist_along_line, values);
180  }
181  }
182  _fe_problem.clearActiveMaterialProperties(_tid);
183 }
184 
185 template <typename T>
186 void
188 {
190 }
191 
Base class for VectorPostprocessors that need to do "sampling" of values in the domain.
Definition: SamplerBase.h:40
virtual void initialize()
Initialize the datastructures.
Definition: SamplerBase.C:76
VectorValue< Real > RealVectorValue
Definition: Assembly.h:31
This class is here to combine the VectorPostprocessor interface and the base class VectorPostprocesso...
LineMaterialSamplerBase(const InputParameters &parameters)
Class constructor Sets up variables for output based on the properties to be output.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
std::vector< const MaterialProperty< T > * > _material_properties
The material properties to be output.
virtual void initialize() override
Initialize Calls through to base class&#39;s initialize()
void mooseError(Args &&... args) const
Definition: MooseObject.h:147
Point _start
The beginning of the line.
Point _end
The end of the line.
virtual void finalize() override
Finalize Calls through to base class&#39;s finalize()
const QBase *const & _qrule
The quadrature rule.
void setupVariables(const std::vector< std::string > &variable_names)
You MUST call this in the constructor of the child class and pass down the name of the variables...
Definition: SamplerBase.C:52
virtual void execute() override
Finds all elements along the user-defined line, loops through them, and samples their material proper...
This is a base class for sampling material properties for the integration points in all elements that...
MooseMesh & _mesh
The mesh.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:65
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:74
virtual void swapBackMaterials(THREAD_ID tid)
virtual Real getScalarFromProperty(const T &property, const Point &curr_point)=0
Reduce the material property to a scalar for output.
virtual void finalize()
Finalize the values.
Definition: SamplerBase.C:92
const MooseArray< Point > & _q_point
The quadrature points.
An interface that restricts an object to subdomains via the &#39;blocks&#39; input parameter.
void elementsIntersectedByLine(const Point &p0, const Point &p1, const MeshBase &mesh, const PointLocatorBase &point_locator, std::vector< Elem *> &intersected_elems, std::vector< LineSegment > &segments)
Find all of the elements intersected by a line.
Definition: RayTracing.C:192
The "SwapBackSentinel" class&#39;s destructor guarantees that FEProblemBase::swapBackMaterials{Face,Neighbor}() is called even when an exception is thrown from FEProblemBase::reinitMaterials{Face,Neighbor}.