Skip to content

Commit

Permalink
Add VectorPostprocessors to get materials along a line closes idahola…
Browse files Browse the repository at this point in the history
…b#4462

This provides an abstract base class, LineMaterialSamplerBase, that
has the basic machinery for extracting the material properties for the
elements along a line, and two derived classes,
LineMaterialRealSampler and LineMaterialSymmTensorSampler, for
extracting data for specific material types.
  • Loading branch information
bwspenc committed Dec 22, 2014
1 parent 684f5bf commit 5d91453
Show file tree
Hide file tree
Showing 14 changed files with 699 additions and 0 deletions.
54 changes: 54 additions & 0 deletions framework/include/vectorpostprocessors/LineMaterialRealSampler.h
@@ -0,0 +1,54 @@
/****************************************************************/
/* DO NOT MODIFY THIS HEADER */
/* MOOSE - Multiphysics Object Oriented Simulation Environment */
/* */
/* (c) 2010 Battelle Energy Alliance, LLC */
/* ALL RIGHTS RESERVED */
/* */
/* Prepared by Battelle Energy Alliance, LLC */
/* Under Contract No. DE-AC07-05ID14517 */
/* With the U. S. Department of Energy */
/* */
/* See COPYRIGHT for full restrictions */
/****************************************************************/

#ifndef LINEMATERIALREALSAMPLER_H
#define LINEMATERIALREALSAMPLER_H

#include "LineMaterialSamplerBase.h"

//Forward Declarations
class LineMaterialRealSampler;

template<>
InputParameters validParams<LineMaterialRealSampler>();

/**
* This class samples Real material properties for the integration points
* in all elements that are intersected by a user-defined line.
*/
class LineMaterialRealSampler :
public LineMaterialSamplerBase<Real>
{
public:
/**
* Class constructor
* Sets up variables for output based on the properties to be output
* @param name The name of the class
* @param parameters The input parameters
*/
LineMaterialRealSampler(const std::string & name, InputParameters parameters);

virtual ~LineMaterialRealSampler() {}

/**
* Reduce the material property to a scalar for output
* In this case, the material property is a Real already, so just return it.
* @param property The material property
* @param curr_point The point corresponding to this material property
* @return A scalar value from this material property to be output
*/
virtual Real getScalarFromProperty(Real &property, const Point * /*curr_point*/);
};

#endif
196 changes: 196 additions & 0 deletions framework/include/vectorpostprocessors/LineMaterialSamplerBase.h
@@ -0,0 +1,196 @@
/****************************************************************/
/* DO NOT MODIFY THIS HEADER */
/* MOOSE - Multiphysics Object Oriented Simulation Environment */
/* */
/* (c) 2010 Battelle Energy Alliance, LLC */
/* ALL RIGHTS RESERVED */
/* */
/* Prepared by Battelle Energy Alliance, LLC */
/* Under Contract No. DE-AC07-05ID14517 */
/* With the U. S. Department of Energy */
/* */
/* See COPYRIGHT for full restrictions */
/****************************************************************/

#ifndef LINEMATERIALSAMPLERBASE_H
#define LINEMATERIALSAMPLERBASE_H

#include "GeneralVectorPostprocessor.h"
#include "RayTracing.h"
#include "SamplerBase.h"
#include "FEProblem.h"
#include "InputParameters.h"

//Forward Declarations
template<typename T>
class LineMaterialSamplerBase;

template<>
InputParameters validParams<LineMaterialSamplerBase<Real> >();

/**
* This is a base class for sampling material properties for the integration points
* in all elements that are intersected by a user-defined line. The positions of
* those points are output in x, y, z coordinates, as well as in terms of the projected
* positions of those points along the line. Derived classes can be created to sample
* arbitrary types of material properties.
*/
template<typename T>
class LineMaterialSamplerBase :
public GeneralVectorPostprocessor,
public SamplerBase
{
public:
/**
* Class constructor
* Sets up variables for output based on the properties to be output
* @param name The name of the class
* @param parameters The input parameters
*/
LineMaterialSamplerBase(const std::string & name, InputParameters parameters);

/**
* Class destructor
*/
virtual ~LineMaterialSamplerBase() {}

/**
* Initialize
* Calls through to base class's initialize()
*/
virtual void initialize();

/**
* Finds all elements along the user-defined line, loops through them, and samples their
* material properties.
*/
virtual void execute();

/**
* Finalize
* Calls through to base class's finalize()
*/
virtual void finalize();

/**
* Thread Join
* Calls through to base class's threadJoin()
* @param sb SamplerBase object to be joint into this object
*/
virtual void threadJoin(const SamplerBase & sb);

/**
* Reduce the material property to a scalar for output
* @param property The material property
* @param curr_point The point corresponding to this material property
* @return A scalar value from this material property to be output
*/
virtual Real getScalarFromProperty(T & property, const Point * curr_point) = 0;

protected:
/// The beginning of the line
Point _start;

/// The end of the line
Point _end;

/// The material properties to be output
std::vector<MaterialProperty<T> *> _material_properties;

/// The mesh
MooseMesh & _mesh;

/// The quadrature rule
QBase * & _qrule;

/// The quadrature points
const MooseArray<Point> & _q_point;
};

template <typename T>
LineMaterialSamplerBase<T>::LineMaterialSamplerBase(const std::string & name, InputParameters parameters) :
GeneralVectorPostprocessor(name, parameters),
SamplerBase(name, parameters, this, _communicator),
_start(getParam<Point>("start")),
_end(getParam<Point>("end")),
_mesh(_subproblem.mesh()),
_qrule(_subproblem.assembly(_tid).qRule()),
_q_point(_subproblem.assembly(_tid).qPoints())
{
std::vector<std::string> material_property_names = getParam<std::vector<std::string> >("property");
for (unsigned int i=0; i<material_property_names.size(); ++i)
{
if (!hasMaterialProperty<T>(material_property_names[i]))
mooseError("In LineMaterialSamplerBase material property: " + material_property_names[i] + " does not exist.");
_material_properties.push_back(&getMaterialProperty<T>(material_property_names[i]));
}

SamplerBase::setupVariables(material_property_names);
}

template <typename T>
void
LineMaterialSamplerBase<T>::initialize()
{
SamplerBase::initialize();
}

template <typename T>
void
LineMaterialSamplerBase<T>::execute()
{
std::vector<Elem *> intersected_elems;
Moose::elementsIntersectedByLine(_start, _end, _fe_problem.mesh(), intersected_elems);

const RealVectorValue line_vec = _end - _start;
const Real line_length(line_vec.size());
const RealVectorValue line_unit_vec = line_vec / line_length;
std::vector<Real> values(_material_properties.size());

for (unsigned int i=0; i<intersected_elems.size(); ++i)
{
const Elem * elem = intersected_elems[i];

if (elem->processor_id() != processor_id())
continue;

_subproblem.prepare(elem, _tid);
_subproblem.reinitElem(elem, _tid);
_fe_problem.reinitMaterials(elem->subdomain_id(), _tid);

for (unsigned int qp=0; qp<_qrule->n_points(); ++qp)
{
const RealVectorValue qp_pos(_q_point[qp]);

const RealVectorValue start_to_qp(qp_pos - _start);
const Real qp_proj_dist_along_line = start_to_qp * line_unit_vec;

if (qp_proj_dist_along_line < 0 || qp_proj_dist_along_line > line_length)
continue;

for (unsigned int j=0; j<_material_properties.size(); ++j)
{
values[j] = getScalarFromProperty((*_material_properties[j])[qp], &_q_point[qp]);
}
addSample(_q_point[qp], qp_proj_dist_along_line, values);
}
_fe_problem.swapBackMaterials(_tid);
}
}

template <typename T>
void
LineMaterialSamplerBase<T>::finalize()
{
SamplerBase::finalize();
}

template <typename T>
void
LineMaterialSamplerBase<T>::threadJoin(const SamplerBase & sb)
{
const LineMaterialSamplerBase<T> & lmsb = static_cast<const LineMaterialSamplerBase<T> &>(sb);
SamplerBase::threadJoin(sb);
}

#endif
2 changes: 2 additions & 0 deletions framework/src/base/Moose.C
Expand Up @@ -193,6 +193,7 @@
#include "LineValueSampler.h"
#include "VectorOfPostprocessors.h"
#include "ElementsAlongLine.h"
#include "LineMaterialRealSampler.h"

// user objects
#include "LayeredIntegral.h"
Expand Down Expand Up @@ -553,6 +554,7 @@ registerObjects(Factory & factory)
registerVectorPostprocessor(LineValueSampler);
registerVectorPostprocessor(VectorOfPostprocessors);
registerVectorPostprocessor(ElementsAlongLine);
registerVectorPostprocessor(LineMaterialRealSampler);

// user objects
registerUserObject(LayeredIntegral);
Expand Down
33 changes: 33 additions & 0 deletions framework/src/vectorpostprocessors/LineMaterialRealSampler.C
@@ -0,0 +1,33 @@
/****************************************************************/
/* DO NOT MODIFY THIS HEADER */
/* MOOSE - Multiphysics Object Oriented Simulation Environment */
/* */
/* (c) 2010 Battelle Energy Alliance, LLC */
/* ALL RIGHTS RESERVED */
/* */
/* Prepared by Battelle Energy Alliance, LLC */
/* Under Contract No. DE-AC07-05ID14517 */
/* With the U. S. Department of Energy */
/* */
/* See COPYRIGHT for full restrictions */
/****************************************************************/

#include "LineMaterialRealSampler.h"

template<>
InputParameters validParams<LineMaterialRealSampler>()
{
InputParameters params = validParams<LineMaterialSamplerBase<Real> >();
return params;
}

LineMaterialRealSampler::LineMaterialRealSampler(const std::string & name, InputParameters parameters) :
LineMaterialSamplerBase<Real>(name, parameters)
{
}

Real
LineMaterialRealSampler::getScalarFromProperty(Real & property, const Point * /*curr_point*/)
{
return property;
}
27 changes: 27 additions & 0 deletions framework/src/vectorpostprocessors/LineMaterialSamplerBase.C
@@ -0,0 +1,27 @@
/****************************************************************/
/* DO NOT MODIFY THIS HEADER */
/* MOOSE - Multiphysics Object Oriented Simulation Environment */
/* */
/* (c) 2010 Battelle Energy Alliance, LLC */
/* ALL RIGHTS RESERVED */
/* */
/* Prepared by Battelle Energy Alliance, LLC */
/* Under Contract No. DE-AC07-05ID14517 */
/* With the U. S. Department of Energy */
/* */
/* See COPYRIGHT for full restrictions */
/****************************************************************/

#include "LineMaterialSamplerBase.h"

template<>
InputParameters validParams<LineMaterialSamplerBase<Real> >()
{
InputParameters params = validParams<GeneralVectorPostprocessor>();
params += validParams<SamplerBase>();
params.addRequiredParam<Point>("start", "The beginning of the line");
params.addRequiredParam<Point>("end", "The end of the line");
params.addRequiredParam<std::vector<std::string> >("property", "Name of the material property to be output along a line");

return params;
}
@@ -0,0 +1,58 @@
/****************************************************************/
/* DO NOT MODIFY THIS HEADER */
/* MOOSE - Multiphysics Object Oriented Simulation Environment */
/* */
/* (c) 2010 Battelle Energy Alliance, LLC */
/* ALL RIGHTS RESERVED */
/* */
/* Prepared by Battelle Energy Alliance, LLC */
/* Under Contract No. DE-AC07-05ID14517 */
/* With the U. S. Department of Energy */
/* */
/* See COPYRIGHT for full restrictions */
/****************************************************************/

#ifndef LINEMATERIALTENSORSAMPLER_H
#define LINEMATERIALTENSORSAMPLER_H

#include "LineMaterialSamplerBase.h"
#include "MaterialTensorCalculator.h"

//Forward Declarations
class LineMaterialSymmTensorSampler;

template<>
InputParameters validParams<LineMaterialSymmTensorSampler>();

/**
* This class samples SymmTensor material properties for the integration points
* in all elements that are intersected by a user-defined line. It provides
* access to the full set of options for reducing the SymmTensor to a scalar.
*/
class LineMaterialSymmTensorSampler :
public LineMaterialSamplerBase<SymmTensor>,
public MaterialTensorCalculator
{
public:
/**
* Class constructor
* Sets up variables for output based on the properties to be output
* @param name The name of the class
* @param parameters The input parameters
*/
LineMaterialSymmTensorSampler(const std::string & name, InputParameters parameters);

virtual ~LineMaterialSymmTensorSampler() {}

/**
* Reduce the material property to a scalar for output
* Call through to getTensorQuantity to access the full set of options for reducing
* the SymmTensor to a scalar quantity.
* @param property The material property
* @param curr_point The point corresponding to this material property
* @return A scalar value from this material property to be output
*/
virtual Real getScalarFromProperty(SymmTensor &property, const Point * curr_point);
};

#endif

0 comments on commit 5d91453

Please sign in to comment.