- pointThe x,y,z coordinates of the point
C++ Type:std::vector<double>
Controllable:No
Description:The x,y,z coordinates of the point
- valueThe value of the point source
C++ Type:double
Controllable:Yes
Description:The value of the point source
- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Controllable:No
Description:The name of the variable that this residual object operates on
ConstantPointSource
Residual contribution of a constant point source term.
This applies a load in a single location in the mesh. The value field is controllable, so the Control system
may be leveraged to control the load during the simulation. Alternatively, a FunctionDiracSource may be used for spatial and temporal variations of the point source using a Function
.
Example input syntax
In this example, three ConstantPointSource
are being applied to variable u
with values 0.1 / -0.1 and -1 at position (0.2 0.3 0.0) / (0.2 0.8 0) and (0.8 0.5 0.8) respectively. u
is solution to a diffusion equation with those three sources.
[DiracKernels]
[./point_source1]
type = ConstantPointSource
variable = u
value = 0.1
point = '0.2 0.3 0.0'
[../]
[./point_source2]
type = ConstantPointSource
variable = u
value = -0.1
point = '0.2 0.8 0.0'
[../]
[./point_source3]
type = ConstantPointSource
variable = u
value = -1.0
point = '0.8 0.5 0.8'
[../]
[]
(test/tests/dirackernels/constant_point_source/3d_point_source.i)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
- point_not_found_behaviorIGNOREBy default (IGNORE), it is ignored if an added point cannot be located in the specified subdomains. If this option is set to ERROR, this situation will result in an error. If this option is set to WARNING, then a warning will be issued.
Default:IGNORE
C++ Type:MooseEnum
Options:ERROR, WARNING, IGNORE
Controllable:No
Description:By default (IGNORE), it is ignored if an added point cannot be located in the specified subdomains. If this option is set to ERROR, this situation will result in an error. If this option is set to WARNING, then a warning will be issued.
- 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
- 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.
- drop_duplicate_pointsTrueBy default points added to a DiracKernel are dropped if a point at the same locationhas been added before. If this option is set to false duplicate points are retainedand contribute to residual and Jacobian.
Default:True
C++ Type:bool
Controllable:No
Description:By default points added to a DiracKernel are dropped if a point at the same locationhas been added before. If this option is set to false duplicate points are retainedand contribute to residual and Jacobian.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- 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
- 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
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
Input Files
- (test/tests/dirackernels/constant_point_source/1d_point_source.i)
- (modules/porous_flow/test/tests/fluidstate/water_vapor_phasechange.i)
- (modules/geochemistry/test/tests/kernels/time_deriv_2.i)
- (test/tests/outputs/oversample/ex02_adapt.i)
- (examples/ex02_kernel/ex02_oversample.i)
- (modules/tensor_mechanics/test/tests/shell/static/pinched_cylinder_symm.i)
- (modules/tensor_mechanics/test/tests/truss/truss_hex.i)
- (test/tests/dirackernels/constant_point_source/2d_point_source.i)
- (modules/tensor_mechanics/test/tests/truss/truss_2d_action.i)
- (test/tests/dirackernels/constant_point_source/1d_point_source_fv.i)
- (modules/tensor_mechanics/test/tests/truss/truss_2d.i)
- (test/tests/outputs/oversample/ex02_oversample.i)
- (test/tests/dirackernels/constant_point_source/3d_point_source.i)
- (modules/tensor_mechanics/test/tests/truss/truss_hex_action.i)
- (test/tests/controls/time_periods/dirackernels/dirac.i)
(test/tests/dirackernels/constant_point_source/3d_point_source.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
ny = 10
nz = 10
[]
[Variables]
active = 'u'
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[DiracKernels]
[./point_source1]
type = ConstantPointSource
variable = u
value = 0.1
point = '0.2 0.3 0.0'
[../]
[./point_source2]
type = ConstantPointSource
variable = u
value = -0.1
point = '0.2 0.8 0.0'
[../]
[./point_source3]
type = ConstantPointSource
variable = u
value = -1.0
point = '0.8 0.5 0.8'
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
file_base = 3d_out
exodus = true
[]
(test/tests/dirackernels/constant_point_source/1d_point_source.i)
###########################################################
# This is test of the Dirac delta function System. The
# ConstantPointSource object is used to apply a constant
# Dirac delta contribution at a specified point in the
# domain.
#
# @Requirement F3.50
###########################################################
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Variables]
active = 'u'
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[DiracKernels]
[./point_source1]
type = ConstantPointSource
variable = u
value = 1.0
point = '0.2 0 0'
[../]
[./point_source2]
type = ConstantPointSource
variable = u
value = -0.5
point = '0.7 0 0'
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
file_base = 1d_out
exodus = true
[]
(modules/porous_flow/test/tests/fluidstate/water_vapor_phasechange.i)
# Tests correct calculation of properties in PorousFlowWaterVapor as a phase change
# from liquid to a two-phase model occurs due to a pressure drop.
# A single 10 m^3 element is used, with constant mass and heat production using
# a Dirac kernel. Initial conditions correspond to just outside the two-phase region in
# the liquid state.
#
# An identical problem can be run using TOUGH2, with the following outputs after 1,000s
# Pressure: 8.58 Mpa
# Temperature: 299.92 K
# Vapor saturation: 0.00637
[Mesh]
type = GeneratedMesh
dim = 3
xmax = 10
ymax = 10
zmax = 10
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pliq]
initial_condition = 9e6
[]
[h]
scaling = 1e-3
[]
[]
[ICs]
[hic]
type = PorousFlowFluidPropertyIC
variable = h
porepressure = pliq
property = enthalpy
temperature = 300
temperature_unit = Celsius
fp = water
[]
[]
[DiracKernels]
[mass]
type = ConstantPointSource
point = '5 5 5'
variable = pliq
value = -1
[]
[heat]
type = ConstantPointSource
point = '5 5 5'
variable = h
value = -1.344269e6
[]
[]
[AuxVariables]
[pressure_gas]
order = CONSTANT
family = MONOMIAL
[]
[pressure_water]
order = CONSTANT
family = MONOMIAL
[]
[enthalpy_gas]
order = CONSTANT
family = MONOMIAL
[]
[enthalpy_water]
order = CONSTANT
family = MONOMIAL
[]
[saturation_gas]
order = CONSTANT
family = MONOMIAL
[]
[saturation_water]
order = CONSTANT
family = MONOMIAL
[]
[density_water]
order = CONSTANT
family = MONOMIAL
[]
[density_gas]
order = CONSTANT
family = MONOMIAL
[]
[viscosity_water]
order = CONSTANT
family = MONOMIAL
[]
[viscosity_gas]
order = CONSTANT
family = MONOMIAL
[]
[temperature]
order = CONSTANT
family = MONOMIAL
[]
[e_gas]
order = CONSTANT
family = MONOMIAL
[]
[e_water]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[enthalpy_water]
type = PorousFlowPropertyAux
variable = enthalpy_water
property = enthalpy
phase = 0
execute_on = 'initial timestep_end'
[]
[enthalpy_gas]
type = PorousFlowPropertyAux
variable = enthalpy_gas
property = enthalpy
phase = 1
execute_on = 'initial timestep_end'
[]
[pressure_water]
type = PorousFlowPropertyAux
variable = pressure_water
property = pressure
phase = 0
execute_on = 'initial timestep_end'
[]
[pressure_gas]
type = PorousFlowPropertyAux
variable = pressure_gas
property = pressure
phase = 1
execute_on = 'initial timestep_end'
[]
[saturation_water]
type = PorousFlowPropertyAux
variable = saturation_water
property = saturation
phase = 0
execute_on = 'initial timestep_end'
[]
[saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = 'initial timestep_end'
[]
[density_water]
type = PorousFlowPropertyAux
variable = density_water
property = density
phase = 0
execute_on = 'initial timestep_end'
[]
[density_gas]
type = PorousFlowPropertyAux
variable = density_gas
property = density
phase = 1
execute_on = 'initial timestep_end'
[]
[viscosity_water]
type = PorousFlowPropertyAux
variable = viscosity_water
property = viscosity
phase = 0
execute_on = 'initial timestep_end'
[]
[viscosity_gas]
type = PorousFlowPropertyAux
variable = viscosity_gas
property = viscosity
phase = 1
execute_on = 'initial timestep_end'
[]
[temperature]
type = PorousFlowPropertyAux
variable = temperature
property = temperature
execute_on = 'initial timestep_end'
[]
[e_water]
type = PorousFlowPropertyAux
variable = e_water
property = internal_energy
phase = 0
execute_on = 'initial timestep_end'
[]
[egas]
type = PorousFlowPropertyAux
variable = e_gas
property = internal_energy
phase = 1
execute_on = 'initial timestep_end'
[]
[]
[Kernels]
[mass]
type = PorousFlowMassTimeDerivative
variable = pliq
[]
[heat]
type = PorousFlowEnergyTimeDerivative
variable = h
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pliq h'
number_fluid_phases = 2
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureBC
pe = 1e5
lambda = 2
pc_max = 1e6
[]
[fs]
type = PorousFlowWaterVapor
water_fp = water
capillary_pressure = pc
[]
[]
[Modules]
[FluidProperties]
[water]
type = Water97FluidProperties
[]
[]
[]
[Materials]
[watervapor]
type = PorousFlowFluidStateSingleComponent
porepressure = pliq
enthalpy = h
temperature_unit = Celsius
capillary_pressure = pc
fluid_state = fs
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-14 0 0 0 1e-14 0 0 0 1e-14'
[]
[relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[]
[relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.2
[]
[internal_energy]
type = PorousFlowMatrixInternalEnergy
density = 2650
specific_heat_capacity = 1000
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 1e3
nl_abs_tol = 1e-12
[TimeStepper]
type = IterationAdaptiveDT
dt = 10
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Postprocessors]
[density_water]
type = ElementAverageValue
variable = density_water
execute_on = 'initial timestep_end'
[]
[density_gas]
type = ElementAverageValue
variable = density_gas
execute_on = 'initial timestep_end'
[]
[viscosity_water]
type = ElementAverageValue
variable = viscosity_water
execute_on = 'initial timestep_end'
[]
[viscosity_gas]
type = ElementAverageValue
variable = viscosity_gas
execute_on = 'initial timestep_end'
[]
[enthalpy_water]
type = ElementAverageValue
variable = enthalpy_water
execute_on = 'initial timestep_end'
[]
[enthalpy_gas]
type = ElementAverageValue
variable = enthalpy_gas
execute_on = 'initial timestep_end'
[]
[sg]
type = ElementAverageValue
variable = saturation_gas
execute_on = 'initial timestep_end'
[]
[sw]
type = ElementAverageValue
variable = saturation_water
execute_on = 'initial timestep_end'
[]
[pwater]
type = ElementAverageValue
variable = pressure_water
execute_on = 'initial timestep_end'
[]
[pgas]
type = ElementAverageValue
variable = pressure_gas
execute_on = 'initial timestep_end'
[]
[temperature]
type = ElementAverageValue
variable = temperature
execute_on = 'initial timestep_end'
[]
[enthalpy]
type = ElementAverageValue
variable = h
execute_on = 'initial timestep_end'
[]
[pliq]
type = ElementAverageValue
variable = pliq
execute_on = 'initial timestep_end'
[]
[liquid_mass]
type = PorousFlowFluidMass
phase = 0
execute_on = 'initial timestep_end'
[]
[vapor_mass]
type = PorousFlowFluidMass
phase = 1
execute_on = 'initial timestep_end'
[]
[liquid_heat]
type = PorousFlowHeatEnergy
phase = 0
execute_on = 'initial timestep_end'
[]
[vapor_heat]
type = PorousFlowHeatEnergy
phase = 1
execute_on = 'initial timestep_end'
[]
[e_water]
type = ElementAverageValue
variable = e_water
execute_on = 'initial timestep_end'
[]
[e_gas]
type = ElementAverageValue
variable = e_gas
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
csv = true
perf_graph = false
[]
(modules/geochemistry/test/tests/kernels/time_deriv_2.i)
# A point-source is added to fluid in a material with spatially-varying porosity
# porosity * d(concentration)/dt = 3.0 * delta(x - 1.0)
# where delta is the Dirac delta function (a ConstantPointSource DiracKernel)
# The solution, at x = 1.0 is
# concentration = concentration_old + 3 * dt / porosity
# while concentration is unchanged elsewhere.
# Note that since GeochemistryTimeDerivative is mass-lumped, it produces this solution.
# If mass lumping had not been used, concentration would have decreased at x != 1.0
[Mesh]
type = GeneratedMesh
dim = 1
nx = 2
xmax = 2
[]
[Variables]
[conc]
[]
[]
[Kernels]
[dot]
type = GeochemistryTimeDerivative
porosity = porosity
variable = conc
[]
[]
[DiracKernels]
[source]
type = ConstantPointSource
point = '1.0 0 0'
variable = conc
value = 12.0
[]
[]
[ICs]
[conc]
type = FunctionIC
function = 'x * x'
variable = conc
[]
[]
[AuxVariables]
[porosity]
[]
[expected]
[]
[should_be_zero]
[]
[]
[AuxKernels]
[porosity]
type = FunctionAux
function = '6.0 + x'
variable = porosity
[]
[expected]
type = FunctionAux
function = 'if(x > 0.5 & x < 1.5, x * x + 2.0 * 12.0 / (6.0 + x), x * x)'
variable = expected
[]
[should_be_zero]
type = ParsedAux
args = 'expected conc'
function = 'expected - conc'
variable = should_be_zero
[]
[]
[Postprocessors]
[error]
type = NodalL2Norm
variable = should_be_zero
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 2
end_time = 2
[]
[Outputs]
csv = true
[]
(test/tests/outputs/oversample/ex02_adapt.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
nz = 0
zmax = 0
elem_type = QUAD9
[]
[Variables]
[./diffused]
order = SECOND
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = diffused
[../]
[]
[DiracKernels]
[./foo]
variable = diffused
type = ConstantPointSource
value = 1
point = '0.3 0.3 0.0'
[../]
[]
[BCs]
[./all]
type = DirichletBC
variable = diffused
boundary = 'bottom left right top'
value = 0.0
[../]
[]
[Executioner]
type = Steady
solve_type = PJFNK
[]
[Adaptivity]
max_h_level = 2
initial_steps = 2
marker = marker
steps = 2
initial_marker = marker
[./Indicators]
[./indicator]
type = GradientJumpIndicator
variable = diffused
[../]
[../]
[./Markers]
[./marker]
type = ErrorFractionMarker
indicator = indicator
refine = 0.5
[../]
[../]
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[./os2]
type = Exodus
refinements = 2
[../]
[./os4]
type = Exodus
refinements = 4
[../]
[]
(examples/ex02_kernel/ex02_oversample.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
nz = 0
zmax = 0
elem_type = QUAD9
[]
[Variables]
[./diffused]
order = SECOND
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = diffused
[../]
[]
[DiracKernels]
[./foo]
variable = diffused
type = ConstantPointSource
value = 1
point = '0.3 0.3 0.0'
[../]
[]
[BCs]
[./all]
type = DirichletBC
variable = diffused
boundary = 'bottom left right top'
value = 0.0
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[./refine_2]
type = Exodus
file_base = oversample_2
refinements = 2
[../]
[./refine_4]
type = Exodus
file_base = oversample_4
refinements = 4
[../]
[]
(modules/tensor_mechanics/test/tests/shell/static/pinched_cylinder_symm.i)
# Test for displacement of pinched cylinder
# Ref: Figure 10 and Table 6 from Dvorkin and Bathe, Eng. Comput., Vol. 1, 1984.
# A cylinder of radius 1 m and length 2 m (along Z axis) with clamped ends
# (at z = 0 and 2 m) is pinched at mid-length by placing point loads of 10 N
# at (1, 0, 1) and (-1, 0, 1). Due to the symmetry of the problem, only 1/8th
# of the cylinder needs to be modeled.
# The normalized series solution for the displacement at the loading point is
# w = Wc E t / P = 164.24; where Wc is the displacement in m, E is the Young's
# modulus, t is the thickness and P is the point load.
# For this problem, E = 1e6 Pa, L = 2 m, R = 1 m, t = 0.01 m, P = 10 N and
# Poisson's ratio = 0.3. FEM results from different mesh discretizations are
# presented below. Only the 10x10 mesh is included as a test.
# Mesh of 1/8 cylinder | FEM/analytical (Moose) | FEM/analytical (Dvorkin)
# |ratio of normalized disp.| ratio of normalized disp.
#----------------------|-------------------------|-------------------------
# 10 x 10 | 0.806 | 0.83
# 20 x 20 | 1.06 | 0.96
# 40 x 40 | 0.95 | -
# 80 x 160 | 0.96 | -
# The results from FEM analysis matches well with the series solution and with
# the solution presented by Dvorkin and Bathe (1984).
[Mesh]
[./mesh]
type = FileMeshGenerator
file = pinched_cyl_10_10.msh
[../]
[./block_100]
type = ParsedSubdomainMeshGenerator
input = mesh
combinatorial_geometry = 'x > -1.1 & x < 1.1 & y > -1.1 & y < 1.1 & z > -0.1 & z < 2.1'
block_id = 100
[../]
[./nodeset_1]
type = BoundingBoxNodeSetGenerator
input = block_100
top_right = '1.1 1.1 0'
bottom_left = '-1.1 -1.1 0'
new_boundary = 'CD' #CD
[../]
[./nodeset_2]
type = BoundingBoxNodeSetGenerator
input = nodeset_1
top_right = '1.1 1.1 1.0'
bottom_left = '-1.1 -1.1 1.0'
new_boundary = 'AB' #AB
[../]
[./nodeset_3]
type = BoundingBoxNodeSetGenerator
input = nodeset_2
top_right = '0.02 1.1 1.0'
bottom_left = '-0.1 0.98 0.0'
new_boundary = 'AD' #AD
[../]
[./nodeset_4]
type = BoundingBoxNodeSetGenerator
input = nodeset_3
top_right = '1.1 0.02 1.0'
bottom_left = '0.98 -0.1 0.0'
new_boundary = 'BC' #BC
[../]
[]
[Variables]
[./disp_x]
order = FIRST
family = LAGRANGE
[../]
[./disp_y]
order = FIRST
family = LAGRANGE
[../]
[./disp_z]
order = FIRST
family = LAGRANGE
[../]
[./rot_x]
order = FIRST
family = LAGRANGE
[../]
[./rot_y]
order = FIRST
family = LAGRANGE
[../]
[]
[BCs]
[./simply_support_x]
type = DirichletBC
variable = disp_x
boundary = 'CD AD'
value = 0.0
[../]
[./simply_support_y]
type = DirichletBC
variable = disp_y
boundary = 'CD BC'
value = 0.0
[../]
[./simply_support_z]
type = DirichletBC
variable = disp_z
boundary = 'CD AB'
value = 0.0
[../]
[./simply_support_rot_x]
type = DirichletBC
variable = rot_x
boundary = 'CD BC'
value = 0.0
[../]
[./simply_support_rot_y]
type = DirichletBC
variable = rot_y
boundary = 'CD AD'
value = 0.0
[../]
[]
[DiracKernels]
[./point1]
type = ConstantPointSource
variable = disp_x
point = '1 0 1'
value = -2.5 # P = 10
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
line_search = 'none'
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
dt = 1.0
dtmin = 1.0
end_time = 1.0
[]
[Kernels]
[./solid_disp_x]
type = ADStressDivergenceShell
block = '100'
component = 0
variable = disp_x
through_thickness_order = SECOND
[../]
[./solid_disp_y]
type = ADStressDivergenceShell
block = '100'
component = 1
variable = disp_y
through_thickness_order = SECOND
[../]
[./solid_disp_z]
type = ADStressDivergenceShell
block = '100'
component = 2
variable = disp_z
through_thickness_order = SECOND
[../]
[./solid_rot_x]
type = ADStressDivergenceShell
block = '100'
component = 3
variable = rot_x
through_thickness_order = SECOND
[../]
[./solid_rot_y]
type = ADStressDivergenceShell
block = '100'
component = 4
variable = rot_y
through_thickness_order = SECOND
[../]
[]
[Materials]
[./elasticity]
type = ADComputeIsotropicElasticityTensorShell
youngs_modulus = 1e6
poissons_ratio = 0.3
block = '100'
through_thickness_order = SECOND
[../]
[./strain]
type = ADComputeIncrementalShellStrain
block = '100'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y'
thickness = 0.01
through_thickness_order = SECOND
[../]
[./stress]
type = ADComputeShellStress
block = '100'
through_thickness_order = SECOND
[../]
[]
[Postprocessors]
[./disp_z2]
type = PointValue
point = '1 0 1'
variable = disp_x
[../]
[]
[Outputs]
exodus = true
[]
(modules/tensor_mechanics/test/tests/truss/truss_hex.i)
# This test is designed to check
# whether truss element works well with other multi-dimensional element
# e.g. the hex element in this case, by assigning different brock number
# to different types of elements.
[Mesh]
type = FileMesh
file = truss_hex.e
displacements = 'disp_x disp_y disp_z'
[]
[Variables]
[./disp_x]
order = FIRST
family = LAGRANGE
[../]
[./disp_y]
order = FIRST
family = LAGRANGE
[../]
[./disp_z]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxVariables]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./e_over_l]
order = CONSTANT
family = MONOMIAL
[../]
[./area]
order = CONSTANT
family = MONOMIAL
# initial_condition = 1.0
[../]
[./react_x]
order = FIRST
family = LAGRANGE
[../]
[./react_y]
order = FIRST
family = LAGRANGE
[../]
[./react_z]
order = FIRST
family = LAGRANGE
[../]
[]
[Functions]
[./x2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 .5 1 1'
[../]
[./y2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 0 .5 1'
[../]
[]
[BCs]
[./fixx1]
type = DirichletBC
variable = disp_x
boundary = 1
value = 0.0
[../]
[./fixy1]
type = DirichletBC
variable = disp_y
boundary = 1
value = 0
[../]
[./fixz1]
type = DirichletBC
variable = disp_z
boundary = 1
value = 0
[../]
[./fixx2]
type = DirichletBC
variable = disp_x
boundary = 2
value = 0
[../]
[./fixz2]
type = DirichletBC
variable = disp_z
boundary = 2
value = 0
[../]
[./fixDummyHex_x]
type = DirichletBC
variable = disp_x
boundary = 1000
value = 0
[../]
[./fixDummyHex_y]
type = DirichletBC
variable = disp_y
boundary = 1000
value = 0
[../]
[./fixDummyHex_z]
type = DirichletBC
variable = disp_z
boundary = 1000
value = 0
[../]
[]
[DiracKernels]
[./pull]
type = ConstantPointSource
value = -25
point = '0 -2 0'
variable = disp_y
[../]
[]
[AuxKernels]
[./axial_stress]
type = MaterialRealAux
block = '1 2'
property = axial_stress
variable = axial_stress
[../]
[./e_over_l]
type = MaterialRealAux
block = '1 2'
property = e_over_l
variable = e_over_l
[../]
[./area1]
type = ConstantAux
block = 1
variable = area
value = 1.0
execute_on = 'initial timestep_begin'
[../]
[./area2]
type = ConstantAux
block = 2
variable = area
value = 0.25
execute_on = 'initial timestep_begin'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = PJFNK
petsc_options_iname = '-pc_type -ksp_gmres_restart'
petsc_options_value = 'jacobi 101'
nl_max_its = 15
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1
num_steps = 1
end_time = 1
[]
[Kernels]
[./truss_x]
type = StressDivergenceTensorsTruss
block = '1 2'
variable = disp_x
displacements = 'disp_x disp_y disp_z'
component = 0
area = area
save_in = react_x
[../]
[./truss_y]
type = StressDivergenceTensorsTruss
block = '1 2'
variable = disp_y
component = 1
displacements = 'disp_x disp_y disp_z'
area = area
save_in = react_y
[../]
[./truss_z]
type = StressDivergenceTensorsTruss
block = '1 2'
variable = disp_z
component = 2
displacements = 'disp_x disp_y disp_z'
area = area
save_in = react_z
[../]
[./TensorMechanics]
block = 1000
displacements = 'disp_x disp_y disp_z'
[../]
# [./hex_x]
# type = StressDivergenceTensors
# block = 1000
# variable = disp_x
# component = 0
# displacements = 'disp_x disp_y disp_z'
# [../]
# [./hex_y]
# type = StressDivergenceTensors
# block = 1000
# variable = disp_y
# component = 1
# displacements = 'disp_x disp_y disp_z'
# [../]
# [./hex_z]
# type = StressDivergenceTensors
# block = 1000
# variable = disp_z
# component = 2
# displacements = 'disp_x disp_y disp_z'
# [../]
[]
[Materials]
[./elasticity_tensor]
type = ComputeIsotropicElasticityTensor
block = 1000
youngs_modulus = 1e6
poissons_ratio = 0
[../]
[./strain]
type = ComputeSmallStrain
block = 1000
displacements = 'disp_x disp_y disp_z'
[../]
[./stress]
type = ComputeLinearElasticStress
block = 1000
[../]
[./linelast]
type = LinearElasticTruss
block = '1 2'
displacements = 'disp_x disp_y disp_z'
youngs_modulus = 1e6
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/dirackernels/constant_point_source/2d_point_source.i)
[Mesh]
[./square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[../]
uniform_refine = 4
[]
[Variables]
active = 'u'
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
active = 'diff'
[./diff]
type = Diffusion
variable = u
[../]
[]
[DiracKernels]
active = 'point_source1 point_source2'
[./point_source1]
type = ConstantPointSource
variable = u
value = 1.0
point = '0.2 0.3'
[../]
[./point_source2]
type = ConstantPointSource
variable = u
value = -0.5
point = '0.2 0.8'
[../]
[]
[BCs]
active = 'left right'
[./left]
type = DirichletBC
variable = u
boundary = 3
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = 1
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
file_base = 2d_out
exodus = true
[]
(modules/tensor_mechanics/test/tests/truss/truss_2d_action.i)
#
# Truss in two dimensional space
#
# The truss is made of five equilateral triangles supported at each end.
# The truss starts at (0,0). At (1,0), there is a point load of 25.
# The reactions are therefore
# Ryleft = 2/3 * 25 = 16.7
# Ryright = 1/3 * 25 = 8.33
# The area of each member is 0.8.
# Statics gives the stress in each member. For example, for element 6 (from
# (0,0) to (1/2,sqrt(3)/2)), the force is
# f = 2/3 * 25 * 2/sqrt(3) = 100/3/sqrt(3) (compressive)
# and the stress is
# s = -100/3/sqrt(3)/0.8 = -24.06
#
[Mesh]
type = FileMesh
file = truss_2d.e
displacements = 'disp_x disp_y'
[]
[AuxVariables]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./e_over_l]
order = CONSTANT
family = MONOMIAL
[../]
[./area]
order = CONSTANT
family = MONOMIAL
[../]
[./react_x]
order = FIRST
family = LAGRANGE
[../]
[./react_y]
order = FIRST
family = LAGRANGE
[../]
[./react_z]
order = FIRST
family = LAGRANGE
[../]
[]
[Functions]
[./x2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 .5 1 1'
[../]
[./y2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 0 .5 1'
[../]
[]
[BCs]
[./fixx1]
type = DirichletBC
variable = disp_x
boundary = 1
value = 0
[../]
[./fixy1]
type = DirichletBC
variable = disp_y
boundary = 1
value = 0
[../]
[./fixy4]
type = DirichletBC
variable = disp_y
boundary = 4
value = 0
[../]
[]
[DiracKernels]
[./pull]
type = ConstantPointSource
value = -25
point = '1 0 0'
variable = disp_y
[../]
[]
[AuxKernels]
[./axial_stress]
type = MaterialRealAux
block = 1
property = axial_stress
variable = axial_stress
[../]
[./e_over_l]
type = MaterialRealAux
block = 1
property = e_over_l
variable = e_over_l
[../]
[./area]
type = ConstantAux
block = 1
variable = area
value = 0.8
execute_on = 'initial timestep_begin'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = PJFNK
petsc_options_iname = '-pc_type -ksp_gmres_restart'
petsc_options_value = 'jacobi 101'
nl_max_its = 15
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1
num_steps = 1
end_time = 1
[]
[Modules/TensorMechanics/LineElementMaster]
[./block]
truss = true
add_variables = true
displacements = 'disp_x disp_y'
# area = area # commented out for error check
block = 1
save_in = 'react_x react_y'
[../]
[]
[Materials]
[./linelast]
type = LinearElasticTruss
block = 1
youngs_modulus = 1e6
displacements = 'disp_x disp_y'
[../]
[]
[Outputs]
file_base = 'truss_2d_out'
exodus = true
[]
(test/tests/dirackernels/constant_point_source/1d_point_source_fv.i)
[Mesh]
[mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 10
[]
[]
[Variables]
[u]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[]
[FVKernels]
[diff]
type = FVDiffusion
variable = u
coeff = coeff
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '1'
[]
[]
[DiracKernels]
[point_source1]
type = ConstantPointSource
variable = u
value = 1.0
point = '0.15 0 0'
[]
[point_source2]
type = ConstantPointSource
variable = u
value = -0.5
point = '0.65 0 0'
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = FVDirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
file_base = 1d_fv_out
exodus = true
[]
(modules/tensor_mechanics/test/tests/truss/truss_2d.i)
#
# Truss in two dimensional space
#
# The truss is made of five equilateral triangles supported at each end.
# The truss starts at (0,0). At (1,0), there is a point load of 25.
# The reactions are therefore
# Ryleft = 2/3 * 25 = 16.7
# Ryright = 1/3 * 25 = 8.33
# The area of each member is 0.8.
# Statics gives the stress in each member. For example, for element 6 (from
# (0,0) to (1/2,sqrt(3)/2)), the force is
# f = 2/3 * 25 * 2/sqrt(3) = 100/3/sqrt(3) (compressive)
# and the stress is
# s = -100/3/sqrt(3)/0.8 = -24.06
#
[Mesh]
type = FileMesh
file = truss_2d.e
displacements = 'disp_x disp_y'
[]
[Variables]
[./disp_x]
order = FIRST
family = LAGRANGE
[../]
[./disp_y]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxVariables]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./e_over_l]
order = CONSTANT
family = MONOMIAL
[../]
[./area]
order = CONSTANT
family = MONOMIAL
# initial_condition = 1.0
[../]
[./react_x]
order = FIRST
family = LAGRANGE
[../]
[./react_y]
order = FIRST
family = LAGRANGE
[../]
[./react_z]
order = FIRST
family = LAGRANGE
[../]
[]
[Functions]
[./x2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 .5 1 1'
[../]
[./y2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 0 .5 1'
[../]
[]
[BCs]
[./fixx1]
type = DirichletBC
variable = disp_x
boundary = 1
value = 0
[../]
[./fixy1]
type = DirichletBC
variable = disp_y
boundary = 1
value = 0
[../]
[./fixy4]
type = DirichletBC
variable = disp_y
boundary = 4
value = 0
[../]
[]
[DiracKernels]
[./pull]
type = ConstantPointSource
value = -25
point = '1 0 0'
variable = disp_y
[../]
[]
[AuxKernels]
[./axial_stress]
type = MaterialRealAux
block = 1
property = axial_stress
variable = axial_stress
[../]
[./e_over_l]
type = MaterialRealAux
block = 1
property = e_over_l
variable = e_over_l
[../]
[./area]
type = ConstantAux
block = 1
variable = area
value = 0.8
execute_on = 'initial timestep_begin'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = PJFNK
petsc_options_iname = '-pc_type -ksp_gmres_restart'
petsc_options_value = 'jacobi 101'
nl_max_its = 15
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1
num_steps = 1
end_time = 1
[]
[Kernels]
[./solid_x]
type = StressDivergenceTensorsTruss
block = 1
displacements = 'disp_x disp_y'
component = 0
variable = disp_x
area = area
save_in = react_x
[../]
[./solid_y]
type = StressDivergenceTensorsTruss
block = 1
displacements = 'disp_x disp_y'
component = 1
variable = disp_y
area = area
save_in = react_y
[../]
[]
[Materials]
[./linelast]
type = LinearElasticTruss
block = 1
youngs_modulus = 1e6
displacements = 'disp_x disp_y'
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/outputs/oversample/ex02_oversample.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
nz = 0
zmax = 0
elem_type = QUAD9
[]
[Variables]
[./diffused]
order = SECOND
[../]
[]
[Kernels]
active = 'diff'
[./diff]
type = Diffusion
variable = diffused
[../]
[]
[DiracKernels]
[./foo]
variable = diffused
type = ConstantPointSource
value = 1
point = '0.3 0.3 0.0'
[../]
[]
[BCs]
active = 'all'
[./all]
type = DirichletBC
variable = diffused
boundary = 'bottom left right top'
value = 0.0
[../]
[]
[Executioner]
type = Steady
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[./os2]
type = Exodus
refinements = 2
[../]
[./os4]
type = Exodus
refinements = 4
[../]
[]
(test/tests/dirackernels/constant_point_source/3d_point_source.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
ny = 10
nz = 10
[]
[Variables]
active = 'u'
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[DiracKernels]
[./point_source1]
type = ConstantPointSource
variable = u
value = 0.1
point = '0.2 0.3 0.0'
[../]
[./point_source2]
type = ConstantPointSource
variable = u
value = -0.1
point = '0.2 0.8 0.0'
[../]
[./point_source3]
type = ConstantPointSource
variable = u
value = -1.0
point = '0.8 0.5 0.8'
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
file_base = 3d_out
exodus = true
[]
(modules/tensor_mechanics/test/tests/truss/truss_hex_action.i)
# This test is designed to check
# whether truss element works well with other multi-dimensional element
# e.g. the hex element in this case, by assigning different block number
# to different types of elements.
[Mesh]
type = FileMesh
file = truss_hex.e
displacements = 'disp_x disp_y disp_z'
[]
[Variables]
[./disp_x]
order = FIRST
family = LAGRANGE
[../]
[./disp_y]
order = FIRST
family = LAGRANGE
[../]
[./disp_z]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxVariables]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./e_over_l]
order = CONSTANT
family = MONOMIAL
[../]
[./area]
order = CONSTANT
family = MONOMIAL
[../]
[./react_x]
order = FIRST
family = LAGRANGE
[../]
[./react_y]
order = FIRST
family = LAGRANGE
[../]
[./react_z]
order = FIRST
family = LAGRANGE
[../]
[]
[Functions]
[./x2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 .5 1 1'
[../]
[./y2]
type = PiecewiseLinear
x = '0 1 2 3'
y = '0 0 .5 1'
[../]
[]
[BCs]
[./fixx1]
type = DirichletBC
variable = disp_x
boundary = 1
value = 0.0
[../]
[./fixy1]
type = DirichletBC
variable = disp_y
boundary = 1
value = 0
[../]
[./fixz1]
type = DirichletBC
variable = disp_z
boundary = 1
value = 0
[../]
[./fixx2]
type = DirichletBC
variable = disp_x
boundary = 2
value = 0
[../]
[./fixz2]
type = DirichletBC
variable = disp_z
boundary = 2
value = 0
[../]
[./fixDummyHex_x]
type = DirichletBC
variable = disp_x
boundary = 1000
value = 0
[../]
[./fixDummyHex_y]
type = DirichletBC
variable = disp_y
boundary = 1000
value = 0
[../]
[./fixDummyHex_z]
type = DirichletBC
variable = disp_z
boundary = 1000
value = 0
[../]
[]
[DiracKernels]
[./pull]
type = ConstantPointSource
value = -25
point = '0 -2 0'
variable = disp_y
[../]
[]
[AuxKernels]
[./axial_stress]
type = MaterialRealAux
block = '1 2'
property = axial_stress
variable = axial_stress
[../]
[./e_over_l]
type = MaterialRealAux
block = '1 2'
property = e_over_l
variable = e_over_l
[../]
[./area1]
type = ConstantAux
block = 1
variable = area
value = 1.0
execute_on = 'initial timestep_begin'
[../]
[./area2]
type = ConstantAux
block = 2
variable = area
value = 0.25
execute_on = 'initial timestep_begin'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = PJFNK
petsc_options_iname = '-pc_type -ksp_gmres_restart'
petsc_options_value = 'jacobi 101'
nl_max_its = 15
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1
num_steps = 1
end_time = 1
[]
[Kernels]
[./TensorMechanics]
block = 1000
displacements = 'disp_x disp_y disp_z'
[../]
[]
[Modules/TensorMechanics/LineElementMaster]
[./block]
truss = true
displacements = 'disp_x disp_y disp_z'
area = area
block = '1 2'
save_in = 'react_x react_y react_z'
[../]
[]
[Materials]
[./elasticity_tensor]
type = ComputeIsotropicElasticityTensor
block = 1000
youngs_modulus = 1e6
poissons_ratio = 0
[../]
[./strain]
type = ComputeSmallStrain
block = 1000
displacements = 'disp_x disp_y disp_z'
[../]
[./stress]
type = ComputeLinearElasticStress
block = 1000
[../]
[./linelast]
type = LinearElasticTruss
block = '1 2'
displacements = 'disp_x disp_y disp_z'
youngs_modulus = 1e6
[../]
[]
[Outputs]
file_base = 'truss_hex_out'
exodus = true
[]
(test/tests/controls/time_periods/dirackernels/dirac.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
uniform_refine = 2
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = CoefDiffusion
variable = u
coef = 0.5
[../]
[./time]
type = TimeDerivative
variable = u
[../]
[]
[DiracKernels]
[./point_source]
type = ConstantPointSource
variable = u
value = 1
point = '0.25 0.25'
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 0.1
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
[Controls]
[./point_source]
type = TimePeriod
disable_objects = 'DiracKernel::point_source'
start_time = '0.15'
end_time = '0.35'
execute_on = 'initial timestep_begin'
[../]
[]