- functionName of function to integrate
C++ Type:FunctionName
Controllable:No
Description:Name of function to integrate
FunctionElementIntegral
Integrates a function over elements
This post-processor computes the integral of some function over a domain :
Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- execute_onTIMESTEP_ENDThe list of flag(s) indicating when this object should be executed, the available options include FORWARD, ADJOINT, HOMOGENEOUS_FORWARD, ADJOINT_TIMESTEP_BEGIN, ADJOINT_TIMESTEP_END, NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, FINAL, CUSTOM, ALWAYS.
Default:TIMESTEP_END
C++ Type:ExecFlagEnum
Options:FORWARD, ADJOINT, HOMOGENEOUS_FORWARD, ADJOINT_TIMESTEP_BEGIN, ADJOINT_TIMESTEP_END, NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, FINAL, CUSTOM, ALWAYS, TRANSFER
Controllable:No
Description:The list of flag(s) indicating when this object should be executed, the available options include FORWARD, ADJOINT, HOMOGENEOUS_FORWARD, ADJOINT_TIMESTEP_BEGIN, ADJOINT_TIMESTEP_END, NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, FINAL, CUSTOM, ALWAYS.
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
Optional Parameters
- allow_duplicate_execution_on_initialFalseIn the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
Default:False
C++ Type:bool
Controllable:No
Description:In the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- execution_order_group0Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group.
Default:0
C++ Type:int
Controllable:No
Description:Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group.
- force_postauxFalseForces the UserObject to be executed in POSTAUX
Default:False
C++ Type:bool
Controllable:No
Description:Forces the UserObject to be executed in POSTAUX
- force_preauxFalseForces the UserObject to be executed in PREAUX
Default:False
C++ Type:bool
Controllable:No
Description:Forces the UserObject to be executed in PREAUX
- force_preicFalseForces the UserObject to be executed in PREIC during initial setup
Default:False
C++ Type:bool
Controllable:No
Description:Forces the UserObject to be executed in PREIC during initial setup
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- outputsVector of output names where you would like to restrict the output of variables(s) associated with this object
C++ Type:std::vector<OutputName>
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_element/ins/energy-conservation/q1q1.i)
- (test/tests/kernels/scalar_kernel_constraint/scalar_constraint_kernel.i)
- (test/tests/outputs/debug/show_execution_userobjects.i)
- (test/tests/kernels/scalar_kernel_constraint/diffusion_bipass_scalar.i)
- (test/tests/kernels/ad_scalar_kernel_constraint/scalar_constraint_together.i)
- (test/tests/kernels/ad_scalar_kernel_constraint/diffusion_override_scalar.i)
- (test/tests/outputs/debug/show_execution_ics.i)
- (test/tests/ics/postprocessor_interface/postprocessor_interface.i)
- (test/tests/ics/integral_preserving_function_ic/sinusoidal_z.i)
- (test/tests/postprocessors/function_element_integral/function_element_integral.i)
- (test/tests/kernels/scalar_kernel_constraint/diffusion_override_scalar.i)
- (test/tests/kernels/ad_scalar_kernel_constraint/scalar_constraint_kernel_RJ.i)
- (modules/navier_stokes/test/tests/finite_element/ins/energy-conservation/q2q1.i)
- (test/tests/kernels/ad_scalar_kernel_constraint/diffusion_bipass_scalar.i)
Child Objects
References
No citations exist within this document.(modules/navier_stokes/test/tests/finite_element/ins/energy-conservation/q1q1.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
nx = 10
ny = 10
dim = 2
[]
[subdomain]
type = SubdomainBoundingBoxGenerator
bottom_left = '0.5 0 0'
top_right = '1 1 0'
block_id = 1
input = gen
[]
[break_boundary]
input = subdomain
type = BreakBoundaryOnSubdomainGenerator
boundaries = 'bottom top'
[]
[sideset]
type = SideSetsBetweenSubdomainsGenerator
input = break_boundary
primary_block = '1'
paired_block = '0'
new_boundary = 'fluid_left'
[]
coord_type = RZ
[]
[Variables]
[T][]
[velocity]
family = LAGRANGE_VEC
block = 1
[]
[pressure]
block = 1
[]
[]
[ICs]
[vel]
type = VectorConstantIC
variable = velocity
x_value = 1e-15
y_value = 1e-15
block = 1
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = pressure
block = 1
[]
[pspg]
type = INSADMassPSPG
variable = pressure
block = 1
[]
[momentum_convection]
type = INSADMomentumAdvection
variable = velocity
block = 1
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = velocity
block = 1
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = pressure
integrate_p_by_parts = true
block = 1
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
block = 1
[]
[temperature_advection]
type = INSADEnergyAdvection
variable = T
block = 1
[]
[temperature_supg]
type = INSADEnergySUPG
variable = T
velocity = velocity
block = 1
[]
[temperature_conduction]
type = ADHeatConduction
variable = T
thermal_conductivity = 'k'
[]
[heat_source]
type = BodyForce
variable = T
block = 0
function = 'x + y'
[]
[]
[BCs]
[velocity_inlet]
type = VectorFunctionDirichletBC
variable = velocity
function_y = 1
boundary = 'bottom_to_1'
[]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'fluid_left right'
[]
[convective_heat_transfer]
type = ConvectiveHeatFluxBC
variable = T
T_infinity = 0
heat_transfer_coefficient = 1
boundary = 'right'
[]
[]
[Materials]
[constant]
type = ADGenericConstantMaterial
prop_names = 'cp rho k mu'
prop_values = '1 1 1 1'
[]
[ins]
type = INSADStabilized3Eqn
pressure = pressure
velocity = velocity
temperature = T
block = 1
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[Outputs]
csv = true
[]
[Postprocessors]
[convective_heat_transfer]
type = ConvectiveHeatTransferSideIntegral
T_solid = T
T_fluid = 0
htc = 1
boundary = 'right'
[]
[advection]
type = INSADElementIntegralEnergyAdvection
temperature = T
velocity = velocity
cp = cp
rho = rho
block = 1
[]
[source]
type = FunctionElementIntegral
function = 'x + y'
block = 0
[]
[energy_balance]
type = ParsedPostprocessor
function = 'convective_heat_transfer + advection - source'
pp_names = 'convective_heat_transfer advection source'
[]
[]
(test/tests/kernels/scalar_kernel_constraint/scalar_constraint_kernel.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[ffnk]
type = BodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
exodus = true
hide = lambda
[]
(test/tests/outputs/debug/show_execution_userobjects.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '1.5 2.4'
dy = '1.3 0.9'
ix = '3 2'
iy = '2 3'
subdomain_id = '0 1
1 0'
[]
[add_interface]
type = SideSetsBetweenSubdomainsGenerator
input = 'cmg'
primary_block = 0
paired_block = 1
new_boundary = 'interface'
[]
second_order = true
[]
[Functions]
[forcing_fnu]
type = ParsedFunction
value = -5.8*(x+y)+x*x*x-x+y*y*y-y
[]
[forcing_fnv]
type = ParsedFunction
value = -4
[]
[slnu]
type = ParsedGradFunction
value = x*x*x-x+y*y*y-y
grad_x = 3*x*x-1
grad_y = 3*y*y-1
[]
[slnv]
type = ParsedGradFunction
value = x*x+y*y
grad_x = 2*x
grad_y = 2*y
[]
# NeumannBC functions
[bc_fnut]
type = ParsedFunction
value = 3*y*y-1
[]
[bc_fnub]
type = ParsedFunction
value = -3*y*y+1
[]
[bc_fnul]
type = ParsedFunction
value = -3*x*x+1
[]
[bc_fnur]
type = ParsedFunction
value = 3*x*x-1
[]
[]
[Variables]
[u]
order = SECOND
family = HIERARCHIC
[]
[v]
order = SECOND
family = LAGRANGE
initial_condition = 1
[]
[]
[AuxVariables]
[v_elem]
order = CONSTANT
family = MONOMIAL
[]
[]
[Kernels]
active = 'diff1 diff2 test1 forceu forcev react'
[diff1]
type = Diffusion
variable = u
[]
[test1]
type = CoupledConvection
variable = u
velocity_vector = v
[]
[diff2]
type = Diffusion
variable = v
[]
[react]
type = Reaction
variable = u
[]
[forceu]
type = BodyForce
variable = u
function = forcing_fnu
[]
[forcev]
type = BodyForce
variable = v
function = forcing_fnv
[]
[]
[AuxKernels]
[set_v_elem]
type = FunctionAux
variable = v_elem
# selected not to be the solution for no particular reason
function = forcing_fnv
[]
[]
[BCs]
[bc_v]
type = FunctionDirichletBC
variable = v
function = slnv
boundary = 'left right top bottom'
[]
[bc_u_tb]
type = CoupledKernelGradBC
variable = u
var2 = v
vel = '0.1 0.1'
boundary = 'top bottom left right'
[]
[bc_ul]
type = FunctionNeumannBC
variable = u
function = bc_fnul
boundary = 'left'
[]
[bc_ur]
type = FunctionNeumannBC
variable = u
function = bc_fnur
boundary = 'right'
[]
[bc_ut]
type = FunctionNeumannBC
variable = u
function = bc_fnut
boundary = 'top'
[]
[bc_ub]
type = FunctionNeumannBC
variable = u
function = bc_fnub
boundary = 'bottom'
[]
[]
[Postprocessors]
# Global user objects
[dofs]
type = NumDOFs
[]
[h]
type = AverageElementSize
[]
# Elemental user objects
[L2u]
type = ElementL2Error
variable = u
function = slnu
# Testing an option
force_preic = true
[]
[L2v]
type = ElementL2Error
variable = v
function = slnv
# Testing an option
force_preaux = true
[]
[H1error]
type = ElementH1Error
variable = u
function = slnu
[]
[H1Semierror]
type = ElementH1SemiError
variable = u
function = slnu
[]
[L2v_elem]
type = ElementL2Error
variable = v_elem
function = slnv
[]
[f_integral]
type = FunctionElementIntegral
function = slnv
[]
[int_v]
type = ElementIntegralVariablePostprocessor
variable = v
block = 1
execute_on = 'TIMESTEP_END transfer'
[]
[int_v_elem]
type = ElementIntegralVariablePostprocessor
variable = v_elem
block = 1
execute_on = 'TIMESTEP_END transfer'
[]
# Side user objects
[integral_v]
type = SideIntegralVariablePostprocessor
variable = v
boundary = 0
[]
[]
[VectorPostprocessors]
# General UOs
[memory]
type = VectorMemoryUsage
[]
[line]
type = LineValueSampler
variable = v
num_points = 10
start_point = '0 0 0'
end_point = '0.5 0.5 0'
sort_by = 'x'
[]
# Nodal UOs
[nodal_sampler_y]
type = NodalValueSampler
variable = v
sort_by = 'y'
[]
[nodal_sampler_x]
type = NodalValueSampler
variable = v
sort_by = 'x'
[]
# Element UO
[elem_sample]
type = ElementValueSampler
variable = v_elem
sort_by = 'x'
[]
[]
[UserObjects]
# Nodal user objects
[find_node]
type = NearestNodeNumberUO
point = '0.5 0.5 0'
[]
# Side user objects
[side_int]
type = LayeredSideIntegral
variable = v
boundary = 0
direction = y
num_layers = 4
[]
[side_int_2]
type = NearestPointLayeredSideIntegral
variable = v
boundary = 0
direction = x
num_layers = 3
points = '1 1 0'
[]
# Interface user objects
[values]
type = InterfaceQpValueUserObject
var = v
boundary = interface
[]
inactive = 'prime_1 prime_2'
# Threaded general user objects
[prime_2]
type = PrimeProductUserObject
[]
[prime_1]
type = PrimeProductUserObject
[]
# Domain user objects
[domain_2]
type = InterfaceDomainUserObject
u = u
v = v
block = '0'
robin_boundaries = 'left'
interface_boundaries = 'interface'
interface_penalty = 1e-10
nl_abs_tol = 1e1
[]
[domain_1]
type = InterfaceDomainUserObject
u = u
v = v
block = '0 1'
robin_boundaries = 'left'
interface_boundaries = 'interface'
interface_penalty = 1e-10
nl_abs_tol = 1e1
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
nl_rel_tol = 1e-10
nl_abs_tol = 1e-10
l_tol = 1e-5
[]
[Problem]
kernel_coverage_check = false
[]
[MultiApps]
active = ''
[full_solve]
type = FullSolveMultiApp
execute_on = 'initial timestep_end final'
input_files = show_execution_userobjects.i
cli_args = 'Problem/solve=false'
[]
[]
[Transfers]
active = ''
[conservative]
type = MultiAppNearestNodeTransfer
from_multi_app = full_solve
source_variable = v
variable = v_elem
from_postprocessors_to_be_preserved = int_v
to_postprocessors_to_be_preserved = int_v_elem
[]
[]
[Debug]
show_execution_order = 'ALWAYS INITIAL NONLINEAR LINEAR TIMESTEP_BEGIN TIMESTEP_END FINAL'
[]
(test/tests/kernels/scalar_kernel_constraint/diffusion_bipass_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = DiffusionNoScalar
variable = u
[]
[ffnk]
type = BodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/kernels/ad_scalar_kernel_constraint/scalar_constraint_together.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[ffnk]
type = ADBodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = ADFunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = ADFunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = ADFunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = ADFunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Executioner]
type = Steady
residual_and_jacobian_together = true
nl_rel_tol = 1e-9
l_tol = 1.e-10
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
solve_type = NEWTON
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/kernels/ad_scalar_kernel_constraint/diffusion_override_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = ADDiffusionNoScalar
variable = u
scalar_variable = lambda
[]
[ffnk]
type = ADBodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem. ILU(0) seems to do OK in both serial and parallel in my testing,
# I have not seen any zero pivot issues.
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'bjacobi ilu'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/outputs/debug/show_execution_ics.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
[]
[Debug]
show_execution_order = ALWAYS
[]
[AuxVariables]
[a]
[]
[b]
[]
[]
[Variables]
[u]
[]
[v]
[]
[]
# From dependency ics test
[ICs]
[u_ic]
type = ConstantIC
variable = u
value = -1
[]
[v_ic]
type = MTICSum
variable = v
var1 = u
var2 = a
[]
[a_ic]
type = ConstantIC
variable = a
value = 10
[]
[b_ic]
type = MTICMult
variable = b
var1 = v
factor = 2
[]
[]
[AuxKernels]
[a_ak]
type = ConstantAux
variable = a
value = 256
[]
[b_ak]
type = ConstantAux
variable = b
value = 42
[]
[]
[Kernels]
[diff_u]
type = Diffusion
variable = u
[]
[diff_v]
type = Diffusion
variable = v
[]
[]
# From depend on uo test
[AuxVariables]
[ghost]
family = MONOMIAL
order = CONSTANT
[]
[]
[ICs]
[ghost_ic]
type = ElementUOIC
variable = ghost
element_user_object = ghost_uo
field_name = "ghosted"
field_type = long
[]
[]
[UserObjects]
[ghost_uo]
type = ElemSideNeighborLayersTester
execute_on = initial
element_side_neighbor_layers = 1
[]
[]
# From postprocessor interface ICs test
[Functions]
# The integral of this function is 2*3 + 3*6 + 5*2 = 34
[test_fn]
type = PiecewiseConstant
axis = x
x = '0 2 5'
y = '3 6 2'
[]
[]
[Postprocessors]
[integral_pp]
type = FunctionElementIntegral
function = test_fn
execute_on = 'INITIAL'
[]
[pp2]
type = FunctionValuePostprocessor
function = 6
execute_on = 'INITIAL'
[]
[]
[AuxVariables]
[test_var]
order = CONSTANT
family = MONOMIAL
[]
[]
[ICs]
[test_var_ic]
type = PostprocessorIC
variable = test_var
pp1 = integral_pp
[]
[]
# From integral preserving test
[AuxVariables]
[power]
family = MONOMIAL
order = CONSTANT
[]
[]
[ICs]
[power]
type = IntegralPreservingFunctionIC
variable = power
magnitude = 550.0
function = 'sin(pi * z / 1.9)'
integral = vol
[]
[]
[Postprocessors]
[vol]
type = FunctionElementIntegral
function = 'sin(pi * x / 1.9)'
execute_on = 'initial'
[]
[]
[BCs]
[left_u]
type = DirichletBC
variable = u
boundary = left
value = 0
[]
[right_u]
type = DirichletBC
variable = u
boundary = right
value = 1
[]
[left_v]
type = DirichletBC
variable = v
boundary = left
value = 2
[]
[right_v]
type = DirichletBC
variable = v
boundary = right
value = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-10
[]
(test/tests/ics/postprocessor_interface/postprocessor_interface.i)
[Mesh]
type = GeneratedMesh
dim = 1
xmin = 0
xmax = 10
nx = 10
[]
[Functions]
# The integral of this function is 2*3 + 3*6 + 5*2 = 34
[test_fn]
type = PiecewiseConstant
axis = x
x = '0 2 5'
y = '3 6 2'
[]
[]
[Postprocessors]
[integral_pp]
type = FunctionElementIntegral
function = test_fn
execute_on = 'INITIAL'
[]
[pp2]
type = FunctionValuePostprocessor
function = 6
execute_on = 'INITIAL'
[]
[]
[AuxVariables]
[test_var]
order = CONSTANT
family = MONOMIAL
[]
[]
[ICs]
[test_var_ic]
type = PostprocessorIC
variable = test_var
pp1 = integral_pp
[]
[]
[Problem]
solve = false
[]
[Executioner]
type = Steady
[]
[Postprocessors]
# This PP should have the sum of the other two PPs: 34 + 6 = 40
[test_var_pp]
type = ElementAverageValue
variable = test_var
execute_on = 'INITIAL'
[]
[]
[Outputs]
csv = true
[]
(test/tests/ics/integral_preserving_function_ic/sinusoidal_z.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 5
ny = 5
nz = 20
xmax = 1.5
ymax = 1.7
zmax = 1.9
xmin = 0.0
ymin = 0.0
zmin = 0.0
[]
[Problem]
type = FEProblem
solve = false
[]
[AuxVariables]
[power]
family = MONOMIAL
order = CONSTANT
[]
[]
[ICs]
[power]
type = IntegralPreservingFunctionIC
variable = power
magnitude = 550.0
function = 'sin(pi * z / 1.9)'
integral = vol
[]
[]
[Postprocessors]
[vol]
type = FunctionElementIntegral
function = 'sin(pi * z / 1.9)'
execute_on = 'initial'
[]
[integrated_power] # should equal 550
type = ElementIntegralVariablePostprocessor
variable = power
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]
(test/tests/postprocessors/function_element_integral/function_element_integral.i)
dx = 2
y1 = 3
y2 = 6
y3 = 8
integral = ${fparse dx * ((y1 + y2) * 0.5 + (y2 + y3) * 0.5)}
[Mesh]
type = GeneratedMesh
dim = 1
nx = 2
xmax = 4
[]
[Functions]
[./function]
type = PiecewiseLinear
axis = x
x = '0 2 4'
y = '${y1} ${y2} ${y3}'
[../]
[]
[Postprocessors]
[./integral_pp]
type = FunctionElementIntegral
function = function
execute_on = 'initial'
[../]
[./integral_rel_err]
type = RelativeDifferencePostprocessor
value1 = integral_pp
value2 = ${integral}
execute_on = 'initial'
[../]
[]
[Problem]
solve = false
[]
[Executioner]
type = Steady
[]
[Outputs]
csv = true
show = 'integral_rel_err'
[]
(test/tests/kernels/scalar_kernel_constraint/diffusion_override_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = DiffusionNoScalar
variable = u
scalar_variable = lambda
[]
[ffnk]
type = BodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem. ILU(0) seems to do OK in both serial and parallel in my testing,
# I have not seen any zero pivot issues.
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'bjacobi ilu'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/kernels/ad_scalar_kernel_constraint/scalar_constraint_kernel_RJ.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = ADDirichletBC
variable = u
boundary = 'bottom'
value = 0
[]
[right]
type = ADDirichletBC
variable = u
boundary = 'right'
value = 0
[]
[top]
type = ADDirichletBC
variable = u
boundary = 'top'
value = 0
[]
[left]
type = ADDirichletBC
variable = u
boundary = 'left'
value = 0
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[]
# Force LU decomposition, nonlinear iterations, to check Jacobian terms with single factorization
[Executioner]
type = Steady
residual_and_jacobian_together = true
nl_rel_tol = 1e-9
l_tol = 1.e-10
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
solve_type = NEWTON
[]
[Outputs]
exodus = true
hide = lambda
[]
(modules/navier_stokes/test/tests/finite_element/ins/energy-conservation/q2q1.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
nx = 10
ny = 10
dim = 2
[]
[subdomain]
type = SubdomainBoundingBoxGenerator
bottom_left = '0.5 0 0'
top_right = '1 1 0'
block_id = 1
input = gen
[]
[break_boundary]
input = subdomain
type = BreakBoundaryOnSubdomainGenerator
boundaries = 'bottom top'
[]
[sideset]
type = SideSetsBetweenSubdomainsGenerator
input = break_boundary
primary_block = '1'
paired_block = '0'
new_boundary = 'fluid_left'
[]
coord_type = RZ
second_order = true
[]
[Variables]
[T]
order = SECOND
[]
[velocity]
family = LAGRANGE_VEC
order = SECOND
block = 1
[]
[pressure]
block = 1
[]
[]
[ICs]
[vel]
type = VectorConstantIC
variable = velocity
x_value = 1e-15
y_value = 1e-15
block = 1
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = pressure
block = 1
[]
[momentum_convection]
type = INSADMomentumAdvection
variable = velocity
block = 1
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = velocity
block = 1
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = pressure
integrate_p_by_parts = true
block = 1
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
block = 1
[]
[temperature_advection]
type = INSADEnergyAdvection
variable = T
block = 1
[]
[temperature_supg]
type = INSADEnergySUPG
variable = T
velocity = velocity
block = 1
[]
[temperature_conduction]
type = ADHeatConduction
variable = T
thermal_conductivity = 'k'
[]
[heat_source]
type = BodyForce
variable = T
block = 0
function = 'x + y'
[]
[]
[BCs]
[velocity_inlet]
type = VectorFunctionDirichletBC
variable = velocity
function_y = 1
boundary = 'bottom_to_1'
[]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'fluid_left right'
[]
[convective_heat_transfer]
type = ConvectiveHeatFluxBC
variable = T
T_infinity = 0
heat_transfer_coefficient = 1
boundary = 'right'
[]
[]
[Materials]
[constant]
type = ADGenericConstantMaterial
prop_names = 'cp rho k mu'
prop_values = '1 1 1 1'
[]
[ins]
type = INSADStabilized3Eqn
pressure = pressure
velocity = velocity
temperature = T
block = 1
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[Outputs]
csv = true
[]
[Postprocessors]
[convective_heat_transfer]
type = ConvectiveHeatTransferSideIntegral
T_solid = T
T_fluid = 0
htc = 1
boundary = 'right'
[]
[advection]
type = INSADElementIntegralEnergyAdvection
temperature = T
velocity = velocity
cp = cp
rho = rho
block = 1
[]
[source]
type = FunctionElementIntegral
function = 'x + y'
block = 0
[]
[energy_balance]
type = ParsedPostprocessor
function = 'convective_heat_transfer + advection - source'
pp_names = 'convective_heat_transfer advection source'
[]
[]
(test/tests/kernels/ad_scalar_kernel_constraint/diffusion_bipass_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = ADDiffusionNoScalar
variable = u
[]
[ffnk]
type = ADBodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = ADFunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = ADFunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = ADFunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = ADFunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
residual_and_jacobian_together = true
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(modules/thermal_hydraulics/include/postprocessors/FunctionElementIntegralRZ.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FunctionElementIntegral.h"
#include "RZSymmetry.h"
/**
* Integrates a function over elements for RZ geometry modeled by XY domain
*/
class FunctionElementIntegralRZ : public FunctionElementIntegral, public RZSymmetry
{
public:
FunctionElementIntegralRZ(const InputParameters & parameters);
protected:
virtual Real computeQpIntegral() override;
public:
static InputParameters validParams();
};
(framework/include/postprocessors/FunctionElementAverage.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "FunctionElementIntegral.h"
/**
* Computes the average of a function over a volume.
*/
class FunctionElementAverage : public FunctionElementIntegral
{
public:
static InputParameters validParams();
FunctionElementAverage(const InputParameters & parameters);
virtual void initialize() override;
virtual void execute() override;
virtual Real getValue() override;
virtual void finalize() override;
virtual void threadJoin(const UserObject & y) override;
protected:
Real _volume;
};