Postprocessor System

A postprocessor is an object that computes a single scalar (Real) value, such as a value sampled from the solution at a point in the domain, or an integral/average over some subdomain or boundary. This value may be used purely for output purposes, or it may be retrieved by other systems via the getPostprocessorValue method, which is available in most MOOSE objects. Furthermore, postprocessors are also functors, which allows them to be retrieved into various objects via the getFunctor<Real> method.

MOOSE includes a large number of postprocessors within the framework, the complete list is provided in Available Objects list section.

commentnote

The Reporter System is a newer, more flexible system for computing aggregate values. It is recommended that new objects for aggregate calculations use the Reporter system.

Example Input File

The following input file snippet demonstrates the use of the ElementExtremeValue to compute the minimum and maximum of the solution variable "u".

[Postprocessors]
  [max]
    type = ElementExtremeValue
    variable = u
  []
  [min]
    type = ElementExtremeValue
    variable = u
    value_type = min
  []
[]
(moose/test/tests/postprocessors/element_extreme_value/element_extreme_value.i)

This snippet is a part of a test that may be executed using the MOOSE test application as follows.


cd ~/projects/moose/test
make -j8
cd tests/postprocessors/element_extreme_value
~/projects/moose/test/moose_test-opt -i element_extreme_value.i

The data from this calculation is reported in the terminal output by default and if Exodus output is enabled the values will automatically be included in the output file. It is also possible to export the data to a comma separated value (csv) file by enabling the CSV object within the Outputs block.


Postprocessor Values:
+----------------+----------------+----------------+
| time           | max            | min            |
+----------------+----------------+----------------+
|   0.000000e+00 |   0.000000e+00 |   0.000000e+00 |
|   1.000000e+00 |   9.788675e-01 |   2.113249e-02 |
+----------------+----------------+----------------+

Coupling Example Code

The values computed within a Postprocessor object may be used within other objects that inherit from the PostprocessorInterface, which is nearly every system within MOOSE. For example, the PostprocessorNeumannBC object allows for a Neumann boundary condition to be set to a value computed from a postprocessor; this object will be used as example to demonstrate how coupling is performed.

To understand how the coupling is coded it is beneficial to first see how the coupling is defined via the input file. The following input file snippet shows that a PointValue postprocessor is created and named "right_pp" and the PostprocessorNeumannBC uses this value to set the boundary condition.

[Postprocessors]
  [right_pp]
    type = PointValue
    point = '0.5 0.5 0'
    variable = aux
    execute_on = 'initial'
  []
[]

[BCs]
  [left]
    type = DirichletBC
    variable = u
    boundary = left
    value = 0
  []
  [right]
    type = PostprocessorNeumannBC
    variable = u
    boundary = right
    postprocessor = right_pp
  []
[]
(moose/test/tests/bcs/pp_neumann/pp_neumann.i)

This first step of coding this type of coupling begins by adding the necessary input file syntax to the object that requires a postprocessor value, PostprocessorNeumannBC in this example. In all MOOSE objects input file syntax is governed by the validParams function of an object. To add the ability to couple a postprocessor, simply add a new parameter using the PostprocessorName type, as shown below. Notice, that the add parameters call includes a default value that makes the use of the postprocessor optional.


#include "PostprocessorNeumannBC.h"

registerMooseObject("MooseApp", PostprocessorNeumannBC);

InputParameters
(moose/framework/src/bcs/PostprocessorNeumannBC.C)

The actual postprocessor value must be assigned to a member variable of the class, thus in the header a member variable must be created, which should always be a constant reference to a PostprocessorValue type. Since this is a reference it must be initialized, this occurs in the source file by calling the getPostprocessorValue method and providing the name used in the validParams function. The following snippets show declaration of the reference in the header and the initialization of this reference in the source file. The _value member variable is then available for use anywhere inside the object, for the case of the boundary condition it is utilized in the computation of the residual.

  const PostprocessorValue & _value;
(moose/framework/include/bcs/PostprocessorNeumannBC.h)
PostprocessorNeumannBC::validParams()
{
  InputParameters params = IntegratedBC::validParams();
  params.addClassDescription(
      "Neumann boundary condition with value prescribed by a Postprocessor value.");
  params.addParam<PostprocessorName>(
      "postprocessor", 0.0, "The postprocessor to use for value of the gradient on the boundary.");
  return params;
}
(moose/framework/src/bcs/PostprocessorNeumannBC.C)

Coupling to other values

Just as Postprocessor values can be used in other objects, Postprocessors themselves can couple to Functions and Scalar Variables. See the following example that couples a scalar variable into a Postprocessor:

[Postprocessors]
  [./totalFlux]
    type = ScalarCoupledPostprocessor
    variable = u
    coupled_scalar = scalar_variable
    boundary = left
  [../]
[]
(moose/test/tests/postprocessors/scalar_coupled_postprocessor/scalar_coupled_postprocessor_test.i)

Creating a Postprocessor Object

In general, every Postprocessor object has two methods that must be defined "execute" and "getValue".

First, consider the execute method. This method is called by MOOSE at different time depending on the type of postprocessor object. Therefore, when creating a Postprocessor object the new object should inherit from one of the following C++ classes:

  • GeneralPostprocessor: "execute" is called once on each execution flag.

  • NodalPostprocessor: "execute" is called for each node within the mesh on each execution flag.

  • ElementalPostprocessor: "execute" is called for each element within the mesh on each execution flag.

  • InternalSidePostprocessor: "execute" is called for each side, that is not on a boundary, within the mesh on each execution flag.

  • SidePostprocessor: "execute" is called for each side, that is on a boundary, within the mesh on each execution flag.

The use of execution flags is discussed in the Execute On section.

The getValue method is responsible for returning the value of the postprocessor object, this value is what is used by all objects that are coupled to the postprocessor. In some cases the necessary communication is performed within this method, but in general this following is preferred.

Parallel Considerations

When creating a postprocessor it is often necessary to perform some parallel communication to ensure that the value being computed is correct across processes and threads. Three additional methods exists for making this process as simple as possible.

  • initialize: This is called prior to the execution of the postprocessor and should be used to setup the object to be in a known state. It is important to point out that execution in this context includes all calls to the execute method. For example, for a NodalPostprocessor object the initialize method is called and then the execute method is called for all nodes.

  • finalize: This is called after the execution of the postprocessor and is intended to perform communication to prepare the object for the call to the getValue method.

  • threadJoin: This is called after the execution of the postprocessor and is intended to perform aggregation for shared memory parallelism.

To understand the use of these methods the AverageNodalVariableValue postprocessor shall be used as an example. As the name suggests this postprocessor computes the average of the value of a variable at the nodes. To perform this calculation the variable values from each node are summed as is the number of values within the execute method. Then the getValue method returns the average by returning the sum divided by the count. The following snippet shows the these two methods: the _u[_qp] is the value of the variable at the current node that comes from a shared base class and _sum and _n are a member variables defined within class for performing the calculation.

void
AverageNodalVariableValue::execute()
{
  _sum += _u[_qp];
  _n++;
}

Real
AverageNodalVariableValue::getValue() const
{
  return _sum / _n;
}
(moose/framework/src/postprocessors/AverageNodalVariableValue.C)

In parallel, the calls to the execute method occur on each process or thread on a subset of the domain, in this case nodes. Therefore, the computed values must be combined to get the actual summations required to compute the average value. The first step is to setup the state of this calculation within the initialize method, which in this example sets the _sum and _n member variables to zero.

void
AverageNodalVariableValue::initialize()
{
  _sum = 0;
  _n = 0;
}
(moose/framework/src/postprocessors/AverageNodalVariableValue.C)

After the aforementioned execute method is called for each node the computed values for _sum and _n must be aggregated from across processes to the root processes. For this problem a gather operation is required to collect the values computed on all processes to the root process. This is accomplished via the gatherSum method.

void
AverageNodalVariableValue::finalize()
{
  gatherSum(_sum);
  gatherSum(_n);
}
(moose/framework/src/postprocessors/AverageNodalVariableValue.C)

Of course, the type of communication necessary depends on the calculation being performed. The UserObject base class includes helper methods for common parallel communications functions.

The initialize and finalize methods are utilized to aggregate for message passing (MPI) based parallelism. For shared memory parallelism the threadJoin method is used. This method is called, like finalize, after execution is complete and includes a single argument. This argument is a reference to a UserObject, which is a base class of Postprocessor objects. The purpose of this method is to enable the aggregation for the Postprocessor objects that were executed on other threads to the object on the root thread. For the AverageNodalVariableValue postprocessor the values for _sum and _n on the root process object are updated to include the these same values from the other threads.

void
AverageNodalVariableValue::threadJoin(const UserObject & y)
{
  const auto & pps = static_cast<const AverageNodalVariableValue &>(y);
  _sum += pps._sum;
  _n += pps._n;
}
(moose/framework/src/postprocessors/AverageNodalVariableValue.C)

Execute On...

Postprocessor objects inherit from the SetupInterface (execute_on) that allows the objects to execute and varying and multiple times during a simulation, such as during initialization and at the end of each time step. Refer to the SetupInterface (execute_on) for additional information.

Using Old and Older values

MOOSE maintains previously computed values in the postprocessor system for using lagged information in a computation. Both the previous time step's value and the value computed two steps back may be retrieved. One reason you might use older values is to break cyclic dependencies. MOOSE does not consider a dependence on an old value when considering the order of evaluation among objects with dependencies.

Available Objects

  • Moose App
  • ADElementAverageMaterialPropertyComputes the average of a material property over a volume.
  • ADElementExtremeFunctorValueFinds either the min or max elemental value of a variable over the domain.
  • ADElementExtremeMaterialPropertyDetermines the minimum or maximum of a material property over a volume.
  • ADElementIntegralFunctorPostprocessorComputes a volume integral of the specified functor
  • ADElementIntegralMaterialPropertyCompute the integral of the material property over the domain
  • ADElementL2FunctorErrorComputes L2 error between an 'approximate' functor and an 'exact' functor
  • ADInterfaceDiffusiveFluxAverageComputes the diffusive flux on the interface.
  • ADInterfaceDiffusiveFluxIntegralComputes the diffusive flux on the interface.
  • ADSideAdvectiveFluxIntegralComputes the volumetric advected quantity through a sideset.
  • ADSideAverageMaterialPropertyComputes the average of a material property over a side set.
  • ADSideDiffusiveFluxAverageComputes the integral of the diffusive flux over the specified boundary
  • ADSideDiffusiveFluxIntegralComputes the integral of the diffusive flux over the specified boundary
  • ADSideFluxAverageComputes the integral of the diffusive flux over the specified boundary
  • ADSideFluxIntegralComputes the integral of the diffusive flux over the specified boundary
  • ADSideIntegralFunctorPostprocessorComputes a surface integral of the specified functor, using the single-sided face argument, which usually means that the functor will be evaluated from a single side of the surface, not interpolated between both sides.
  • ADSideIntegralMaterialPropertyCompute the integral of a scalar material property component over the domain.
  • ADSideVectorDiffusivityFluxIntegralComputes the integral of the diffusive flux over the specified boundary
  • AreaPostprocessorComputes the "area" or dimension - 1 "volume" of a given boundary or boundaries in your mesh.
  • AverageElementSizeComputes the average element size.
  • AverageNodalVariableValueComputes the average value of a field by sampling all nodal solutions on the domain or within a subdomain
  • AxisymmetricCenterlineAverageValueComputes the average value of a variable on a sideset located along the centerline of an axisymmetric model.
  • ChangeOverFixedPointPostprocessorComputes the change or relative change in a post-processor value over a single or multiple fixed point iterations
  • ChangeOverTimePostprocessorComputes the change or relative change in a post-processor value over a timestep or the entire transient
  • ChangeOverTimestepPostprocessorComputes the change or relative change in a post-processor value over a timestep or the entire transient
  • ConstantPostprocessorPostprocessor that holds a constant value
  • CumulativeValuePostprocessorCreates a cumulative sum of a Postprocessor value with time.
  • DifferencePostprocessorComputes the difference between two postprocessors
  • ElementArrayL2NormEvaluates L2-norm of a component of an array variable
  • ElementAverageMaterialPropertyComputes the average of a material property over a volume.
  • ElementAverageSecondTimeDerivativeComputes the element averaged second derivative of variable
  • ElementAverageTimeDerivativeComputes a volume integral of the time derivative of a given variable
  • ElementAverageValueComputes the volumetric average of a variable
  • ElementExtremeFunctorValueFinds either the min or max elemental value of a variable over the domain.
  • ElementExtremeMaterialPropertyDetermines the minimum or maximum of a material property over a volume.
  • ElementExtremeValueFinds either the min or max elemental value of a variable over the domain.
  • ElementH1ErrorComputes the H1 error between a variable and a function
  • ElementH1SemiErrorReturns the gradient difference norm part of the H1 error
  • ElementHCurlErrorReturns the H(curl)-norm of the difference between a pair of computed and analytical vector-valued solutions.
  • ElementHCurlSemiErrorReturns the H(curl)-seminorm of the difference between a pair of computed and analytical vector-valued solutions.
  • ElementHDivErrorReturns the H(div)-norm of the difference between a pair of computed and analytical vector-valued solutions.
  • ElementHDivSemiErrorReturns the H(div)-seminorm of the difference between a pair of computed and analytical vector-valued solutions.
  • ElementIntegralArrayVariablePostprocessorIntegral of one component of an array variable.
  • ElementIntegralFunctorPostprocessorComputes a volume integral of the specified functor
  • ElementIntegralMaterialPropertyCompute the integral of the material property over the domain
  • ElementIntegralVariablePostprocessorComputes a volume integral of the specified variable
  • ElementL1ErrorComputes L1 error between an elemental field variable and an analytical function.
  • ElementL2DifferenceComputes the element-wise L2 difference between the current variable and a coupled variable.
  • ElementL2ErrorComputes L2 error between a field variable and an analytical function
  • ElementL2FunctorErrorComputes L2 error between an 'approximate' functor and an 'exact' functor
  • ElementL2NormComputes a volume integral of the specified variable
  • ElementSidesL2NormComputes the L2 norm of a variable over element sides.
  • ElementVectorL2ErrorReturns the L2-norm of the difference between a pair of computed and analytical vector-valued solutions.
  • ElementW1pErrorComputes the W1p norm of the difference between a variable and an analytic solution, as a function
  • ElementalVariableValueOutputs an elemental variable value at a particular location
  • EmptyPostprocessorA postprocessor object that returns a value of zero.
  • FindValueOnLineFind a specific target value along a sampling line. The variable values along the line should change monotonically. The target value is searched using a bisection algorithm.
  • FunctionElementAverageComputes the average of a function over a volume.
  • FunctionElementIntegralIntegrates a function over elements
  • FunctionSideAverageComputes the average of a function over a boundary.
  • FunctionSideIntegralComputes the integral of a function over a boundary.
  • FunctionValuePostprocessorComputes the value of a supplied function at a single point (scalable)
  • GreaterThanLessThanPostprocessorCount number of DOFs of a non-linear variable that are greater than or less than a given threshold
  • InterfaceAverageVariableValuePostprocessorComputes the average value of a variable on an interface. Note that this cannot be used on the centerline of an axisymmetric model.
  • InterfaceDiffusiveFluxAverageComputes the diffusive flux on the interface.
  • InterfaceDiffusiveFluxIntegralComputes the diffusive flux on the interface.
  • InterfaceIntegralVariableValuePostprocessorAdd access to variables and their gradient on an interface.
  • InternalSideIntegralVariablePostprocessorComputes an integral on internal sides of the specified variable
  • IsMatrixSymmetricReport whether the system matrix is symmetric or not.
  • LinearCombinationPostprocessorComputes a linear combination between an arbitrary number of post-processors
  • MemoryUsageMemory usage statistics for the running simulation.
  • NearestNodeNumberOutputs the nearest node number to a point
  • NodalExtremeValueFinds either the min or max elemental value of a variable over the domain.
  • NodalL2ErrorThe L2-norm of the difference between a variable and a function computed at nodes.
  • NodalL2NormComputes the nodal L2-norm of the coupled variable, which is defined by summing the square of its value at every node and taking the square root.
  • NodalMaxValueComputes the maximum (over all the nodal values) of a variable.
  • NodalMaxValueIdFinds the node id with the maximum nodal value across all postprocessors.
  • NodalProxyMaxValueFinds the node id with the maximum nodal value across all postprocessors.
  • NodalSumComputes the sum of all of the nodal values of the specified variable. Note: This object sets the default "unique_node_execute" flag to true to avoid double counting nodes between shared blocks.
  • NodalVariableValueOutputs values of a nodal variable at a particular location
  • NumDOFsReturn the number of Degrees of freedom from either the NL, Aux or both systems.
  • NumElementsReturn the number of active or total elements in the simulation.
  • NumElemsReturn the number of active or total elements in the simulation.
  • NumFailedTimeStepsCollects the number of failed time steps from the time stepper.
  • NumFixedPointIterationsReturns the number of fixed point iterations taken by the executioner.
  • NumLinearIterationsCompute the number of linear iterations.
  • NumMeshDivisionsReturn the number of divisions/regions from a MeshDivision object.
  • NumNodesReturns the total number of nodes in a simulation (works with DistributedMesh)
  • NumNonlinearIterationsOutputs the number of nonlinear iterations
  • NumPicardIterationsReturns the number of fixed point iterations taken by the executioner.
  • NumPositionsReturn the number of Positions from a Positions object.
  • NumRelationshipManagersReturn the number of relationship managers active.
  • NumResidualEvaluationsReturns the total number of residual evaluations performed.
  • NumVarsReturn the number of variables from either the NL, Aux, or both systems.
  • ParsedPostprocessorComputes a parsed expression with post-processors
  • PercentChangePostprocessorComputes the percent change of a postprocessor value compared to the value at the previous timestep.
  • PerfGraphDataRetrieves performance information about a section from the PerfGraph.
  • PointValueCompute the value of a variable at a specified location
  • PostprocessorComparisonCompares two post-processors and produces a boolean value
  • PseudoTimestepComputes pseudo-time steps for obtaining steady-state solutions through a pseudo transient process.
  • ReceiverReports the value stored in this processor, which is usually filled in by another object. The Receiver does not compute its own value.
  • RelativeDifferencePostprocessorComputes the absolute value of the relative difference between 2 post-processor values.
  • RelativeSolutionDifferenceNormComputes the relative norm of the solution difference of two consecutive time steps.
  • ResidualReport the non-linear residual.
  • ScalarL2ErrorCompute L2 error of a scalar variable using analytic function.
  • ScalarVariableReturns the value of a scalar variable as a postprocessor value.
  • ScalePostprocessorScales a post-processor by a value
  • SideAdvectiveFluxIntegralComputes the volumetric advected quantity through a sideset.
  • SideAverageFunctorPostprocessorComputes the average of a functor over a side set.
  • SideAverageMaterialPropertyComputes the average of a material property over a side set.
  • SideAverageValueComputes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.
  • SideDiffusiveFluxAverageComputes the integral of the diffusive flux over the specified boundary
  • SideDiffusiveFluxIntegralComputes the integral of the diffusive flux over the specified boundary
  • SideExtremeValueFinds either the min or max variable value of a variable over a boundary.
  • SideFluxAverageComputes the integral of the diffusive flux over the specified boundary
  • SideFluxIntegralComputes the integral of the diffusive flux over the specified boundary
  • SideIntegralFunctorPostprocessorComputes a surface integral of the specified functor, using the single-sided face argument, which usually means that the functor will be evaluated from a single side of the surface, not interpolated between both sides.
  • SideIntegralMaterialPropertyCompute the integral of a scalar material property component over the domain.
  • SideIntegralVariablePostprocessorComputes a surface integral of the specified variable
  • SideVectorDiffusivityFluxIntegralComputes the integral of the diffusive flux over the specified boundary
  • TagVectorSumComputes the sum of components of the requested tagged vector
  • TimeExtremeValueA postprocessor for reporting the extreme value of another postprocessor over time.
  • TimeIntegratedPostprocessorIntegrate a Postprocessor value over time using trapezoidal rule.
  • TimePostprocessorReports the current time
  • TimestepSizeReports the timestep size
  • TotalVariableValueIntegrate a Postprocessor value over time using trapezoidal rule.
  • VariableInnerProductComputes a volume integral of the specified variable
  • VariableResidualComputes the L2 norm of the residual of a single variable in the solution vector.
  • VectorPostprocessorComparisonCompares two vector post-processors of equal size and produces a boolean value
  • VectorPostprocessorComponentReturns the value of the specified component of a VectorPostprocessor
  • VectorPostprocessorReductionValueTakes a VectorPostprocessor and performs a reduction operation on it (max, min, sum, average) and stores as postprocessor.
  • VolumePostprocessorComputes the volume of a specified block
  • Phase Field App
  • AverageGrainVolumeCalculate average grain area in a polycrystal
  • DiscreteNucleationDataOutput diagnostic data on a DiscreteNucleationInserter
  • DiscreteNucleationTimeStepReturn a time step limit for nucleation event to be used by IterationAdaptiveDT
  • FauxGrainTrackerFake grain tracker object for cases where the number of grains is equal to the number of order parameters.
  • FauxPolycrystalVoronoiRandom Voronoi tessellation polycrystal when the number of order parameters equal to the number of grains
  • FeatureFloodCountThe object is able to find and count "connected components" in any solution field or number of solution fields. A primary example would be to count "bubbles".
  • FeatureVolumeFraction
  • GrainBoundaryAreaCalculate total grain boundary length in 2D and area in 3D
  • GrainTrackerGrain Tracker object for running reduced order parameter simulations without grain coalescence.
  • GrainTrackerElasticityGrain Tracker object for running reduced order parameter simulations without grain coalescence.
  • PFCElementEnergyIntegral
  • PolycrystalCirclesPolycrystal circles generated from a vector input or read from a file
  • PolycrystalEBSDObject for setting up a polycrystal structure from an EBSD Datafile
  • PolycrystalHexPerturbed hexagonal polycrystal
  • PolycrystalVoronoiRandom Voronoi tessellation polycrystal (used by PolycrystalVoronoiAction)
  • WeightedVariableAverageAverage a variable value using a weight mask given by a material property.
  • Ray Tracing App
  • RayDataValueObtains a value from the data or aux data of a Ray after tracing has been completed.
  • RayIntegralValueObtains the integrated value accumulated into a Ray from an IntegralRayKernel-derived class.
  • RayTracingStudyResultGets a result from a RayTracingStudy.
  • raccoon App
  • ExternalWorkThis class computes the total external work. The power expenditure (rate of external work) is defined as . The power expenditure is integrated in time to get the total work.
  • KineticEnergyThis class computes the total kinetic energy of the form .
  • LargeDeformationJIntegralThis class computes the J integral for a phase-field model of fracture
  • PhaseFieldJIntegralCompute the J integral for a phase-field model of fracture
  • SolutionChangeNormThis class computes the solution change (L2) norm of selected variables.
  • Heat Transfer App
  • ADConvectiveHeatTransferSideIntegralComputes the total convective heat transfer across a boundary.
  • ConvectiveHeatTransferSideIntegralComputes the total convective heat transfer across a boundary.
  • ExposedSideAverageValueComputes the average value of a variable on the exposed portion of a sideset. Note that this cannot be used on the centerline of an axisymmetric model.
  • GrayLambertSurfaceRadiationPPThis postprocessor allows to extract radiosity, heat flux density, and temperature from the GrayLambertSurfaceRadiationBase object.
  • HomogenizedThermalConductivityPostprocessor for asymptotic expansion homogenization for thermal conductivity
  • ThermalConductivityComputes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.
  • ViewFactorPPThis postprocessor allows to extract view factors from ViewFactor userobjects.
  • Solid Mechanics App
  • ADMassComputes the mass of the solid as the integral of the density material property
  • ADMaterialTensorAverageComputes the average of a RankTwoTensor component over a volume.
  • ADMaterialTensorIntegralThis postprocessor computes an element integral of a component of a material tensor as specified by the user-supplied indices
  • ADSidesetReactionComputes the integrated reaction force in a user-specified direction on a sideset from the surface traction
  • AsymptoticExpansionHomogenizationElasticConstantsPostprocessor for asymptotic expansion homogenization for elasticity
  • CavityPressurePostprocessorInterfaces with the CavityPressureUserObject to store the initial number of moles of a gas contained within an internal volume.
  • CrackFrontDataDetermines which nodes are along the crack front
  • CriticalTimeStepComputes and reports the critical time step for the explicit solver.
  • MassComputes the mass of the solid as the integral of the density material property
  • MaterialTensorAverageComputes the average of a RankTwoTensor component over a volume.
  • MaterialTensorIntegralThis postprocessor computes an element integral of a component of a material tensor as specified by the user-supplied indices
  • MaterialTimeStepPostprocessorThis postprocessor estimates a timestep that reduces the increment change in a material property below a given threshold.
  • NormalBoundaryDisplacementThis postprocessor computes the normal displacement on a given set of boundaries.
  • PolarMomentOfInertiaCompute the polar moment of inertia of a sideset w.r.t. a point and a direction
  • SidesetReactionComputes the integrated reaction force in a user-specified direction on a sideset from the surface traction
  • TorqueReactionTorqueReaction calculates the torque in 2D and 3Dabout a user-specified axis of rotation centeredat a user-specified origin.

Available Actions