- vThe coupled variable whose components are coupled to AuxVariable
C++ Type:std::vector<VariableName>
Controllable:No
Description:The coupled variable whose components are coupled to AuxVariable
- variableThe name of the variable that this object applies to
C++ Type:AuxVariableName
Controllable:No
Description:The name of the variable that this object applies to
- vector_tagTag Name this Aux works on
C++ Type:TagName
Controllable:No
Description:Tag Name this Aux works on
TagVectorAux
The value of a tagged vector for a given node and a given variable is coupled to the current AuxVariable. TagVectorAux return the coupled nodal value. AuxVariable then is written out in an exodus file.
Couple a tag vector, and return its nodal value
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
- boundaryThe list of boundaries (ids or names) from the mesh where this boundary condition applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this boundary condition applies
- check_boundary_restrictedTrueWhether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
Default:True
C++ Type:bool
Controllable:No
Description:Whether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
- execute_onTIMESTEP_ENDThe list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM, PRE_DISPLACE, ALWAYS.
Default:TIMESTEP_END
C++ Type:ExecFlagEnum
Options:NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM, PRE_DISPLACE, ALWAYS
Controllable:No
Description:The list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM, PRE_DISPLACE, 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.
- scaledTrueReturn value depending on the variable scaling/autoscaling. Set this to false to obtain unscaled physical reaction forces.
Default:True
C++ Type:bool
Controllable:No
Description:Return value depending on the variable scaling/autoscaling. Set this to false to obtain unscaled physical reaction forces.
Optional Parameters
- 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.
- 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/tensor_mechanics/test/tests/plane_stress/3D_finite_tension_pull.i)
- (modules/tensor_mechanics/test/tests/generalized_plane_strain/generalized_plane_strain_ref_resid.i)
- (test/tests/tag/2d_diffusion_tag_vector.i)
- (test/tests/tag/tag_interface_kernels.i)
- (test/tests/tag/2d_diffusion_dg_tag.i)
- (test/tests/tag/tag_nodal_kernels.i)
- (modules/heat_conduction/test/tests/convective_flux_function/convective_flux_function.i)
- (test/tests/tag/eigen_tag.i)
- (test/tests/tag/2d_diffusion_vector_tag_test.i)
- (modules/tensor_mechanics/test/tests/torque_reaction/torque_reaction.i)
- (test/tests/tag/tag_neumann.i)
- (modules/tensor_mechanics/test/tests/plane_stress/weak_plane_stress_finite_tension_pull.i)
- (modules/tensor_mechanics/test/tests/torque_reaction/torque_reaction_3D.i)
- (modules/tensor_mechanics/test/tests/torque_reaction/torque_reaction_cylinder.i)
- (test/tests/tag/tag_dirac_kernels.i)
References
No citations exist within this document.(modules/tensor_mechanics/test/tests/plane_stress/3D_finite_tension_pull.i)
[GlobalParams]
order = FIRST
family = LAGRANGE
displacements = 'disp_x disp_y disp_z'
[]
[Problem]
extra_tag_vectors = 'ref'
[]
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 3
nx = 1
ny = 1
nz = 1
[]
[]
[AuxVariables]
[react_x]
[]
[]
[Postprocessors]
[react_x]
type = NodalSum
variable = 'react_x'
boundary = 'right'
[]
[stress_xx]
type = ElementalVariableValue
variable = 'stress_xx'
elementid = 0
[]
[strain_zz]
type = ElementalVariableValue
variable = 'strain_zz'
elementid = 0
[]
[]
[Modules/TensorMechanics/Master]
[plane_stress]
strain = FINITE
extra_vector_tags = 'ref'
generate_output = 'stress_xx stress_xy stress_yy stress_zz strain_xx strain_xy strain_yy strain_zz'
add_variables = true
[]
[]
[AuxKernels]
[react_x]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_x'
variable = 'react_x'
[]
[]
[BCs]
[leftx]
type = DirichletBC
boundary = left
variable = disp_x
value = 0.0
[]
[bottomy]
type = DirichletBC
boundary = bottom
variable = disp_y
value = 0.0
[]
[backz]
type = DirichletBC
boundary = back
variable = disp_z
value = 0.0
[]
[rightx]
type = FunctionDirichletBC
boundary = right
variable = disp_x
function = 't'
[]
[]
[Materials]
[elasticity_tensor]
type = ComputeIsotropicElasticityTensor
poissons_ratio = 0.3
youngs_modulus = 1e6
[]
[stress]
type = ComputeFiniteStrainElasticStress
[]
[]
[Executioner]
type = Transient
solve_type = PJFNK
line_search = none
nl_rel_tol = 1e-12
nl_abs_tol = 1e-10
# time control
start_time = 0.0
dt = 0.01
dtmin = 0.01
end_time = 0.2
[]
[Outputs]
csv = true
[]
(modules/tensor_mechanics/test/tests/generalized_plane_strain/generalized_plane_strain_ref_resid.i)
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[./square]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 2
[../]
[]
[Problem]
type = ReferenceResidualProblem
reference_vector = 'ref'
extra_tag_vectors = 'ref'
group_variables = 'disp_x disp_y'
[]
[Variables]
[./scalar_strain_zz]
order = FIRST
family = SCALAR
[../]
[]
[AuxVariables]
[./temp]
order = FIRST
family = LAGRANGE
[../]
[./saved_x]
order = FIRST
family = LAGRANGE
[../]
[./saved_y]
order = FIRST
family = LAGRANGE
[../]
[]
[Postprocessors]
[./react_z]
type = MaterialTensorIntegral
rank_two_tensor = stress
index_i = 2
index_j = 2
[../]
[]
[Modules/TensorMechanics/Master]
[./all]
strain = SMALL
add_variables = true
displacements = 'disp_x disp_y'
generate_output = 'stress_xx stress_xy stress_yy stress_zz strain_xx strain_xy strain_yy strain_zz'
planar_formulation = GENERALIZED_PLANE_STRAIN
eigenstrain_names = eigenstrain
scalar_out_of_plane_strain = scalar_strain_zz
temperature = temp
extra_vector_tags = 'ref'
[../]
[]
[AuxKernels]
[./tempfuncaux]
type = FunctionAux
variable = temp
function = tempfunc
use_displaced_mesh = false
[../]
[./saved_x]
type = TagVectorAux
variable = 'saved_x'
vector_tag = 'ref'
v = 'disp_x'
execute_on = timestep_end
[../]
[./saved_y]
type = TagVectorAux
variable = 'saved_y'
vector_tag = 'ref'
execute_on = timestep_end
v = 'disp_y'
[../]
[]
[Functions]
[./tempfunc]
type = ParsedFunction
value = '(1-x)*t'
[../]
[]
[BCs]
[./bottomx]
type = DirichletBC
boundary = 0
variable = disp_x
value = 0.0
[../]
[./bottomy]
type = DirichletBC
boundary = 0
variable = disp_y
value = 0.0
[../]
[]
[Materials]
[./elastic_tensor]
type = ComputeIsotropicElasticityTensor
poissons_ratio = 0.3
youngs_modulus = 1e6
[../]
[./thermal_strain]
type = ComputeThermalExpansionEigenstrain
temperature = temp
thermal_expansion_coeff = 0.02
stress_free_temperature = 0.5
eigenstrain_name = eigenstrain
[../]
[./stress]
type = ComputeLinearElasticStress
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = PJFNK
line_search = none
# controls for linear iterations
l_max_its = 100
l_tol = 1e-4
# controls for nonlinear iterations
nl_max_its = 15
nl_rel_tol = 1e-10
nl_abs_tol = 1e-5
# time control
start_time = 0.0
dt = 1.0
dtmin = 1.0
end_time = 2.0
num_steps = 5000
[]
[Outputs]
exodus = true
[]
(test/tests/tag/2d_diffusion_tag_vector.i)
[Mesh]
[./square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[../]
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = LAGRANGE
[../]
[./tag_variable2]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag1
[../]
[./TagVectorAux2]
type = TagVectorAux
variable = tag_variable2
v = u
vector_tag = vec_tag2
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = 3
value = 0
preset = false
extra_vector_tags = vec_tag1
[../]
[./right]
type = DirichletBC
variable = u
boundary = 1
value = 1
preset = false
extra_vector_tags = vec_tag2
[../]
[]
[Problem]
type = FEProblem
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
file_base = tag_vector_out
exodus = true
[]
(test/tests/tag/tag_interface_kernels.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 3
nx = 2
xmax = 2
ny = 2
ymax = 2
nz = 2
zmax = 2
[]
[./subdomain1]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0 0'
top_right = '1 1 1'
block_id = 1
[../]
[./break_boundary]
input = subdomain1
type = BreakBoundaryOnSubdomainGenerator
[../]
[./interface]
type = SideSetsBetweenSubdomainsGenerator
input = break_boundary
primary_block = '0'
paired_block = '1'
new_boundary = 'primary0_interface'
[../]
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
block = 0
[../]
[./v]
order = FIRST
family = LAGRANGE
block = 1
[../]
[]
[Kernels]
[./diff_u]
type = CoeffParamDiffusion
variable = u
D = 4
block = 0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./diff_v]
type = CoeffParamDiffusion
variable = v
D = 2
block = 1
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./source_u]
type = BodyForce
variable = u
value = 1
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[InterfaceKernels]
[./interface]
type = PenaltyInterfaceDiffusion
variable = u
neighbor_var = v
boundary = primary0_interface
penalty = 1e6
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[BCs]
[./u]
type = VacuumBC
variable = u
boundary = 'left_to_0 bottom_to_0 back_to_0 right top front'
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./v]
type = VacuumBC
variable = v
boundary = 'left_to_1 bottom_to_1 back_to_1'
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = LAGRANGE
block = 0
[../]
[./tag_variable2]
order = FIRST
family = LAGRANGE
block = 1
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = u
block = 0
vector_tag = vec_tag2
execute_on = timestep_end
[../]
[./TagVectorAux2]
type = TagMatrixAux
variable = tag_variable2
v = v
block = 1
matrix_tag = mat_tag2
execute_on = timestep_end
[../]
[]
[Postprocessors]
[./u_int]
type = ElementIntegralVariablePostprocessor
variable = u
block = 0
[../]
[./v_int]
type = ElementIntegralVariablePostprocessor
variable = v
block = 1
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'nontime residual vec_tag1 vec_tag2'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
exodus = true
[]
(test/tests/tag/2d_diffusion_dg_tag.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[./InitialCondition]
type = ConstantIC
value = 1
[../]
[../]
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = MONOMIAL
[../]
[./tag_variable2]
order = FIRST
family = MONOMIAL
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag2
execute_on = timestep_end
[../]
[./TagVectorAux2]
type = TagMatrixAux
variable = tag_variable2
v = u
matrix_tag = mat_tag2
execute_on = timestep_end
[../]
[]
[Functions]
[./forcing_fn]
type = ParsedFunction
value = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
value = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[./abs]
type = Reaction
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[./forcing]
type = BodyForce
variable = u
function = forcing_fn
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[DGKernels]
[./dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'nontime residual vec_tag1 vec_tag2'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-10
[]
[Postprocessors]
[./h]
type = AverageElementSize
[../]
[./dofs]
type = NumDOFs
[../]
[./l2_err]
type = ElementL2Error
variable = u
function = exact_fn
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/tag/tag_nodal_kernels.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[./u]
[../]
[./nodal_ode]
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./time]
type = TimeDerivative
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[NodalKernels]
[./td]
type = TimeDerivativeNodalKernel
variable = nodal_ode
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./constant_rate]
type = ConstantRate
variable = nodal_ode
rate = 1.0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 10
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'time nontime residual vec_tag1 vec_tag2'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = LAGRANGE
[../]
[./tag_variable2]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = nodal_ode
vector_tag = vec_tag2
execute_on = timestep_end
[../]
[./TagVectorAux2]
type = TagMatrixAux
variable = tag_variable2
v = u
matrix_tag = mat_tag2
execute_on = timestep_end
[../]
[]
[Executioner]
type = Transient
num_steps = 10
nl_rel_tol = 1e-08
dt = 0.01
[]
[Outputs]
exodus = true
[]
(modules/heat_conduction/test/tests/convective_flux_function/convective_flux_function.i)
# This is a test of the ConvectiveFluxFunction BC.
# There is a single 1x1 element with a prescribed temperature
# on the left side and a convective flux BC on the right side.
# The temperature on the left is 100, and the far-field temp is 200.
# The conductance of the body (conductivity * length) is 10
#
# If the conductance in the BC is also 10, the temperature on the
# right side of the solid element should be 150 because half of the
# temperature drop should occur over the body and half in the BC.
#
# The integrated flux is deltaT * conductance, or -50 * 10 = -500.
# The negative sign indicates that heat is going into the body.
#
# The conductance is defined multiple ways using this input, and
# as long as it evaluates to 10, the result described above will
# be obtained.
[Mesh]
type = GeneratedMesh
dim = 2
nx = 1
ny = 1
[]
[Problem]
extra_tag_vectors = 'bcs'
[]
[Variables]
[temp]
initial_condition = 100.0
[]
[]
[AuxVariables]
[flux]
[]
[]
[AuxKernels]
[flux]
type = TagVectorAux
variable = flux
v = temp
vector_tag = 'bcs'
execute_on = timestep_end
[]
[]
[Kernels]
[heat_conduction]
type = HeatConduction
variable = temp
[]
[]
[Materials]
[thermal]
type = HeatConductionMaterial
thermal_conductivity = 10.0
[]
[]
[BCs]
[left]
type = DirichletBC
variable = temp
boundary = left
value = 100.0
[]
[right]
type = ConvectiveFluxFunction
variable = temp
boundary = right
T_infinity = 200.0
coefficient = 10.0 #This will behave as described in the header of this file if this evaluates to 10
extra_vector_tags = 'bcs'
[]
[]
[Postprocessors]
[integrated_flux]
type = NodalSum
variable = flux
boundary = right
[]
[]
[Executioner]
type = Transient
start_time = 0.0
end_time = 1.0
dt = 1.0
nl_rel_tol=1e-12
[]
[Outputs]
csv = true
[]
(test/tests/tag/eigen_tag.i)
[Mesh/gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 5
ny = 5
[]
[Variables/u]
[]
[AuxVariables]
[vec_tag_diff]
order = FIRST
family = LAGRANGE
[]
[vec_tag_rhs]
order = FIRST
family = LAGRANGE
[]
[mat_tag_diff]
order = FIRST
family = LAGRANGE
[]
[mat_tag_rhs]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
extra_vector_tags = 'tag_diff'
extra_matrix_tags = 'tag_diff'
[]
[rhs]
type = CoefReaction
variable = u
extra_vector_tags = 'eigen tag_rhs'
extra_matrix_tags = 'tag_rhs'
[]
[]
[AuxKernels]
[vec_tag_diff]
type = TagVectorAux
variable = vec_tag_diff
v = u
vector_tag = tag_diff
[]
[vec_tag_rhs]
type = TagVectorAux
variable = vec_tag_rhs
v = u
vector_tag = tag_rhs
[]
[mat_tag_diff]
type = TagVectorAux
variable = mat_tag_diff
v = u
vector_tag = tag_diff
[]
[mat_tag_rhs]
type = TagVectorAux
variable = mat_tag_diff
v = u
vector_tag = tag_rhs
[]
[]
[BCs/homogeneous]
type = DirichletBC
boundary = 'top right bottom left'
variable = u
value = 0
[]
[Problem]
extra_tag_vectors = 'tag_diff tag_rhs'
extra_tag_matrices = 'tag_diff tag_rhs'
[]
[Executioner]
type = Eigenvalue
solve_type = NEWTON
eigen_problem_type = GEN_NON_HERMITIAN
[]
[Outputs]
exodus = true
[]
(test/tests/tag/2d_diffusion_vector_tag_test.i)
[Mesh]
[./square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[../]
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
initial_condition = 1.0
[../]
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = LAGRANGE
[../]
[./tag_variable2]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./reaction1]
type = Reaction
variable = u
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[./reaction2]
type = Reaction
variable = u
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[./reaction3]
type = Reaction
variable = u
[../]
[./reaction4]
type = Reaction
variable = u
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag1
execute_on = timestep_end
[../]
[./TagVectorAux2]
type = TagVectorAux
variable = tag_variable2
v = u
vector_tag = vec_tag2
execute_on = timestep_end
[../]
[]
[BCs]
active = 'left right'
[./left]
type = DirichletBC
variable = u
preset = false
boundary = 3
value = 10
extra_vector_tags = vec_tag1
[../]
[./right]
type = DirichletBC
variable = u
preset = false
boundary = 1
value = 100
extra_vector_tags = vec_tag2
[../]
[./right1]
type = DirichletBC
variable = u
preset = false
boundary = 1
value = 100
[../]
[./right2]
type = DirichletBC
variable = u
preset = false
boundary = 1
value = 100
[../]
[]
[Problem]
type = TagTestProblem
extra_tag_vectors = 'vec_tag1 vec_tag2'
test_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
file_base = vector_tag_test_out
exodus = true
[]
(modules/tensor_mechanics/test/tests/torque_reaction/torque_reaction.i)
# Scalar torque reaction
# This test computes the sum of the torques acting on a ten element 2D bar mesh
# and is intended to replicate the classical wrench problem from statics.
# A displacement in the y along the right face is applied to the bar end to create
# a shear force along the bar end. The rotation origin default (the global origin)
# and the axis of rotation direction vector used to compute the torque reaction
# is set to (0, 0, 1) out of the plane.
# Torque is calculated for the two nodes on the left of the bar. For the bottom
# node on the right, the torque/ moment lever is the x coordinate value, and for
# the top node on the right the torque lever is the hypotenuse of the x and y
# coordinates. The expected sum of the torque reaction is just over 37.
[GlobalParams]
order = FIRST
family = LAGRANGE
displacements = 'disp_x disp_y'
[]
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 1
[]
[Problem]
extra_tag_vectors = 'ref'
[]
[AuxVariables]
[./saved_x]
[../]
[./saved_y]
[../]
[]
[AuxKernels]
[saved_x]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_x'
variable = 'saved_x'
[]
[saved_y]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_y'
variable = 'saved_y'
[]
[]
[Modules/TensorMechanics/Master]
[master]
strain = SMALL
generate_output = 'stress_xx stress_yy'
add_variables = true
extra_vector_tags = 'ref'
[]
[]
[BCs]
[./left_x]
type = DirichletBC
variable = disp_x
boundary = left
value = 0.0
[../]
[./left_y]
type = DirichletBC
variable = disp_y
boundary = left
value = 0.0
[../]
[./right_shear_y]
type = FunctionDirichletBC
variable = disp_y
boundary = right
function = '0.001*t'
[../]
[]
[Materials]
[./elasticity_tensor]
type = ComputeIsotropicElasticityTensor
youngs_modulus = 207000
poissons_ratio = 0.3
[../]
[./elastic_stress]
type = ComputeLinearElasticStress
[../]
[]
[Executioner]
type = Transient
line_search = 'none'
l_max_its = 30
nl_max_its = 20
nl_abs_tol = 1e-12
nl_rel_tol = 1e-10
l_tol = 1e-8
start_time = 0.0
dt = 0.5
end_time = 1
num_steps = 2
[]
[Postprocessors]
[./_dt]
type = TimestepSize
[../]
[./torque]
type = TorqueReaction
boundary = right
reaction_force_variables = 'saved_x saved_y'
direction_vector = '0. 0. 1.'
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/tag/tag_neumann.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
construct_side_list_from_node_list = true
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = LAGRANGE
[../]
[./tag_variable2]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag2
execute_on = timestep_end
[../]
[./TagVectorAux2]
type = TagMatrixAux
variable = tag_variable2
v = u
matrix_tag = mat_tag2
execute_on = timestep_end
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./right]
type = NeumannBC
variable = u
boundary = right
value = 2
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'nontime residual vec_tag1 vec_tag2'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/tensor_mechanics/test/tests/plane_stress/weak_plane_stress_finite_tension_pull.i)
[GlobalParams]
order = FIRST
family = LAGRANGE
displacements = 'disp_x disp_y'
out_of_plane_strain = strain_zz
[]
[Problem]
extra_tag_vectors = 'ref'
[]
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 1
ny = 1
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[strain_zz]
[]
[]
[AuxVariables]
[react_x]
[]
[]
[Postprocessors]
[react_x]
type = NodalSum
variable = 'react_x'
boundary = 'right'
[]
[stress_xx]
type = ElementalVariableValue
variable = 'stress_xx'
elementid = 0
[]
[strain_zz]
type = ElementalVariableValue
variable = 'strain_zz'
elementid = 0
[]
[]
[Modules/TensorMechanics/Master]
[plane_stress]
strain = FINITE
planar_formulation = WEAK_PLANE_STRESS
extra_vector_tags = 'ref'
generate_output = 'stress_xx stress_xy stress_yy stress_zz strain_xx strain_xy strain_yy'
[]
[]
[AuxKernels]
[react_x]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_x'
variable = 'react_x'
[]
[]
[BCs]
[leftx]
type = DirichletBC
boundary = left
variable = disp_x
value = 0.0
[]
[bottomy]
type = DirichletBC
boundary = bottom
variable = disp_y
value = 0.0
[]
[rightx]
type = FunctionDirichletBC
boundary = right
variable = disp_x
function = 't'
[]
[]
[Materials]
[elasticity_tensor]
type = ComputeIsotropicElasticityTensor
poissons_ratio = 0.3
youngs_modulus = 1e6
[]
[stress]
type = ComputeFiniteStrainElasticStress
[]
[]
[Executioner]
type = Transient
solve_type = PJFNK
line_search = none
nl_rel_tol = 1e-12
nl_abs_tol = 1e-10
# time control
start_time = 0.0
dt = 0.01
dtmin = 0.01
end_time = 0.2
[]
[Outputs]
csv = true
[]
(modules/tensor_mechanics/test/tests/torque_reaction/torque_reaction_3D.i)
# Scalar torque reaction
# This test computes the sum of the torques acting on a single element cube mesh.
# Equal displacements in the x and the z are applied along the cube top to
# create a shear force along the (1, 0, 1) direction. The rotation origin is
# set to the middle of the bottom face of the cube (0.5, 0, 0.5), and the axis of
# rotation direction vector used to compute the torque reaction is set to (-1, 0, 1).
# Torque is calculated for the four nodes on the top of the cube. The projection
# of the node coordinates is zero for nodes 3 and 6, +1 for node 7, and -1 for
# node 2 from the selection of the direction vector and the rotation axis origin.
[GlobalParams]
order = FIRST
family = LAGRANGE
displacements = 'disp_x disp_y disp_z'
[]
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
[]
[Problem]
extra_tag_vectors = 'ref'
[]
[GlobalParams]
volumetric_locking_correction=true
[]
[AuxVariables]
[./saved_x]
[../]
[./saved_y]
[../]
[./saved_z]
[../]
[]
[AuxKernels]
[saved_x]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_x'
variable = 'saved_x'
[]
[saved_y]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_y'
variable = 'saved_y'
[]
[saved_z]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_z'
variable = 'saved_z'
[]
[]
[Modules/TensorMechanics/Master]
[master]
strain = SMALL
generate_output = 'stress_xx stress_yy stress_zz'
add_variables = true
extra_vector_tags = 'ref'
[]
[]
[BCs]
[./bottom_x]
type = DirichletBC
variable = disp_x
boundary = bottom
value = 0.0
[../]
[./bottom_y]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0.0
[../]
[./bottom_z]
type = DirichletBC
variable = disp_z
boundary = bottom
value = 0.0
[../]
[./top_shear_z]
type = FunctionDirichletBC
variable = disp_z
boundary = top
function = '0.01*t'
[../]
[./top_shear_x]
type = FunctionDirichletBC
variable = disp_x
boundary = top
function = '0.01*t'
[../]
[]
[Materials]
[./elasticity_tensor]
type = ComputeIsotropicElasticityTensor
youngs_modulus = 207000
poissons_ratio = 0.3
[../]
[./elastic_stress]
type = ComputeLinearElasticStress
[../]
[]
[Executioner]
type = Transient
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-ksp_gmres_restart'
petsc_options_value = '101'
line_search = 'none'
l_max_its = 30
nl_max_its = 20
nl_abs_tol = 1e-14
nl_rel_tol = 1e-12
l_tol = 1e-8
start_time = 0.0
dt = 0.5
end_time = 1
num_steps = 2
[]
[Postprocessors]
[./_dt]
type = TimestepSize
[../]
[./torque]
type = TorqueReaction
boundary = top
reaction_force_variables = 'saved_x saved_y saved_z'
axis_origin = '0.5 0. 0.5'
direction_vector = '-1. 0. 1.'
[../]
[]
[Outputs]
exodus = true
[]
(modules/tensor_mechanics/test/tests/torque_reaction/torque_reaction_cylinder.i)
# This test uses the DisplacementAboutAxis boundary condition to twist the top
# of a cylinder while the bottom face of the cylinder remains fixed. The
# TorqueReaction postprocessor is used to calculate the applied torque acting
# on the cylinder at the top face. This test can be extended, with a new mesh,
# to model a crack in the center of the cylinder face under type III loading.
[GlobalParams]
order = FIRST
family = LAGRANGE
displacements = 'disp_x disp_y disp_z'
[]
[Mesh]
file = cylinder.e
[]
[Problem]
extra_tag_vectors = 'ref'
[]
[GlobalParams]
volumetric_locking_correction=true
[]
[AuxVariables]
[./saved_x]
[../]
[./saved_y]
[../]
[./saved_z]
[../]
[]
[AuxKernels]
[saved_x]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_x'
variable = 'saved_x'
[]
[saved_y]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_y'
variable = 'saved_y'
[]
[saved_z]
type = TagVectorAux
vector_tag = 'ref'
v = 'disp_z'
variable = 'saved_z'
[]
[]
[Modules/TensorMechanics/Master]
[master]
strain = FINITE
generate_output = 'stress_xx'
add_variables = true
extra_vector_tags = 'ref'
[]
[]
[BCs]
[./bottom_x]
type = DirichletBC
variable = disp_x
boundary = 1
value = 0.0
[../]
[./bottom_y]
type = DirichletBC
variable = disp_y
boundary = 1
value = 0.0
[../]
[./bottom_z]
type = DirichletBC
variable = disp_z
boundary = 1
value = 0.0
[../]
[./top_x]
type = DisplacementAboutAxis
boundary = 2
function = '0.1*t'
angle_units = degrees
axis_origin = '10. 10. 10.'
axis_direction = '0 -1.0 1.0'
component = 0
variable = disp_x
[../]
[./top_y]
type = DisplacementAboutAxis
boundary = 2
function = '0.1*t'
angle_units = degrees
axis_origin = '10. 10. 10.'
axis_direction = '0 -1.0 1.0'
component = 1
variable = disp_y
[../]
[./top_z]
type = DisplacementAboutAxis
boundary = 2
function = '0.1*t'
angle_units = degrees
axis_origin = '10. 10. 10.'
axis_direction = '0 -1.0 1.0'
component = 2
variable = disp_z
[../]
[]
[Materials]
[./elasticity_tensor]
type = ComputeIsotropicElasticityTensor
youngs_modulus = 207000
poissons_ratio = 0.3
[../]
[./elastic_stress]
type = ComputeFiniteStrainElasticStress
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-ksp_gmres_restart'
petsc_options_value = '101'
line_search = 'none'
l_max_its = 50
nl_max_its = 20
nl_abs_tol = 1e-12
nl_rel_tol = 1e-11
l_tol = 1e-10
start_time = 0.0
dt = 0.25
end_time = 0.5
[]
[Postprocessors]
[./torque]
type = TorqueReaction
boundary = 2
reaction_force_variables = 'saved_x saved_y saved_z'
axis_origin = '10. 10. 10.'
direction_vector = '0 -1.0 1.0'
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/tag/tag_dirac_kernels.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = QUAD4
uniform_refine = 4
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./ddt_u]
type = TimeDerivative
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./diff_u]
type = Diffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./ddt_v]
type = TimeDerivative
variable = v
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./diff_v]
type = Diffusion
variable = v
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[DiracKernels]
[./nonlinear_source]
type = NonlinearSource
variable = u
coupled_var = v
scale_factor = 1000
point = '0.2 0.3 0'
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[BCs]
[./left_u]
type = DirichletBC
variable = u
boundary = 3
value = 0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./right_u]
type = DirichletBC
variable = u
boundary = 1
value = 1
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./left_v]
type = DirichletBC
variable = v
boundary = 3
value = 1
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./right_v]
type = DirichletBC
variable = v
boundary = 1
value = 0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[Preconditioning]
[./precond]
type = SMP
full = true
[../]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'time nontime residual vec_tag1 vec_tag2'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = LAGRANGE
[../]
[./tag_variable2]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag2
execute_on = timestep_end
[../]
[./TagVectorAux2]
type = TagMatrixAux
variable = tag_variable2
v = u
matrix_tag = mat_tag2
execute_on = timestep_end
[../]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON' # NEWTON provides a more stringent test of off-diagonal Jacobians
num_steps = 5
dt = 1
dtmin = 1
l_max_its = 100
nl_max_its = 6
nl_abs_tol = 1.e-13
[]
[Postprocessors]
[./point_value]
type = PointValue
variable = u
point = '0.2 0.3 0'
[../]
[]
[Outputs]
exodus = true
[]