- long_term_viscosityViscosity corresponding to the long-term part of the deformation
C++ Type:double
Description:Viscosity corresponding to the long-term part of the deformation
- poissons_ratioInitial poisson ratio of the material
C++ Type:double
Description:Initial poisson ratio of the material
- youngs_modulusInitial elastic modulus of the material
C++ Type:double
Description:Initial elastic modulus of the material
ConcreteLogarithmicCreepModel
Logarithmic viscoelastic model for cementitious materials.
Description
This material represents a logarithmic viscoelastic behavior for concrete and other cementitious materials.
The strain is decomposed into three components: an elastic strain , a recoverable viscoelastic strain and an irrecoverable creep strain . This corresponds to a Burgers type of material with an elastic spring, a Kelvin-Voigt module, and a dashpot placed in series.
The implementation uses the linear viscoelasticity framework defined in MOOSE tensor_mechanics module.
Elastic Strain
The elastic strain is directly related to the stress with the usual Hooke's law, where is the fourth-order elasticity tensor of the material. In the current model, is calculated using the youngs_modulus
and poissons_ratio
parameters provided by the user.
Irrecoverable Creep Strain
This strain corresponds to the long-term visco-elastic strain of the material. It is linear with the logarithm of time, and cannot be recovered upon unloading. It is calculated using an aging dashpot: is the same elasticity tensor used for the calculation of the elastic strain.
is defined using the long_term_viscosity
parameter, and controls the slope of the creep curve in the logarithmic space.
is defined using the long_term_characteristic_time
parameter (default value of 1 day) and controls when the logarithmic behavior starts.
Recoverable Creep Strain
This strain corresponds to the short-term visco-elastic strain of the material, which can be partially recovered upon unloading. It is calculated using a single Kelvin-Voigt module, where is the fourth-order elasticity tensor of the spring and the characteristic time of the dashpot.
is calculated using the recoverable_youngs_modulus
and recoverable_poissons_ratio
parameters provided by the user, or using youngs_modulus
and poissons_ratio
parameters if they are undefined.
is defined using the recoverable_viscosity
parameter, or with the long_term_characteristic_time
parameter if it is undefined.
To disable the recoverable strain, the user must set the use_recovery
parameter to FALSE.
Creep At Elevated Temperature
The model can include the effect on temperature on creep.
The creep strain rates are increased using an Arrhenius-type law, with the activation_temperature
parameter controlling the increase of strain rate with temperature.
The model assumes that temperature
and reference_temperature
are in Celsius, while activation_temperature
is in Kelvin.
Creep Under Low Relative Humidity
The model can include the effect of a constant low relative humidity on creep.
In this case, the irrecoverable creep strain rate is decreased linearly with the current relative humidity.
This effect is automatically activated when the user defines a humidity coupled variable.
Drying Creep
The model can include the effect of drying on the concrete creep, also known as the Pickett effect.
In this case, the characteristic time of the irrecoverable creep strain is modified to account for the increase of the irrecoverable creep strain rate on drying. The new characteristic time depends on , with the drying_creep_viscosity
parameter provided by the user, and the rate of drying.
This effect if accounted in the simulation when the user defines a humidity coupled variable, and .
!!! The calculation of the drying creep requires either humidity
to be computed on timestep_begin
(for example using an AuxKernel or a MultiApp transfer), or the force_recompute_properties
parameter to be set to TRUE
.
Time Integration
The visco-elastic strains are integrated in time using a time-stepping procedure defined for linear viscoelastic materials in the tensor_mechanics module.
For the time-stepping scheme to be properly updated, a LinearViscoelasticityManager object must be included in the input file, and linked to the ConcreteLogarithmicCreepModel material
The ConcreteLogarithmicCreepModel is compatible with either the total small strain approximation, or either of the incremental strain approximation (incremental small strains or finite strains). It needs the following stress calculators:
Strain | Stress | Additional Materials |
---|---|---|
ComputeSmallStrain | ComputeLinearViscoelasticStress | none |
ComputeIncrementalSmallStrain | ComputeMultipleInelasticStress | LinearViscoelasticStressUpdate |
ComputeFiniteStrain | ComputeMultipleInelasticStress | LinearViscoelasticStressUpdate |
The stress calculators use the actual elasticity tensor of the material , which is provided by the ConcreteLogarithmicCreepModel itself.
Input Parameters
- activation_temperatureActivation temperature for the creep [in Kelvin]
C++ Type:double
Options:
Description:Activation temperature for the creep [in Kelvin]
- base_nameOptional parameter that allows the user to define multiple mechanics material systems on the same block, i.e. for multiple phases
C++ Type:std::string
Options:
Description:Optional parameter that allows the user to define multiple mechanics material systems on the same block, i.e. for multiple phases
- blockThe list of block ids (SubdomainID) that this object will be applied
C++ Type:std::vector<SubdomainName>
Options:
Description:The list of block ids (SubdomainID) that this object will be applied
- boundaryThe list of boundary IDs from the mesh where this boundary condition applies
C++ Type:std::vector<BoundaryName>
Options:
Description:The list of boundary IDs from the mesh where this boundary condition applies
- computeTrueWhen false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
Default:True
C++ Type:bool
Options:
Description:When false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
- constant_onNONEWhen ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
Default:NONE
C++ Type:MooseEnum
Options:NONE, ELEMENT, SUBDOMAIN
Description:When ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
- creep_strain_namecreep_strainname of the true creep strain of the material(computed by LinearViscoelasticStressUpdate orComputeLinearViscoelasticStress)
Default:creep_strain
C++ Type:std::string
Options:
Description:name of the true creep strain of the material(computed by LinearViscoelasticStressUpdate orComputeLinearViscoelasticStress)
- driving_eigenstrainname of the eigenstrain that increases the creep strains
C++ Type:std::string
Options:
Description:name of the eigenstrain that increases the creep strains
- drying_creep_viscosityViscosity corresponding to the drying creep
C++ Type:double
Options:
Description:Viscosity corresponding to the drying creep
- elastic_strain_nameelastic_strainname of the true elastic strain of the material
Default:elastic_strain
C++ Type:std::string
Options:
Description:name of the true elastic strain of the material
- force_recompute_propertiesFalseforces the computation of the viscoelastic properties at each step ofthe solver (default: false)
Default:False
C++ Type:bool
Options:
Description:forces the computation of the viscoelastic properties at each step ofthe solver (default: false)
- humidityHumidity variable
C++ Type:std::vector<VariableName>
Options:
Description:Humidity variable
- integration_rulebackward-eulerdescribes how the viscoelastic behavior is integrated through time
Default:backward-euler
C++ Type:MooseEnum
Options:backward-euler, mid-point, newmark, zienkiewicz
Description:describes how the viscoelastic behavior is integrated through time
- long_term_characteristic_time1Rate at which the long_term viscosity increases
Default:1
C++ Type:double
Options:
Description:Rate at which the long_term viscosity increases
- recoverable_poissons_ratioPoisson coefficient of the recoverable part of the deformation
C++ Type:double
Options:
Description:Poisson coefficient of the recoverable part of the deformation
- recoverable_viscosityViscosity corresponding to the recoverable part of the deformation
C++ Type:double
Options:
Description:Viscosity corresponding to the recoverable part of the deformation
- recoverable_youngs_modulusModulus corresponding to the recoverable part of the deformation
C++ Type:double
Options:
Description:Modulus corresponding to the recoverable part of the deformation
- reference_temperature20Reference temperature [in Celsius]
Default:20
C++ Type:double
Options:
Description:Reference temperature [in Celsius]
- temperatureTemperature variable [in Celsius]
C++ Type:std::vector<VariableName>
Options:
Description:Temperature variable [in Celsius]
- theta1coefficient for Newmark integration rule (between 0 and 1)
Default:1
C++ Type:double
Options:
Description:coefficient for Newmark integration rule (between 0 and 1)
- use_recoveryTrueEnables or disables creep recovery
Default:True
C++ Type:bool
Options:
Description:Enables or disables creep recovery
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
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
Options:
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
Options:
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
Options:
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
- output_propertiesList of material properties, from this material, to output (outputs must also be defined to an output type)
C++ Type:std::vector<std::string>
Options:
Description:List of material properties, from this material, to output (outputs must also be defined to an output type)
- outputsnone Vector of output names were you would like to restrict the output of variables(s) associated with this object
Default:none
C++ Type:std::vector<OutputName>
Options:
Description:Vector of output names were you would like to restrict the output of variables(s) associated with this object
Outputs Parameters
Input Files
- (assessment/asr_validation/wald2017b/analysis/A1-unreinforced.i)
- (assessment/asr_validation/wald2017b/analysis/A1-biaxial.i)
- (test/tests/concrete_logcreep/concrete_logcreep_norec.i)
- (assessment/asr_validation/wald2017b/analysis/A3-biaxial.i)
- (test/tests/concrete_logcreep/concrete_logcreep.i)
- (assessment/asr_validation/wald2017b/analysis/A1-triaxial.i)
- (test/tests/concrete_logcreep/concrete_logcreep_temp.i)
- (test/tests/concrete_logcreep/concrete_logcreep_humidity.i)
- (test/tests/concrete_logcreep/concrete_logcreep_drying.i)
- (assessment/asr_validation/wald2017b/analysis/A1-uniaxial.i)
(assessment/asr_validation/wald2017b/analysis/A1-unreinforced.i)
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
block = 1
[]
[Problem]
type = ReferenceResidualProblem
reference_vector = 'ref'
extra_tag_vectors = 'ref'
group_variables = 'disp_x disp_y disp_z'
[]
[Mesh]
file = A1-unreinforced.e
[]
[Variables]
[./T]
order = FIRST
family = LAGRANGE
initial_condition = 10.6
[../]
[./rh]
order = FIRST
family = LAGRANGE
initial_condition = 0.8
block = 1
[../]
[]
[AuxVariables]
[./resid_x]
[../]
[./resid_y]
[../]
[./resid_z]
[../]
[./ASR_ex]
order = CONSTANT
family = MONOMIAL
[../]
[./ASR_vstrain]
order = CONSTANT
family = MONOMIAL
[../]
[./ASR_strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./ASR_strain_yy]
order = CONSTANT
family = MONOMIAL
[../]
[./ASR_strain_zz]
order = CONSTANT
family = MONOMIAL
[../]
[./ASR_strain_xy]
order = CONSTANT
family = MONOMIAL
[../]
[./ASR_strain_yz]
order = CONSTANT
family = MONOMIAL
[../]
[./ASR_strain_zx]
order = CONSTANT
family = MONOMIAL
[../]
[./volumetric_strain]
order = CONSTANT
family = MONOMIAL
[../]
[./thermal_strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./thermal_strain_yy]
order = CONSTANT
family = MONOMIAL
[../]
[./thermal_strain_zz]
order = CONSTANT
family = MONOMIAL
[../]
[./thermal_conductivity]
order = CONSTANT
family = Monomial
[../]
[./thermal_capacity]
order = CONSTANT
family = Monomial
[../]
[./moisture_capacity]
order = CONSTANT
family = Monomial
[../]
[./humidity_diffusivity]
order = CONSTANT
family = Monomial
[../]
[./water_content]
order = CONSTANT
family = Monomial
[../]
[./water_hydrated]
order = CONSTANT
family = Monomial
[../]
[damage_index]
order = CONSTANT
family = MONOMIAL
[]
[]
[Modules/TensorMechanics/Master]
[./concrete]
strain = FINITE
add_variables = true
eigenstrain_names = 'asr_expansion thermal_expansion'
generate_output = 'stress_xx stress_yy stress_zz stress_xy stress_yz stress_zx vonmises_stress hydrostatic_stress elastic_strain_xx elastic_strain_yy elastic_strain_zz strain_xx strain_yy strain_zz'
save_in = 'resid_x resid_y resid_z'
[../]
[]
[Kernels]
[./T_td]
type = ConcreteThermalTimeIntegration
variable = T
extra_vector_tags = 'ref'
[../]
[./T_diff]
type = ConcreteThermalConduction
variable = T
extra_vector_tags = 'ref'
[../]
[./T_conv]
type = ConcreteThermalConvection
variable = T
relative_humidity = rh
extra_vector_tags = 'ref'
[../]
[./T_adsorption]
type = ConcreteLatentHeat
variable = T
H = rh
extra_vector_tags = 'ref'
[../]
[./rh_td]
type = ConcreteMoistureTimeIntegration
variable = rh
extra_vector_tags = 'ref'
[../]
[./rh_diff]
type = ConcreteMoistureDiffusion
variable = rh
temperature = T
extra_vector_tags = 'ref'
[../]
[]
[AuxKernels]
[./ASR_ex]
type = MaterialRealAux
variable = ASR_ex
property = ASR_extent
execute_on = 'timestep_end'
[../]
[./ASR_vstrain]
type = MaterialRealAux
variable = ASR_vstrain
property = ASR_volumetric_strain
execute_on = 'timestep_end'
[../]
[./ASR_strain_xx]
type = RankTwoAux
rank_two_tensor = asr_expansion
variable = ASR_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./ASR_strain_yy]
type = RankTwoAux
rank_two_tensor = asr_expansion
variable = ASR_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_zz]
type = RankTwoAux
rank_two_tensor = asr_expansion
variable = ASR_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_xy]
type = RankTwoAux
rank_two_tensor = asr_expansion
variable = ASR_strain_xy
index_i = 0
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_yz]
type = RankTwoAux
rank_two_tensor = asr_expansion
variable = ASR_strain_yz
index_i = 1
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_zx]
type = RankTwoAux
rank_two_tensor = asr_expansion
variable = ASR_strain_zx
index_i = 0
index_j = 2
execute_on = 'timestep_end'
[../]
[./thermal_strain_xx]
type = RankTwoAux
rank_two_tensor = thermal_expansion
variable = thermal_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./thermal_strain_yy]
type = RankTwoAux
rank_two_tensor = thermal_expansion
variable = thermal_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./thermal_strain_zz]
type = RankTwoAux
rank_two_tensor = thermal_expansion
variable = thermal_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./volumetric_strain]
type = RankTwoScalarAux
scalar_type = VolumetricStrain
rank_two_tensor = total_strain
variable = volumetric_strain
[../]
[./k]
type = MaterialRealAux
variable = thermal_conductivity
property = thermal_conductivity
execute_on = 'timestep_end'
[../]
[./capacity]
type = MaterialRealAux
variable = thermal_capacity
property = thermal_capacity
execute_on = 'timestep_end'
[../]
[./rh_capacity]
type = MaterialRealAux
variable = moisture_capacity
property = moisture_capacity
execute_on = 'timestep_end'
[../]
[./rh_duff]
type = MaterialRealAux
variable = humidity_diffusivity
property = humidity_diffusivity
execute_on = 'timestep_end'
[../]
[./wc_duff]
type = MaterialRealAux
variable = water_content
property = moisture_content
execute_on = 'timestep_end'
[../]
[./hydrw_duff]
type = MaterialRealAux
variable = water_hydrated
property = hydrated_water
execute_on = 'timestep_end'
[../]
[damage_index]
type = MaterialRealAux
variable = damage_index
property = damage_index
execute_on = timestep_end
[]
[]
[Functions]
[./ramp_temp]
type = PiecewiseLinear
data_file = temperature_history.csv
format = columns
[../]
[./ramp_humidity]
type = PiecewiseLinear
data_file = humidity_history.csv
format = columns
[../]
[]
[Materials]
[./concrete]
type = ConcreteThermalMoisture
# setup thermal property models and parameters
# options available: CONSTANT ASCE-1992 KODUR-2004 EUROCODE-2004 KIM-2003
thermal_conductivity_model = KODUR-2004
thermal_capacity_model = KODUR-2004
aggregate_type = Siliceous #options: Siliceous Carbonate
ref_density_of_concrete = 2231.0 # in kg/m^3
ref_specific_heat_of_concrete = 1100.0 # in J/(Kg.0C)
ref_thermal_conductivity_of_concrete = 3 # in W/(m.0C)
# setup moisture capacity and humidity diffusivity models
aggregate_pore_type = dense #options: dense porous
aggregate_mass = 1877.0 #mass of aggregate (kg) per m^3 of concrete
cement_type = 1 #options: 1 2 3 4
cement_mass = 354.0 #mass of cement (kg) per m^3 of concrete
water_to_cement_ratio = 0.5
concrete_cure_time = 28.0 #curing time in (days)
# options available for humidity diffusivity:
moisture_diffusivity_model = Bazant #options: Bazant Mensi
D1 = 3.0e-8
coupled_moisture_diffusivity_factor = 1.0e-2 # factor for mositure diffusivity due to heat
# coupled nonlinear variables
relative_humidity = rh
temperature = T
[../]
[./creep]
type = LinearViscoelasticStressUpdate
# block = 1
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
poissons_ratio = 0.22
youngs_modulus = 37.3e9
recoverable_youngs_modulus = 37.3e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
humidity = rh
temperature = T
activation_temperature = 23.0
[../]
[ASR_expansion]
type = ConcreteASREigenstrain
expansion_type = Anisotropic
reference_temperature = 35.0
temperature_unit = Celsius
max_volumetric_expansion = 2.2e-2
characteristic_time = 18.9
latency_time = 18.0
characteristic_activation_energy = 5400.0
latency_activation_energy = 9400.0
stress_latency_factor = 1.0
compressive_strength = 31.0e6
compressive_stress_exponent = 0.0
expansion_stress_limit = 8.0e6
tensile_strength = 3.2e6
tensile_retention_factor = 1.0
tensile_absorption_factor = 1.0
ASR_dependent_tensile_strength = false
residual_tensile_strength_fraction = 1.0
temperature = T
relative_humidity = rh
rh_exponent = 1.0
eigenstrain_name = asr_expansion
absolute_tolerance = 1e-10
output_iteration_info_on_error = true
[]
[thermal_strain_concrete]
type = ComputeThermalExpansionEigenstrain
temperature = T
thermal_expansion_coeff = 8.0e-6
stress_free_temperature = 10.6
eigenstrain_name = thermal_expansion
[]
[ASR_damage_concrete]
type = ConcreteASRMicrocrackingDamage
residual_youngs_modulus_fraction = 0.1
[]
[./stress]
type = ComputeMultipleInelasticStress
inelastic_models = 'creep'
damage_model = ASR_damage_concrete
[../]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
viscoelastic_model = logcreep
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = disp_x
boundary = '2000 2005'
value = 0.0
[../]
[./bottom]
type = DirichletBC
variable = disp_y
boundary = '2000 2001'
value = 0.0
[../]
[./back]
type = DirichletBC
variable = disp_z
boundary = '2000 2005'
value = 0.0
[../]
[./T]
type = FunctionDirichletBC
variable = T
boundary = '101 102 103 104 105 106'
function = ramp_temp
[../]
[./rh]
type = FunctionDirichletBC
variable = rh
boundary = '101 102 103 104 105 106'
function = ramp_humidity
[../]
[]
[Postprocessors]
[./nelem]
type = NumElems
[../]
[./ndof]
type = NumDOFs
[../]
[./ASR_strain]
type = ElementAverageValue
variable = ASR_vstrain
[../]
[./ASR_strain_xx]
type = ElementAverageValue
variable = ASR_strain_xx
[../]
[./ASR_strain_yy]
type = ElementAverageValue
variable = ASR_strain_yy
[../]
[./ASR_strain_zz]
type = ElementAverageValue
variable = ASR_strain_zz
[../]
[ASR_ext]
type = ElementAverageValue
variable = ASR_ex
[]
[./vonmises]
type = ElementAverageValue
variable = vonmises_stress
[../]
[./vstrain]
type = ElementAverageValue
variable = volumetric_strain
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
[../]
[./strain_yy]
type = ElementAverageValue
variable = strain_yy
[../]
[./strain_zz]
type = ElementAverageValue
variable = strain_zz
[../]
[./temp]
type = ElementAverageValue
variable = T
[../]
[./temp_bc]
type = SideAverageValue
variable = T
boundary = 102
[../]
[./humidity]
type = ElementAverageValue
variable = rh
[../]
[./humidity_bc]
type = SideAverageValue
variable = rh
boundary = 102
[../]
[./thermal_strain_xx]
type = ElementAverageValue
variable = thermal_strain_xx
[../]
[./thermal_strain_yy]
type = ElementAverageValue
variable = thermal_strain_yy
[../]
[./thermal_strain_zz]
type = ElementAverageValue
variable = thermal_strain_zz
[../]
[./disp_x_101]
type = SideAverageValue
variable = disp_x
boundary = 101
[../]
[./disp_x_102]
type = SideAverageValue
variable = disp_x
boundary = 102
[../]
[./disp_x_103]
type = SideAverageValue
variable = disp_x
boundary = 103
[../]
[./disp_x_104]
type = SideAverageValue
variable = disp_x
boundary = 104
[../]
[./disp_x_105]
type = SideAverageValue
variable = disp_x
boundary = 105
[../]
[./disp_x_106]
type = SideAverageValue
variable = disp_x
boundary = 106
[../]
[./disp_y_101]
type = SideAverageValue
variable = disp_y
boundary = 101
[../]
[./disp_y_102]
type = SideAverageValue
variable = disp_y
boundary = 102
[../]
[./disp_y_103]
type = SideAverageValue
variable = disp_y
boundary = 103
[../]
[./disp_y_104]
type = SideAverageValue
variable = disp_y
boundary = 104
[../]
[./disp_y_105]
type = SideAverageValue
variable = disp_y
boundary = 105
[../]
[./disp_y_106]
type = SideAverageValue
variable = disp_y
boundary = 106
[../]
[./disp_z_101]
type = SideAverageValue
variable = disp_z
boundary = 101
[../]
[./disp_z_102]
type = SideAverageValue
variable = disp_z
boundary = 102
[../]
[./disp_z_103]
type = SideAverageValue
variable = disp_z
boundary = 103
[../]
[./disp_z_104]
type = SideAverageValue
variable = disp_z
boundary = 104
[../]
[./disp_z_105]
type = SideAverageValue
variable = disp_z
boundary = 105
[../]
[./disp_z_106]
type = SideAverageValue
variable = disp_z
boundary = 106
[../]
[disp_x_p1_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 -0.08'
[../]
[disp_x_p1_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 -0.08'
[../]
[disp_x_p2_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 0.08'
[../]
[disp_x_p2_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 0.08'
[../]
[disp_x_p3_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.08'
[../]
[disp_x_p3_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.08'
[../]
[disp_x_p4_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.08'
[../]
[disp_x_p4_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.08'
[../]
[disp_x_p5_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.235'
[../]
[disp_x_p5_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.235'
[../]
[disp_x_p6_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.235'
[../]
[disp_x_p6_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.235'
[../]
[disp_y_p1_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 -0.08'
[../]
[disp_y_p1_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 -0.08'
[../]
[disp_y_p2_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 0.08'
[../]
[disp_y_p2_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 0.08'
[../]
[disp_y_p3_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.08'
[../]
[disp_y_p3_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.08'
[../]
[disp_y_p4_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.08'
[../]
[disp_y_p4_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.08'
[../]
[disp_y_p5_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.235'
[../]
[disp_y_p5_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.235'
[../]
[disp_y_p6_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.235'
[../]
[disp_y_p6_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.235'
[../]
[disp_y_p7_pos]
type = PointValue
variable = disp_y
point = '-0.235 0.24 0.08'
[../]
[disp_y_p7_neg]
type = PointValue
variable = disp_y
point = '-0.235 -0.24 0.08'
[../]
[disp_y_p8_pos]
type = PointValue
variable = disp_y
point = '0.235 0.24 0.08'
[../]
[disp_y_p8_neg]
type = PointValue
variable = disp_y
point = '0.235 -0.24 0.08'
[../]
[disp_z_p1_pos]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 0.24'
[../]
[disp_z_p1_neg]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 -0.24'
[../]
[disp_z_p2_pos]
type = PointValue
variable = disp_z
point = '-0.08 0.08 0.24'
[../]
[disp_z_p2_neg]
type = PointValue
variable = disp_z
point = '-0.08 0.08 -0.24'
[../]
[disp_z_p3_pos]
type = PointValue
variable = disp_z
point = '0.08 -0.08 0.24'
[../]
[disp_z_p3_neg]
type = PointValue
variable = disp_z
point = '0.08 -0.08 -0.24'
[../]
[disp_z_p4_pos]
type = PointValue
variable = disp_z
point = '0.08 0.08 0.24'
[../]
[disp_z_p4_neg]
type = PointValue
variable = disp_z
point = '0.08 0.08 -0.24'
[../]
[disp_z_p5_pos]
type = PointValue
variable = disp_z
point = '0.235 0.08 0.24'
[../]
[disp_z_p5_neg]
type = PointValue
variable = disp_z
point = '0.235 0.08 -0.24'
[../]
[disp_z_p6_pos]
type = PointValue
variable = disp_z
point = '-0.235 0.08 0.24'
[../]
[disp_z_p6_neg]
type = PointValue
variable = disp_z
point = '-0.235 0.08 -0.24'
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
line_search = none
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
start_time = 2419200 #28 days
dt = 86400 #1 day in sec
automatic_scaling = true
end_time = 38880000 #450 days
l_max_its = 20
nl_max_its = 10
nl_rel_tol = 1e-6
# Because this problem is unrestrained, the displacement reference is 0,
# so this controls the displacement convergence:
nl_abs_tol = 1e-10
[]
[Outputs]
perf_graph = true
csv = true
#exodus = true #Turned off to save space
[]
[Debug]
show_var_residual_norms = true
[]
(assessment/asr_validation/wald2017b/analysis/A1-biaxial.i)
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
penalty = 1e12
[]
[Problem]
type = ReferenceResidualProblem
reference_vector = 'ref'
extra_tag_vectors = 'ref'
group_variables = 'disp_x disp_y disp_z'
[]
[Mesh]
file = A1-biaxial.e
[]
[Variables]
[./T]
order = FIRST
family = LAGRANGE
initial_condition = 10.6
[../]
[./rh]
order = FIRST
family = LAGRANGE
initial_condition = 0.8
block = 1
[../]
[]
[AuxVariables]
[./resid_x]
[../]
[./resid_y]
[../]
[./resid_z]
[../]
[./ASR_ex]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_vstrain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./volumetric_strain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_conductivity]
order = CONSTANT
family = Monomial
[../]
[./thermal_capacity]
order = CONSTANT
family = Monomial
[../]
[./moisture_capacity]
order = CONSTANT
family = Monomial
[../]
[./humidity_diffusivity]
order = CONSTANT
family = Monomial
[../]
[./water_content]
order = CONSTANT
family = Monomial
[../]
[./water_hydrated]
order = CONSTANT
family = Monomial
[../]
[damage_index]
order = CONSTANT
family = MONOMIAL
[]
[./area]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_strain]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Modules/TensorMechanics/Master]
[./concrete]
block = 1
strain = FINITE
add_variables = true
eigenstrain_names = 'asr_expansion thermal_expansion'
generate_output = 'stress_xx stress_yy stress_zz stress_xy stress_yz stress_zx vonmises_stress hydrostatic_stress elastic_strain_xx elastic_strain_yy elastic_strain_zz strain_xx strain_yy strain_zz'
extra_vector_tags = 'ref'
[../]
[]
[Modules/TensorMechanics/LineElementMaster]
[./Reinforcement_block]
block = '2 3'
truss = true
area = area
displacements = 'disp_x disp_y disp_z'
#Note: Intentially not including this here to have it give a nonzero
# displacement reference residual since it's an unrestrained problem
#extra_vector_tags = 'ref'
[../]
[]
[Constraints]
[./rebar_x2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[./rebar_x3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[]
[Kernels]
[./T_td]
type = ConcreteThermalTimeIntegration
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_diff]
type = ConcreteThermalConduction
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_conv]
type = ConcreteThermalConvection
variable = T
relative_humidity = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./T_adsorption]
type = ConcreteLatentHeat
variable = T
H = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_td]
type = ConcreteMoistureTimeIntegration
variable = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_diff]
type = ConcreteMoistureDiffusion
variable = rh
temperature = T
block = 1
extra_vector_tags = 'ref'
[../]
[./heat_dt]
type = TimeDerivative
variable = T
block = '2 3'
extra_vector_tags = 'ref'
[../]
[./heat_conduction]
type = HeatConduction
variable = T
diffusion_coefficient = 53.0
block = '2 3'
extra_vector_tags = 'ref'
[../]
[]
[AuxKernels]
[./ASR_ex]
type = MaterialRealAux
variable = ASR_ex
block = 1
property = ASR_extent
execute_on = 'timestep_end'
[../]
[./ASR_vstrain]
type = MaterialRealAux
block = 1
variable = ASR_vstrain
property = ASR_volumetric_strain
execute_on = 'timestep_end'
[../]
[./ASR_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./ASR_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_xy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xy
index_i = 0
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_yz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yz
index_i = 1
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_zx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zx
index_i = 0
index_j = 2
execute_on = 'timestep_end'
[../]
[./thermal_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./thermal_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./thermal_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./volumetric_strain]
type = RankTwoScalarAux
scalar_type = VolumetricStrain
rank_two_tensor = total_strain
variable = volumetric_strain
block = 1
[../]
[./k]
type = MaterialRealAux
variable = thermal_conductivity
property = thermal_conductivity
execute_on = 'timestep_end'
block = 1
[../]
[./capacity]
type = MaterialRealAux
variable = thermal_capacity
property = thermal_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_capacity]
type = MaterialRealAux
variable = moisture_capacity
property = moisture_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_duff]
type = MaterialRealAux
variable = humidity_diffusivity
property = humidity_diffusivity
execute_on = 'timestep_end'
block = 1
[../]
[./wc_duff]
type = MaterialRealAux
variable = water_content
property = moisture_content
execute_on = 'timestep_end'
block = 1
[../]
[./hydrw_duff]
type = MaterialRealAux
variable = water_hydrated
property = hydrated_water
execute_on = 'timestep_end'
block = 1
[../]
[damage_index]
type = MaterialRealAux
block = 1
variable = damage_index
property = damage_index
execute_on = timestep_end
[]
[./areax]
type = ConstantAux
block = '2'
variable = area
value = 1.33e-4
execute_on = 'initial timestep_begin'
[../]
[./areaz]
type = ConstantAux
block = '3'
variable = area
value = 1.33e-4
execute_on = 'initial timestep_begin'
[../]
[./axial_stress]
type = MaterialRealAux
block = '2 3'
variable = axial_stress
property = axial_stress
[../]
[]
[Functions]
[./ramp_temp]
type = PiecewiseLinear
data_file = temperature_history.csv
format = columns
[../]
[./ramp_humidity]
type = PiecewiseLinear
data_file = humidity_history.csv
format = columns
[../]
[]
[Materials]
[./concrete]
type = ConcreteThermalMoisture
block = 1
# setup thermal property models and parameters
# options available: CONSTANT ASCE-1992 KODUR-2004 EUROCODE-2004 KIM-2003
thermal_conductivity_model = KODUR-2004
thermal_capacity_model = KODUR-2004
aggregate_type = Siliceous #options: Siliceous Carbonate
ref_density_of_concrete = 2231.0 # in kg/m^3
ref_specific_heat_of_concrete = 1100.0 # in J/(Kg.0C)
ref_thermal_conductivity_of_concrete = 3 # in W/(m.0C)
# setup moisture capacity and humidity diffusivity models
aggregate_pore_type = dense #options: dense porous
aggregate_mass = 1877.0 #mass of aggregate (kg) per m^3 of concrete
cement_type = 1 #options: 1 2 3 4
cement_mass = 354.0 #mass of cement (kg) per m^3 of concrete
water_to_cement_ratio = 0.5
concrete_cure_time = 28.0 #curing time in (days)
# options available for humidity diffusivity:
moisture_diffusivity_model = Bazant #options: Bazant Mensi
D1 = 3.0e-8
coupled_moisture_diffusivity_factor = 1.0e-2 # factor for mositure diffusivity due to heat
# coupled nonlinear variables
relative_humidity = rh
temperature = T
[../]
[./creep]
type = LinearViscoelasticStressUpdate
block = 1
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
block = 1
poissons_ratio = 0.22
youngs_modulus = 37.3e9
recoverable_youngs_modulus = 37.3e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
humidity = rh
temperature = T
activation_temperature = 23.0
[../]
[ASR_expansion]
type = ConcreteASREigenstrain
block = 1
expansion_type = Anisotropic
reference_temperature = 35.0
temperature_unit = Celsius
max_volumetric_expansion = 2.2e-2
characteristic_time = 18.9
latency_time = 18.0
characteristic_activation_energy = 5400.0
latency_activation_energy = 9400.0
stress_latency_factor = 1.0
compressive_strength = 31.0e6
compressive_stress_exponent = 0.0
expansion_stress_limit = 8.0e6
tensile_strength = 3.2e6
tensile_retention_factor = 1.0
tensile_absorption_factor = 1.0
ASR_dependent_tensile_strength = false
residual_tensile_strength_fraction = 1.0
temperature = T
relative_humidity = rh
rh_exponent = 1.0
eigenstrain_name = asr_expansion
absolute_tolerance = 1e-10
output_iteration_info_on_error = true
[]
[thermal_strain_concrete]
type = ComputeThermalExpansionEigenstrain
block = 1
temperature = T
thermal_expansion_coeff = 8.0e-6
stress_free_temperature = 10.6
eigenstrain_name = thermal_expansion
[]
[ASR_damage_concrete]
type = ConcreteASRMicrocrackingDamage
residual_youngs_modulus_fraction = 0.1
block = 1
[]
[./stress]
type = ComputeMultipleInelasticStress
block = 1
inelastic_models = 'creep'
damage_model = ASR_damage_concrete
[../]
[truss]
type = LinearElasticTruss
block = '2 3'
youngs_modulus = 2e11
temperature = T
thermal_expansion_coeff = 11.3e-6
temperature_ref = 10.6
[]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
block = 1
viscoelastic_model = logcreep
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = disp_x
boundary = '2000 2005'
value = 0.0
[../]
[./bottom]
type = DirichletBC
variable = disp_y
boundary = '2000 2001'
value = 0.0
[../]
[./back]
type = DirichletBC
variable = disp_z
boundary = '2000 2005'
value = 0.0
[../]
[./T]
type = FunctionDirichletBC
variable = T
boundary = '101 102 103 104 105 106'
function = ramp_temp
[../]
[./rh]
type = FunctionDirichletBC
variable = rh
boundary = '101 102 103 104 105 106'
function = ramp_humidity
[../]
[]
[Postprocessors]
[./nelem]
type = NumElems
[../]
[./ndof]
type = NumDOFs
[../]
[./ASR_strain]
type = ElementAverageValue
variable = ASR_vstrain
block = 1
[../]
[./ASR_strain_xx]
type = ElementAverageValue
variable = ASR_strain_xx
block = 1
[../]
[./ASR_strain_yy]
type = ElementAverageValue
variable = ASR_strain_yy
block = 1
[../]
[./ASR_strain_zz]
type = ElementAverageValue
variable = ASR_strain_zz
block = 1
[../]
[ASR_ext]
type = ElementAverageValue
variable = ASR_ex
block = 1
[]
[./vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 1
[../]
[./vstrain]
type = ElementAverageValue
variable = volumetric_strain
block = 1
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 1
[../]
[./strain_yy]
type = ElementAverageValue
variable = strain_yy
block = 1
[../]
[./strain_zz]
type = ElementAverageValue
variable = strain_zz
block = 1
[../]
[./temp]
type = ElementAverageValue
variable = T
block = 1
[../]
[./humidity]
type = ElementAverageValue
variable = rh
block = 1
[../]
[./thermal_strain_xx]
type = ElementAverageValue
variable = thermal_strain_xx
block = 1
[../]
[./thermal_strain_yy]
type = ElementAverageValue
variable = thermal_strain_yy
block = 1
[../]
[./thermal_strain_zz]
type = ElementAverageValue
variable = thermal_strain_zz
block = 1
[../]
[./disp_x_101]
type = SideAverageValue
variable = disp_x
boundary = 101
[../]
[./disp_x_102]
type = SideAverageValue
variable = disp_x
boundary = 102
[../]
[./disp_x_103]
type = SideAverageValue
variable = disp_x
boundary = 103
[../]
[./disp_x_104]
type = SideAverageValue
variable = disp_x
boundary = 104
[../]
[./disp_x_105]
type = SideAverageValue
variable = disp_x
boundary = 105
[../]
[./disp_x_106]
type = SideAverageValue
variable = disp_x
boundary = 106
[../]
[./disp_y_101]
type = SideAverageValue
variable = disp_y
boundary = 101
[../]
[./disp_y_102]
type = SideAverageValue
variable = disp_y
boundary = 102
[../]
[./disp_y_103]
type = SideAverageValue
variable = disp_y
boundary = 103
[../]
[./disp_y_104]
type = SideAverageValue
variable = disp_y
boundary = 104
[../]
[./disp_y_105]
type = SideAverageValue
variable = disp_y
boundary = 105
[../]
[./disp_y_106]
type = SideAverageValue
variable = disp_y
boundary = 106
[../]
[./disp_z_101]
type = SideAverageValue
variable = disp_z
boundary = 101
[../]
[./disp_z_102]
type = SideAverageValue
variable = disp_z
boundary = 102
[../]
[./disp_z_103]
type = SideAverageValue
variable = disp_z
boundary = 103
[../]
[./disp_z_104]
type = SideAverageValue
variable = disp_z
boundary = 104
[../]
[./disp_z_105]
type = SideAverageValue
variable = disp_z
boundary = 105
[../]
[./disp_z_106]
type = SideAverageValue
variable = disp_z
boundary = 106
[../]
[disp_x_p1_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 -0.08'
[../]
[disp_x_p1_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 -0.08'
[../]
[disp_x_p2_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 0.08'
[../]
[disp_x_p2_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 0.08'
[../]
[disp_x_p3_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.08'
[../]
[disp_x_p3_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.08'
[../]
[disp_x_p4_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.08'
[../]
[disp_x_p4_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.08'
[../]
[disp_x_p5_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.235'
[../]
[disp_x_p5_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.235'
[../]
[disp_x_p6_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.235'
[../]
[disp_x_p6_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.235'
[../]
[disp_y_p1_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 -0.08'
[../]
[disp_y_p1_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 -0.08'
[../]
[disp_y_p2_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 0.08'
[../]
[disp_y_p2_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 0.08'
[../]
[disp_y_p3_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.08'
[../]
[disp_y_p3_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.08'
[../]
[disp_y_p4_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.08'
[../]
[disp_y_p4_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.08'
[../]
[disp_y_p5_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.235'
[../]
[disp_y_p5_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.235'
[../]
[disp_y_p6_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.235'
[../]
[disp_y_p6_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.235'
[../]
[disp_y_p7_pos]
type = PointValue
variable = disp_y
point = '-0.235 0.24 0.08'
[../]
[disp_y_p7_neg]
type = PointValue
variable = disp_y
point = '-0.235 -0.24 0.08'
[../]
[disp_y_p8_pos]
type = PointValue
variable = disp_y
point = '0.235 0.24 0.08'
[../]
[disp_y_p8_neg]
type = PointValue
variable = disp_y
point = '0.235 -0.24 0.08'
[../]
[disp_z_p1_pos]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 0.24'
[../]
[disp_z_p1_neg]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 -0.24'
[../]
[disp_z_p2_pos]
type = PointValue
variable = disp_z
point = '-0.08 0.08 0.24'
[../]
[disp_z_p2_neg]
type = PointValue
variable = disp_z
point = '-0.08 0.08 -0.24'
[../]
[disp_z_p3_pos]
type = PointValue
variable = disp_z
point = '0.08 -0.08 0.24'
[../]
[disp_z_p3_neg]
type = PointValue
variable = disp_z
point = '0.08 -0.08 -0.24'
[../]
[disp_z_p4_pos]
type = PointValue
variable = disp_z
point = '0.08 0.08 0.24'
[../]
[disp_z_p4_neg]
type = PointValue
variable = disp_z
point = '0.08 0.08 -0.24'
[../]
[disp_z_p5_pos]
type = PointValue
variable = disp_z
point = '0.235 0.08 0.24'
[../]
[disp_z_p5_neg]
type = PointValue
variable = disp_z
point = '0.235 0.08 -0.24'
[../]
[disp_z_p6_pos]
type = PointValue
variable = disp_z
point = '-0.235 0.08 0.24'
[../]
[disp_z_p6_neg]
type = PointValue
variable = disp_z
point = '-0.235 0.08 -0.24'
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
line_search = none
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
start_time = 2419200
dt = 1000000
automatic_scaling = true
end_time = 38880000
l_max_its = 20
nl_max_its = 10
nl_rel_tol = 1e-6
nl_abs_tol = 1e-10
[]
[Outputs]
perf_graph = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
(test/tests/concrete_logcreep/concrete_logcreep_norec.i)
# Tests the ConcreteLogarithmicCreepModel without recoverable creep
#
# Expected value for strain_xx
# 0.001 * (1 + log(1 + t))
#
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
elem_type = HEX8
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]
[./stress_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./creep_strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
[./TensorMechanics]
displacements = 'disp_x disp_y disp_z'
use_displaced_mesh = true
[../]
[]
[AuxKernels]
[./stress_xx]
type = RankTwoAux
variable = stress_xx
rank_two_tensor = stress
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./strain_xx]
type = RankTwoAux
variable = strain_xx
rank_two_tensor = total_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./creep_strain_xx]
type = RankTwoAux
variable = creep_strain_xx
rank_two_tensor = creep_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[]
[BCs]
[./symmy]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0
[../]
[./symmx]
type = DirichletBC
variable = disp_x
boundary = left
value = 0
[../]
[./symmz]
type = DirichletBC
variable = disp_z
boundary = back
value = 0
[../]
[./axial_load]
type = NeumannBC
variable = disp_x
boundary = right
value = 10e6
[../]
[]
[Materials]
[./stress]
type = ComputeLinearViscoelasticStress
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
poissons_ratio = 0.2
youngs_modulus = 10e9
long_term_viscosity = 1
long_term_characteristic_time = 1
use_recovery = false
[../]
[./strain]
type = ComputeSmallStrain
displacements = 'disp_x disp_y disp_z'
[../]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
viscoelastic_model = logcreep
[../]
[]
[Postprocessors]
[./stress_xx]
type = ElementAverageValue
variable = stress_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./creep_strain_xx]
type = ElementAverageValue
variable = creep_strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
l_max_its = 50
l_tol = 1e-8
nl_max_its = 20
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
dtmin = 0.01
end_time = 100
[./TimeStepper]
type = LogConstantDT
log_dt = 0.1
first_dt = 0.1
[../]
[]
[Outputs]
file_base = concrete_logcreep_norec_out
exodus = true
csv = true
[]
(assessment/asr_validation/wald2017b/analysis/A3-biaxial.i)
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
penalty = 1e12
[]
[Problem]
type = ReferenceResidualProblem
reference_vector = 'ref'
extra_tag_vectors = 'ref'
group_variables = 'disp_x disp_y disp_z'
[]
[Mesh]
file = A3-biaxial.e
[]
[Variables]
[./T]
order = FIRST
family = LAGRANGE
initial_condition = 10.6
[../]
[./rh]
order = FIRST
family = LAGRANGE
initial_condition = 0.8
block = 1
[../]
[]
[AuxVariables]
[./resid_x]
[../]
[./resid_y]
[../]
[./resid_z]
[../]
[./ASR_ex]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_vstrain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./volumetric_strain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_conductivity]
order = CONSTANT
family = Monomial
[../]
[./thermal_capacity]
order = CONSTANT
family = Monomial
[../]
[./moisture_capacity]
order = CONSTANT
family = Monomial
[../]
[./humidity_diffusivity]
order = CONSTANT
family = Monomial
[../]
[./water_content]
order = CONSTANT
family = Monomial
[../]
[./water_hydrated]
order = CONSTANT
family = Monomial
[../]
[damage_index]
order = CONSTANT
family = MONOMIAL
[]
[./area]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_strain]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Modules/TensorMechanics/Master]
[./concrete]
block = 1
strain = FINITE
add_variables = true
eigenstrain_names = 'asr_expansion thermal_expansion'
generate_output = 'stress_xx stress_yy stress_zz stress_xy stress_yz stress_zx vonmises_stress hydrostatic_stress elastic_strain_xx elastic_strain_yy elastic_strain_zz strain_xx strain_yy strain_zz'
extra_vector_tags = 'ref'
[../]
[]
[Modules/TensorMechanics/LineElementMaster]
[./Reinforcement_block]
block = '2 3'
truss = true
area = area
displacements = 'disp_x disp_y disp_z'
#Note: Intentially not including this here to have it give a nonzero
# displacement reference residual since it's an unrestrained problem
#extra_vector_tags = 'ref'
[../]
[]
[Constraints]
[./rebar_x2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[./rebar_x3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[]
[Kernels]
[./T_td]
type = ConcreteThermalTimeIntegration
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_diff]
type = ConcreteThermalConduction
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_conv]
type = ConcreteThermalConvection
variable = T
relative_humidity = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./T_adsorption]
type = ConcreteLatentHeat
variable = T
H = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_td]
type = ConcreteMoistureTimeIntegration
variable = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_diff]
type = ConcreteMoistureDiffusion
variable = rh
temperature = T
block = 1
extra_vector_tags = 'ref'
[../]
[./heat_dt]
type = TimeDerivative
variable = T
block = '2 3'
extra_vector_tags = 'ref'
[../]
[./heat_conduction]
type = HeatConduction
variable = T
diffusion_coefficient = 53.0
block = '2 3'
extra_vector_tags = 'ref'
[../]
[]
[AuxKernels]
[./ASR_ex]
type = MaterialRealAux
variable = ASR_ex
block = 1
property = ASR_extent
execute_on = 'timestep_end'
[../]
[./ASR_vstrain]
type = MaterialRealAux
block = 1
variable = ASR_vstrain
property = ASR_volumetric_strain
execute_on = 'timestep_end'
[../]
[./ASR_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./ASR_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_xy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xy
index_i = 0
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_yz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yz
index_i = 1
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_zx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zx
index_i = 0
index_j = 2
execute_on = 'timestep_end'
[../]
[./thermal_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./thermal_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./thermal_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./volumetric_strain]
type = RankTwoScalarAux
scalar_type = VolumetricStrain
rank_two_tensor = total_strain
variable = volumetric_strain
block = 1
[../]
[./k]
type = MaterialRealAux
variable = thermal_conductivity
property = thermal_conductivity
execute_on = 'timestep_end'
block = 1
[../]
[./capacity]
type = MaterialRealAux
variable = thermal_capacity
property = thermal_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_capacity]
type = MaterialRealAux
variable = moisture_capacity
property = moisture_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_duff]
type = MaterialRealAux
variable = humidity_diffusivity
property = humidity_diffusivity
execute_on = 'timestep_end'
block = 1
[../]
[./wc_duff]
type = MaterialRealAux
variable = water_content
property = moisture_content
execute_on = 'timestep_end'
block = 1
[../]
[./hydrw_duff]
type = MaterialRealAux
variable = water_hydrated
property = hydrated_water
execute_on = 'timestep_end'
block = 1
[../]
[damage_index]
type = MaterialRealAux
block = 1
variable = damage_index
property = damage_index
execute_on = timestep_end
[]
[./areax]
type = ConstantAux
block = '2'
variable = area
value = 1.33e-4
execute_on = 'initial timestep_begin'
[../]
[./areaz]
type = ConstantAux
block = '3'
variable = area
value = 2.84e-4
execute_on = 'initial timestep_begin'
[../]
[./axial_stress]
type = MaterialRealAux
block = '2 3'
variable = axial_stress
property = axial_stress
[../]
[]
[Functions]
[./ramp_temp]
type = PiecewiseLinear
data_file = temperature_history.csv
format = columns
[../]
[./ramp_humidity]
type = PiecewiseLinear
data_file = humidity_history.csv
format = columns
[../]
[]
[Materials]
[./concrete]
type = ConcreteThermalMoisture
block = 1
# setup thermal property models and parameters
# options available: CONSTANT ASCE-1992 KODUR-2004 EUROCODE-2004 KIM-2003
thermal_conductivity_model = KODUR-2004
thermal_capacity_model = KODUR-2004
aggregate_type = Siliceous #options: Siliceous Carbonate
ref_density_of_concrete = 2231.0 # in kg/m^3
ref_specific_heat_of_concrete = 1100.0 # in J/(Kg.0C)
ref_thermal_conductivity_of_concrete = 3 # in W/(m.0C)
# setup moisture capacity and humidity diffusivity models
aggregate_pore_type = dense #options: dense porous
aggregate_mass = 1877.0 #mass of aggregate (kg) per m^3 of concrete
cement_type = 1 #options: 1 2 3 4
cement_mass = 354.0 #mass of cement (kg) per m^3 of concrete
water_to_cement_ratio = 0.5
concrete_cure_time = 28.0 #curing time in (days)
# options available for humidity diffusivity:
moisture_diffusivity_model = Bazant #options: Bazant Mensi
D1 = 3.0e-8
coupled_moisture_diffusivity_factor = 1.0e-2 # factor for mositure diffusivity due to heat
# coupled nonlinear variables
relative_humidity = rh
temperature = T
[../]
[./creep]
type = LinearViscoelasticStressUpdate
block = 1
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
block = 1
poissons_ratio = 0.22
youngs_modulus = 37.3e9
recoverable_youngs_modulus = 37.3e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
humidity = rh
temperature = T
activation_temperature = 23.0
[../]
[ASR_expansion]
type = ConcreteASREigenstrain
block = 1
expansion_type = Anisotropic
reference_temperature = 35.0
temperature_unit = Celsius
max_volumetric_expansion = 2.2e-2
characteristic_time = 18.9
latency_time = 18.0
characteristic_activation_energy = 5400.0
latency_activation_energy = 9400.0
stress_latency_factor = 1.0
compressive_strength = 31.0e6
compressive_stress_exponent = 0.0
expansion_stress_limit = 8.0e6
tensile_strength = 3.2e6
tensile_retention_factor = 1.0
tensile_absorption_factor = 1.0
ASR_dependent_tensile_strength = false
residual_tensile_strength_fraction = 1.0
temperature = T
relative_humidity = rh
rh_exponent = 1.0
eigenstrain_name = asr_expansion
absolute_tolerance = 1e-10
output_iteration_info_on_error = true
[]
[thermal_strain_concrete]
type = ComputeThermalExpansionEigenstrain
block = 1
temperature = T
thermal_expansion_coeff = 8.0e-6
stress_free_temperature = 10.6
eigenstrain_name = thermal_expansion
[]
[ASR_damage_concrete]
type = ConcreteASRMicrocrackingDamage
residual_youngs_modulus_fraction = 0.1
block = 1
[]
[./stress]
type = ComputeMultipleInelasticStress
block = 1
inelastic_models = 'creep'
damage_model = ASR_damage_concrete
[../]
[truss]
type = LinearElasticTruss
block = '2 3'
youngs_modulus = 2e11
temperature = T
thermal_expansion_coeff = 11.3e-6
temperature_ref = 10.6
[]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
block = 1
viscoelastic_model = logcreep
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = disp_x
boundary = '2000 2005'
value = 0.0
[../]
[./bottom]
type = DirichletBC
variable = disp_y
boundary = '2000 2001'
value = 0.0
[../]
[./back]
type = DirichletBC
variable = disp_z
boundary = '2000 2005'
value = 0.0
[../]
[./T]
type = FunctionDirichletBC
variable = T
boundary = '101 102 103 104 105 106'
function = ramp_temp
[../]
[./rh]
type = FunctionDirichletBC
variable = rh
boundary = '101 102 103 104 105 106'
function = ramp_humidity
[../]
[]
[Postprocessors]
[./nelem]
type = NumElems
[../]
[./ndof]
type = NumDOFs
[../]
[./ASR_strain]
type = ElementAverageValue
variable = ASR_vstrain
block = 1
[../]
[./ASR_strain_xx]
type = ElementAverageValue
variable = ASR_strain_xx
block = 1
[../]
[./ASR_strain_yy]
type = ElementAverageValue
variable = ASR_strain_yy
block = 1
[../]
[./ASR_strain_zz]
type = ElementAverageValue
variable = ASR_strain_zz
block = 1
[../]
[ASR_ext]
type = ElementAverageValue
variable = ASR_ex
block = 1
[]
[./vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 1
[../]
[./vstrain]
type = ElementAverageValue
variable = volumetric_strain
block = 1
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 1
[../]
[./strain_yy]
type = ElementAverageValue
variable = strain_yy
block = 1
[../]
[./strain_zz]
type = ElementAverageValue
variable = strain_zz
block = 1
[../]
[./temp]
type = ElementAverageValue
variable = T
block = 1
[../]
[./humidity]
type = ElementAverageValue
variable = rh
block = 1
[../]
[./thermal_strain_xx]
type = ElementAverageValue
variable = thermal_strain_xx
block = 1
[../]
[./thermal_strain_yy]
type = ElementAverageValue
variable = thermal_strain_yy
block = 1
[../]
[./thermal_strain_zz]
type = ElementAverageValue
variable = thermal_strain_zz
block = 1
[../]
[./disp_x_101]
type = SideAverageValue
variable = disp_x
boundary = 101
[../]
[./disp_x_102]
type = SideAverageValue
variable = disp_x
boundary = 102
[../]
[./disp_x_103]
type = SideAverageValue
variable = disp_x
boundary = 103
[../]
[./disp_x_104]
type = SideAverageValue
variable = disp_x
boundary = 104
[../]
[./disp_x_105]
type = SideAverageValue
variable = disp_x
boundary = 105
[../]
[./disp_x_106]
type = SideAverageValue
variable = disp_x
boundary = 106
[../]
[./disp_y_101]
type = SideAverageValue
variable = disp_y
boundary = 101
[../]
[./disp_y_102]
type = SideAverageValue
variable = disp_y
boundary = 102
[../]
[./disp_y_103]
type = SideAverageValue
variable = disp_y
boundary = 103
[../]
[./disp_y_104]
type = SideAverageValue
variable = disp_y
boundary = 104
[../]
[./disp_y_105]
type = SideAverageValue
variable = disp_y
boundary = 105
[../]
[./disp_y_106]
type = SideAverageValue
variable = disp_y
boundary = 106
[../]
[./disp_z_101]
type = SideAverageValue
variable = disp_z
boundary = 101
[../]
[./disp_z_102]
type = SideAverageValue
variable = disp_z
boundary = 102
[../]
[./disp_z_103]
type = SideAverageValue
variable = disp_z
boundary = 103
[../]
[./disp_z_104]
type = SideAverageValue
variable = disp_z
boundary = 104
[../]
[./disp_z_105]
type = SideAverageValue
variable = disp_z
boundary = 105
[../]
[./disp_z_106]
type = SideAverageValue
variable = disp_z
boundary = 106
[../]
[disp_x_p1_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 -0.08'
[../]
[disp_x_p1_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 -0.08'
[../]
[disp_x_p2_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 0.08'
[../]
[disp_x_p2_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 0.08'
[../]
[disp_x_p3_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.08'
[../]
[disp_x_p3_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.08'
[../]
[disp_x_p4_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.08'
[../]
[disp_x_p4_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.08'
[../]
[disp_x_p5_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.235'
[../]
[disp_x_p5_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.235'
[../]
[disp_x_p6_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.235'
[../]
[disp_x_p6_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.235'
[../]
[disp_x_p7_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 -0.235'
[../]
[disp_x_p7_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 -0.235'
[../]
[disp_x_p8_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 0.235'
[../]
[disp_x_p8_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 0.235'
[../]
[disp_x_p9_pos]
type = PointValue
variable = disp_x
point = '0.24 0.235 -0.235'
[../]
[disp_x_p9_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.235 -0.235'
[../]
[disp_x_p10_pos]
type = PointValue
variable = disp_x
point = '0.24 0.235 0.235'
[../]
[disp_x_p10_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.235 0.235'
[../]
[disp_x_p11_pos]
type = PointValue
variable = disp_x
point = '0.24 0.235 -0.08'
[../]
[disp_x_p11_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.235 -0.08'
[../]
[disp_x_p12_pos]
type = PointValue
variable = disp_x
point = '0.24 0.235 0.08'
[../]
[disp_x_p12_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.235 0.08'
[../]
[disp_x_p13_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.235 -0.235'
[../]
[disp_x_p13_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.235 -0.235'
[../]
[disp_x_p14_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.235 0.235'
[../]
[disp_x_p14_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.235 0.235'
[../]
[disp_x_p15_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.235 -0.08'
[../]
[disp_x_p15_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.235 -0.08'
[../]
[disp_x_p16_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.235 0.08'
[../]
[disp_x_p16_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.235 0.08'
[../]
[disp_y_p1_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 -0.08'
[../]
[disp_y_p1_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 -0.08'
[../]
[disp_y_p2_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 0.08'
[../]
[disp_y_p2_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 0.08'
[../]
[disp_y_p3_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.08'
[../]
[disp_y_p3_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.08'
[../]
[disp_y_p4_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.08'
[../]
[disp_y_p4_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.08'
[../]
[disp_y_p5_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.235'
[../]
[disp_y_p5_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.235'
[../]
[disp_y_p6_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.235'
[../]
[disp_y_p6_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.235'
[../]
[disp_y_p7_pos]
type = PointValue
variable = disp_y
point = '-0.235 0.24 0.08'
[../]
[disp_y_p7_neg]
type = PointValue
variable = disp_y
point = '-0.235 -0.24 0.08'
[../]
[disp_y_p8_pos]
type = PointValue
variable = disp_y
point = '0.235 0.24 0.08'
[../]
[disp_y_p8_neg]
type = PointValue
variable = disp_y
point = '0.235 -0.24 0.08'
[../]
[disp_y_p9_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 -0.235'
[../]
[disp_y_p9_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 -0.235'
[../]
[disp_y_p10_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 0.235'
[../]
[disp_y_p10_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 0.235'
[../]
[disp_y_p11_pos]
type = PointValue
variable = disp_y
point = '-0.235 0.24 -0.08'
[../]
[disp_y_p11_neg]
type = PointValue
variable = disp_y
point = '-0.235 -0.24 -0.08'
[../]
[disp_y_p12_pos]
type = PointValue
variable = disp_y
point = '0.235 0.24 -0.08'
[../]
[disp_y_p12_neg]
type = PointValue
variable = disp_y
point = '0.235 -0.24 -0.08'
[../]
[disp_z_p1_pos]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 0.24'
[../]
[disp_z_p1_neg]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 -0.24'
[../]
[disp_z_p2_pos]
type = PointValue
variable = disp_z
point = '-0.08 0.08 0.24'
[../]
[disp_z_p2_neg]
type = PointValue
variable = disp_z
point = '-0.08 0.08 -0.24'
[../]
[disp_z_p3_pos]
type = PointValue
variable = disp_z
point = '0.08 -0.08 0.24'
[../]
[disp_z_p3_neg]
type = PointValue
variable = disp_z
point = '0.08 -0.08 -0.24'
[../]
[disp_z_p4_pos]
type = PointValue
variable = disp_z
point = '0.08 0.08 0.24'
[../]
[disp_z_p4_neg]
type = PointValue
variable = disp_z
point = '0.08 0.08 -0.24'
[../]
[disp_z_p5_pos]
type = PointValue
variable = disp_z
point = '0.235 0.08 0.24'
[../]
[disp_z_p5_neg]
type = PointValue
variable = disp_z
point = '0.235 0.08 -0.24'
[../]
[disp_z_p6_pos]
type = PointValue
variable = disp_z
point = '-0.235 0.08 0.24'
[../]
[disp_z_p6_neg]
type = PointValue
variable = disp_z
point = '-0.235 0.08 -0.24'
[../]
[disp_z_p7_pos]
type = PointValue
variable = disp_z
point = '-0.235 -0.08 0.24'
[../]
[disp_z_p7_neg]
type = PointValue
variable = disp_z
point = '-0.235 -0.08 -0.24'
[../]
[disp_z_p8_pos]
type = PointValue
variable = disp_z
point = '0.235 -0.08 0.24'
[../]
[disp_z_p8_neg]
type = PointValue
variable = disp_z
point = '0.235 -0.08 -0.24'
[../]
[disp_z_p9_pos]
type = PointValue
variable = disp_z
point = '0.235 0.235 0.24'
[../]
[disp_z_p9_neg]
type = PointValue
variable = disp_z
point = '0.235 0.235 -0.24'
[../]
[disp_z_p10_pos]
type = PointValue
variable = disp_z
point = '-0.235 0.235 0.24'
[../]
[disp_z_p10_neg]
type = PointValue
variable = disp_z
point = '-0.235 0.235 -0.24'
[../]
[disp_z_p11_pos]
type = PointValue
variable = disp_z
point = '0.08 0.235 0.24'
[../]
[disp_z_p11_neg]
type = PointValue
variable = disp_z
point = '0.08 0.235 -0.24'
[../]
[disp_z_p12_pos]
type = PointValue
variable = disp_z
point = '-0.08 0.235 0.24'
[../]
[disp_z_p12_neg]
type = PointValue
variable = disp_z
point = '-0.08 0.235 -0.24'
[../]
[disp_z_p13_pos]
type = PointValue
variable = disp_z
point = '0.235 -0.235 0.24'
[../]
[disp_z_p13_neg]
type = PointValue
variable = disp_z
point = '0.235 -0.235 -0.24'
[../]
[disp_z_p14_pos]
type = PointValue
variable = disp_z
point = '-0.235 -0.235 0.24'
[../]
[disp_z_p14_neg]
type = PointValue
variable = disp_z
point = '-0.235 -0.235 -0.24'
[../]
[disp_z_p15_pos]
type = PointValue
variable = disp_z
point = '0.08 -0.235 0.24'
[../]
[disp_z_p15_neg]
type = PointValue
variable = disp_z
point = '0.08 -0.235 -0.24'
[../]
[disp_z_p16_pos]
type = PointValue
variable = disp_z
point = '-0.08 -0.235 0.24'
[../]
[disp_z_p16_neg]
type = PointValue
variable = disp_z
point = '-0.08 -0.235 -0.24'
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
line_search = none
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
start_time = 14391759
dt = 1000000
automatic_scaling = true
end_time = 38880000
l_max_its = 20
nl_max_its = 10
nl_rel_tol = 1e-6
nl_abs_tol = 1e-10
[]
[Outputs]
perf_graph = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
(test/tests/concrete_logcreep/concrete_logcreep.i)
# Tests the ConcreteLogarithmicCreepModel in its default configuration
#
# Expected value for strain_xx
# 0.001 * (1 + (1 - exp(-t)) + log(1 + t))
#
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
elem_type = HEX8
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]
[./stress_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./creep_strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
[./TensorMechanics]
displacements = 'disp_x disp_y disp_z'
use_displaced_mesh = true
[../]
[]
[AuxKernels]
[./stress_xx]
type = RankTwoAux
variable = stress_xx
rank_two_tensor = stress
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./strain_xx]
type = RankTwoAux
variable = strain_xx
rank_two_tensor = total_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./creep_strain_xx]
type = RankTwoAux
variable = creep_strain_xx
rank_two_tensor = creep_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[]
[BCs]
[./symmy]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0
[../]
[./symmx]
type = DirichletBC
variable = disp_x
boundary = left
value = 0
[../]
[./symmz]
type = DirichletBC
variable = disp_z
boundary = back
value = 0
[../]
[./axial_load]
type = NeumannBC
variable = disp_x
boundary = right
value = 10e6
[../]
[]
[Materials]
[./stress]
type = ComputeLinearViscoelasticStress
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
poissons_ratio = 0.2
youngs_modulus = 10e9
recoverable_youngs_modulus = 10e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
[../]
[./strain]
type = ComputeSmallStrain
displacements = 'disp_x disp_y disp_z'
[../]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
viscoelastic_model = logcreep
[../]
[]
[Postprocessors]
[./stress_xx]
type = ElementAverageValue
variable = stress_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./creep_strain_xx]
type = ElementAverageValue
variable = creep_strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
l_max_its = 50
l_tol = 1e-8
nl_max_its = 20
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
dtmin = 0.01
end_time = 100
[./TimeStepper]
type = LogConstantDT
log_dt = 0.1
first_dt = 0.1
[../]
[]
[Outputs]
file_base = concrete_logcreep_out
exodus = true
csv = true
[]
(assessment/asr_validation/wald2017b/analysis/A1-triaxial.i)
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
penalty = 1e12
[]
[Problem]
type = ReferenceResidualProblem
reference_vector = 'ref'
extra_tag_vectors = 'ref'
group_variables = 'disp_x disp_y disp_z'
[]
[Mesh]
file = A1-triaxial.e
[]
[Variables]
[./T]
order = FIRST
family = LAGRANGE
initial_condition = 10.6
[../]
[./rh]
order = FIRST
family = LAGRANGE
initial_condition = 0.8
block = 1
[../]
[]
[AuxVariables]
[./resid_x]
[../]
[./resid_y]
[../]
[./resid_z]
[../]
[./ASR_ex]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_vstrain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./volumetric_strain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_conductivity]
order = CONSTANT
family = Monomial
[../]
[./thermal_capacity]
order = CONSTANT
family = Monomial
[../]
[./moisture_capacity]
order = CONSTANT
family = Monomial
[../]
[./humidity_diffusivity]
order = CONSTANT
family = Monomial
[../]
[./water_content]
order = CONSTANT
family = Monomial
[../]
[./water_hydrated]
order = CONSTANT
family = Monomial
[../]
[damage_index]
order = CONSTANT
family = MONOMIAL
[]
[./area]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_strain]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Modules/TensorMechanics/Master]
[./concrete]
block = 1
strain = FINITE
add_variables = true
eigenstrain_names = 'asr_expansion thermal_expansion'
generate_output = 'stress_xx stress_yy stress_zz stress_xy stress_yz stress_zx vonmises_stress hydrostatic_stress elastic_strain_xx elastic_strain_yy elastic_strain_zz strain_xx strain_yy strain_zz'
extra_vector_tags = 'ref'
[../]
[]
[Modules/TensorMechanics/LineElementMaster]
[./Reinforcement_block]
block = '2 3 4'
truss = true
area = area
displacements = 'disp_x disp_y disp_z'
#Note: Intentially not including this here to have it give a nonzero
# displacement reference residual since it's an unrestrained problem
#extra_vector_tags = 'ref'
[../]
[]
[Constraints]
[./rebar_x2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[./rebar_x3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T3]
type = EqualValueEmbeddedConstraint
secondary = 3
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[./rebar_x4]
type = EqualValueEmbeddedConstraint
secondary = 4
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y4]
type = EqualValueEmbeddedConstraint
secondary = 4
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z4]
type = EqualValueEmbeddedConstraint
secondary = 4
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T4]
type = EqualValueEmbeddedConstraint
secondary = 4
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[]
[Kernels]
[./T_td]
type = ConcreteThermalTimeIntegration
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_diff]
type = ConcreteThermalConduction
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_conv]
type = ConcreteThermalConvection
variable = T
relative_humidity = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./T_adsorption]
type = ConcreteLatentHeat
variable = T
H = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_td]
type = ConcreteMoistureTimeIntegration
variable = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_diff]
type = ConcreteMoistureDiffusion
variable = rh
temperature = T
block = 1
extra_vector_tags = 'ref'
[../]
[./heat_dt]
type = TimeDerivative
variable = T
block = '2 3 4'
extra_vector_tags = 'ref'
[../]
[./heat_conduction]
type = HeatConduction
variable = T
diffusion_coefficient = 53.0
block = '2 3 4'
extra_vector_tags = 'ref'
[../]
[]
[AuxKernels]
[./ASR_ex]
type = MaterialRealAux
variable = ASR_ex
block = 1
property = ASR_extent
execute_on = 'timestep_end'
[../]
[./ASR_vstrain]
type = MaterialRealAux
block = 1
variable = ASR_vstrain
property = ASR_volumetric_strain
execute_on = 'timestep_end'
[../]
[./ASR_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./ASR_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_xy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xy
index_i = 0
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_yz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yz
index_i = 1
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_zx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zx
index_i = 0
index_j = 2
execute_on = 'timestep_end'
[../]
[./thermal_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./thermal_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./thermal_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./volumetric_strain]
type = RankTwoScalarAux
scalar_type = VolumetricStrain
rank_two_tensor = total_strain
variable = volumetric_strain
block = 1
[../]
[./k]
type = MaterialRealAux
variable = thermal_conductivity
property = thermal_conductivity
execute_on = 'timestep_end'
block = 1
[../]
[./capacity]
type = MaterialRealAux
variable = thermal_capacity
property = thermal_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_capacity]
type = MaterialRealAux
variable = moisture_capacity
property = moisture_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_duff]
type = MaterialRealAux
variable = humidity_diffusivity
property = humidity_diffusivity
execute_on = 'timestep_end'
block = 1
[../]
[./wc_duff]
type = MaterialRealAux
variable = water_content
property = moisture_content
execute_on = 'timestep_end'
block = 1
[../]
[./hydrw_duff]
type = MaterialRealAux
variable = water_hydrated
property = hydrated_water
execute_on = 'timestep_end'
block = 1
[../]
[damage_index]
type = MaterialRealAux
block = 1
variable = damage_index
property = damage_index
execute_on = timestep_end
[]
[./areax]
type = ConstantAux
block = '2'
variable = area
value = 1.33e-4
execute_on = 'initial timestep_begin'
[../]
[./areaz]
type = ConstantAux
block = '3'
variable = area
value = 1.33e-4
execute_on = 'initial timestep_begin'
[../]
[./areay]
type = ConstantAux
block = '4'
variable = area
value = 1.33e-4
execute_on = 'initial timestep_begin'
[../]
[./axial_stress]
type = MaterialRealAux
block = '2 3 4'
variable = axial_stress
property = axial_stress
[../]
[]
[Functions]
[./ramp_temp]
type = PiecewiseLinear
data_file = temperature_history.csv
format = columns
[../]
[./ramp_humidity]
type = PiecewiseLinear
data_file = humidity_history.csv
format = columns
[../]
[]
[Materials]
[./concrete]
type = ConcreteThermalMoisture
block = 1
# setup thermal property models and parameters
# options available: CONSTANT ASCE-1992 KODUR-2004 EUROCODE-2004 KIM-2003
thermal_conductivity_model = KODUR-2004
thermal_capacity_model = KODUR-2004
aggregate_type = Siliceous #options: Siliceous Carbonate
ref_density_of_concrete = 2231.0 # in kg/m^3
ref_specific_heat_of_concrete = 1100.0 # in J/(Kg.0C)
ref_thermal_conductivity_of_concrete = 3 # in W/(m.0C)
# setup moisture capacity and humidity diffusivity models
aggregate_pore_type = dense #options: dense porous
aggregate_mass = 1877.0 #mass of aggregate (kg) per m^3 of concrete
cement_type = 1 #options: 1 2 3 4
cement_mass = 354.0 #mass of cement (kg) per m^3 of concrete
water_to_cement_ratio = 0.5
concrete_cure_time = 28.0 #curing time in (days)
# options available for humidity diffusivity:
moisture_diffusivity_model = Bazant #options: Bazant Mensi
D1 = 3.0e-8
coupled_moisture_diffusivity_factor = 1.0e-2 # factor for mositure diffusivity due to heat
# coupled nonlinear variables
relative_humidity = rh
temperature = T
[../]
[./creep]
type = LinearViscoelasticStressUpdate
block = 1
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
block = 1
poissons_ratio = 0.22
youngs_modulus = 37.3e9
recoverable_youngs_modulus = 37.3e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
humidity = rh
temperature = T
activation_temperature = 23.0
[../]
[ASR_expansion]
type = ConcreteASREigenstrain
block = 1
expansion_type = Anisotropic
reference_temperature = 35.0
temperature_unit = Celsius
max_volumetric_expansion = 2.2e-2
characteristic_time = 18.9
latency_time = 18.0
characteristic_activation_energy = 5400.0
latency_activation_energy = 9400.0
stress_latency_factor = 1.0
compressive_strength = 31.0e6
compressive_stress_exponent = 0.0
expansion_stress_limit = 8.0e6
tensile_strength = 3.2e6
tensile_retention_factor = 1.0
tensile_absorption_factor = 1.0
ASR_dependent_tensile_strength = false
residual_tensile_strength_fraction = 1.0
temperature = T
relative_humidity = rh
rh_exponent = 1.0
eigenstrain_name = asr_expansion
absolute_tolerance = 1e-10
output_iteration_info_on_error = true
[]
[thermal_strain_concrete]
type = ComputeThermalExpansionEigenstrain
block = 1
temperature = T
thermal_expansion_coeff = 8.0e-6
stress_free_temperature = 10.6
eigenstrain_name = thermal_expansion
[]
[ASR_damage_concrete]
type = ConcreteASRMicrocrackingDamage
residual_youngs_modulus_fraction = 0.1
block = 1
[]
[./stress]
type = ComputeMultipleInelasticStress
block = 1
inelastic_models = 'creep'
damage_model = ASR_damage_concrete
[../]
[truss]
type = LinearElasticTruss
block = '2 3 4'
youngs_modulus = 2e11
temperature = T
thermal_expansion_coeff = 11.3e-6
temperature_ref = 10.6
[]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
block = 1
viscoelastic_model = logcreep
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = disp_x
boundary = '2000 2005'
value = 0.0
[../]
[./bottom]
type = DirichletBC
variable = disp_y
boundary = '2000 2001'
value = 0.0
[../]
[./back]
type = DirichletBC
variable = disp_z
boundary = '2000 2005'
value = 0.0
[../]
[./T]
type = FunctionDirichletBC
variable = T
boundary = '101 102 103 104 105 106'
function = ramp_temp
[../]
[./rh]
type = FunctionDirichletBC
variable = rh
boundary = '101 102 103 104 105 106'
function = ramp_humidity
[../]
[]
[Postprocessors]
[./nelem]
type = NumElems
[../]
[./ndof]
type = NumDOFs
[../]
[./ASR_strain]
type = ElementAverageValue
variable = ASR_vstrain
block = 1
[../]
[./ASR_strain_xx]
type = ElementAverageValue
variable = ASR_strain_xx
block = 1
[../]
[./ASR_strain_yy]
type = ElementAverageValue
variable = ASR_strain_yy
block = 1
[../]
[./ASR_strain_zz]
type = ElementAverageValue
variable = ASR_strain_zz
block = 1
[../]
[ASR_ext]
type = ElementAverageValue
variable = ASR_ex
block = 1
[]
[./vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 1
[../]
[./vstrain]
type = ElementAverageValue
variable = volumetric_strain
block = 1
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 1
[../]
[./strain_yy]
type = ElementAverageValue
variable = strain_yy
block = 1
[../]
[./strain_zz]
type = ElementAverageValue
variable = strain_zz
block = 1
[../]
[./temp]
type = ElementAverageValue
variable = T
block = 1
[../]
[./humidity]
type = ElementAverageValue
variable = rh
block = 1
[../]
[./thermal_strain_xx]
type = ElementAverageValue
variable = thermal_strain_xx
block = 1
[../]
[./thermal_strain_yy]
type = ElementAverageValue
variable = thermal_strain_yy
block = 1
[../]
[./thermal_strain_zz]
type = ElementAverageValue
variable = thermal_strain_zz
block = 1
[../]
[./disp_x_101]
type = SideAverageValue
variable = disp_x
boundary = 101
[../]
[./disp_x_102]
type = SideAverageValue
variable = disp_x
boundary = 102
[../]
[./disp_x_103]
type = SideAverageValue
variable = disp_x
boundary = 103
[../]
[./disp_x_104]
type = SideAverageValue
variable = disp_x
boundary = 104
[../]
[./disp_x_105]
type = SideAverageValue
variable = disp_x
boundary = 105
[../]
[./disp_x_106]
type = SideAverageValue
variable = disp_x
boundary = 106
[../]
[./disp_y_101]
type = SideAverageValue
variable = disp_y
boundary = 101
[../]
[./disp_y_102]
type = SideAverageValue
variable = disp_y
boundary = 102
[../]
[./disp_y_103]
type = SideAverageValue
variable = disp_y
boundary = 103
[../]
[./disp_y_104]
type = SideAverageValue
variable = disp_y
boundary = 104
[../]
[./disp_y_105]
type = SideAverageValue
variable = disp_y
boundary = 105
[../]
[./disp_y_106]
type = SideAverageValue
variable = disp_y
boundary = 106
[../]
[./disp_z_101]
type = SideAverageValue
variable = disp_z
boundary = 101
[../]
[./disp_z_102]
type = SideAverageValue
variable = disp_z
boundary = 102
[../]
[./disp_z_103]
type = SideAverageValue
variable = disp_z
boundary = 103
[../]
[./disp_z_104]
type = SideAverageValue
variable = disp_z
boundary = 104
[../]
[./disp_z_105]
type = SideAverageValue
variable = disp_z
boundary = 105
[../]
[./disp_z_106]
type = SideAverageValue
variable = disp_z
boundary = 106
[../]
[disp_x_p1_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 -0.08'
[../]
[disp_x_p1_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 -0.08'
[../]
[disp_x_p2_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 0.08'
[../]
[disp_x_p2_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 0.08'
[../]
[disp_x_p3_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.08'
[../]
[disp_x_p3_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.08'
[../]
[disp_x_p4_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.08'
[../]
[disp_x_p4_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.08'
[../]
[disp_x_p5_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.235'
[../]
[disp_x_p5_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.235'
[../]
[disp_x_p6_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.235'
[../]
[disp_x_p6_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.235'
[../]
[disp_y_p1_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 -0.08'
[../]
[disp_y_p1_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 -0.08'
[../]
[disp_y_p2_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 0.08'
[../]
[disp_y_p2_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 0.08'
[../]
[disp_y_p3_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.08'
[../]
[disp_y_p3_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.08'
[../]
[disp_y_p4_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.08'
[../]
[disp_y_p4_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.08'
[../]
[disp_y_p5_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.235'
[../]
[disp_y_p5_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.235'
[../]
[disp_y_p6_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.235'
[../]
[disp_y_p6_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.235'
[../]
[disp_y_p7_pos]
type = PointValue
variable = disp_y
point = '-0.235 0.24 0.08'
[../]
[disp_y_p7_neg]
type = PointValue
variable = disp_y
point = '-0.235 -0.24 0.08'
[../]
[disp_y_p8_pos]
type = PointValue
variable = disp_y
point = '0.235 0.24 0.08'
[../]
[disp_y_p8_neg]
type = PointValue
variable = disp_y
point = '0.235 -0.24 0.08'
[../]
[disp_z_p1_pos]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 0.24'
[../]
[disp_z_p1_neg]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 -0.24'
[../]
[disp_z_p2_pos]
type = PointValue
variable = disp_z
point = '-0.08 0.08 0.24'
[../]
[disp_z_p2_neg]
type = PointValue
variable = disp_z
point = '-0.08 0.08 -0.24'
[../]
[disp_z_p3_pos]
type = PointValue
variable = disp_z
point = '0.08 -0.08 0.24'
[../]
[disp_z_p3_neg]
type = PointValue
variable = disp_z
point = '0.08 -0.08 -0.24'
[../]
[disp_z_p4_pos]
type = PointValue
variable = disp_z
point = '0.08 0.08 0.24'
[../]
[disp_z_p4_neg]
type = PointValue
variable = disp_z
point = '0.08 0.08 -0.24'
[../]
[disp_z_p5_pos]
type = PointValue
variable = disp_z
point = '0.235 0.08 0.24'
[../]
[disp_z_p5_neg]
type = PointValue
variable = disp_z
point = '0.235 0.08 -0.24'
[../]
[disp_z_p6_pos]
type = PointValue
variable = disp_z
point = '-0.235 0.08 0.24'
[../]
[disp_z_p6_neg]
type = PointValue
variable = disp_z
point = '-0.235 0.08 -0.24'
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
line_search = none
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
start_time = 2419200
dt = 1000000
automatic_scaling = true
end_time = 38880000
l_max_its = 20
nl_max_its = 10
nl_rel_tol = 1e-6
nl_abs_tol = 1e-10
[]
[Outputs]
perf_graph = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
(test/tests/concrete_logcreep/concrete_logcreep_temp.i)
# Tests the ConcreteLogarithmicCreepModel at high temperature
#
# Expected value for strain_xx
# 0.001 * (1 + (1 - exp(-t * A)) + A * log(1 + t))
#
# with A = 1.506173
#
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
elem_type = HEX8
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
[../]
[]
[Functions]
[./temp_function]
type = ParsedFunction
value = 60
[../]
[]
[AuxVariables]
[./T]
order = CONSTANT
family = MONOMIAL
initial_condition = 60
[../]
[./stress_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./creep_strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
[./TensorMechanics]
displacements = 'disp_x disp_y disp_z'
use_displaced_mesh = true
[../]
[]
[AuxKernels]
[./temperature]
type = FunctionAux
variable = T
function = temp_function
[../]
[./stress_xx]
type = RankTwoAux
variable = stress_xx
rank_two_tensor = stress
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./strain_xx]
type = RankTwoAux
variable = strain_xx
rank_two_tensor = total_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./creep_strain_xx]
type = RankTwoAux
variable = creep_strain_xx
rank_two_tensor = creep_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[]
[BCs]
[./symmy]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0
[../]
[./symmx]
type = DirichletBC
variable = disp_x
boundary = left
value = 0
[../]
[./symmz]
type = DirichletBC
variable = disp_z
boundary = back
value = 0
[../]
[./axial_load]
type = NeumannBC
variable = disp_x
boundary = right
value = 10e6
[../]
[]
[Materials]
[./stress]
type = ComputeLinearViscoelasticStress
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
poissons_ratio = 0.2
youngs_modulus = 10e9
recoverable_youngs_modulus = 10e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
temperature = T
activation_temperature = 1000
[../]
[./strain]
type = ComputeSmallStrain
displacements = 'disp_x disp_y disp_z'
[../]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
viscoelastic_model = logcreep
[../]
[]
[Postprocessors]
[./stress_xx]
type = ElementAverageValue
variable = stress_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./creep_strain_xx]
type = ElementAverageValue
variable = creep_strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
l_max_its = 50
l_tol = 1e-8
nl_max_its = 20
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
dtmin = 0.01
end_time = 100
[./TimeStepper]
type = LogConstantDT
log_dt = 0.1
first_dt = 0.1
[../]
[]
[Outputs]
file_base = concrete_logcreep_temp_out
exodus = true
csv = true
[]
(test/tests/concrete_logcreep/concrete_logcreep_humidity.i)
# Tests the ConcreteLogarithmicCreepModel under low relative humidity
#
# Expected value for strain_xx
# 0.001 * (1 + 0.5 * (1 - exp(-t)) + 0.5 * log(1 + t))
#
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
elem_type = HEX8
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
[../]
[]
[Functions]
[./humidity_function]
type = ParsedFunction
value = 0.5
[../]
[]
[AuxVariables]
[./h]
order = CONSTANT
family = MONOMIAL
[../]
[./stress_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./creep_strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
[./TensorMechanics]
displacements = 'disp_x disp_y disp_z'
use_displaced_mesh = true
[../]
[]
[AuxKernels]
[./humidity]
type = FunctionAux
variable = h
function = humidity_function
execute_on = timestep_begin
[../]
[./stress_xx]
type = RankTwoAux
variable = stress_xx
rank_two_tensor = stress
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./strain_xx]
type = RankTwoAux
variable = strain_xx
rank_two_tensor = total_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./creep_strain_xx]
type = RankTwoAux
variable = creep_strain_xx
rank_two_tensor = creep_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[]
[BCs]
[./symmy]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0
[../]
[./symmx]
type = DirichletBC
variable = disp_x
boundary = left
value = 0
[../]
[./symmz]
type = DirichletBC
variable = disp_z
boundary = back
value = 0
[../]
[./axial_load]
type = NeumannBC
variable = disp_x
boundary = right
value = 10e6
[../]
[]
[Materials]
[./stress]
type = ComputeLinearViscoelasticStress
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
poissons_ratio = 0.2
youngs_modulus = 10e9
recoverable_youngs_modulus = 10e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
humidity = h
[../]
[./strain]
type = ComputeSmallStrain
displacements = 'disp_x disp_y disp_z'
[../]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
viscoelastic_model = logcreep
[../]
[]
[Postprocessors]
[./stress_xx]
type = ElementAverageValue
variable = stress_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./creep_strain_xx]
type = ElementAverageValue
variable = creep_strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
l_max_its = 50
l_tol = 1e-8
nl_max_its = 20
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
dtmin = 0.01
end_time = 100
[./TimeStepper]
type = LogConstantDT
log_dt = 0.1
first_dt = 0.1
[../]
[]
[Outputs]
file_base = concrete_logcreep_humidity_out
exodus = true
csv = true
[]
(test/tests/concrete_logcreep/concrete_logcreep_drying.i)
# Tests the ConcreteLogarithmicCreepModel in drying conditions
#
# Expected value for strain_xx
# 0.001 * (1 + 1.02 * (1 - exp(-t)) + 0.01 * (101 * log(1 + t)) + 0.025 * t)
#
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
elem_type = HEX8
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
[../]
[]
[Functions]
[./humidity_function]
type = PiecewiseLinear
x = '0 100 1000'
y = '1 0 0'
[../]
[]
[AuxVariables]
[./h]
order = CONSTANT
family = MONOMIAL
initial_condition = 1
[../]
[./stress_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[./creep_strain_xx]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
[./TensorMechanics]
displacements = 'disp_x disp_y disp_z'
use_displaced_mesh = true
[../]
[]
[AuxKernels]
[./humidity]
type = FunctionAux
variable = h
function = humidity_function
execute_on = timestep_begin
[../]
[./stress_xx]
type = RankTwoAux
variable = stress_xx
rank_two_tensor = stress
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./strain_xx]
type = RankTwoAux
variable = strain_xx
rank_two_tensor = total_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[./creep_strain_xx]
type = RankTwoAux
variable = creep_strain_xx
rank_two_tensor = creep_strain
index_j = 0
index_i = 0
execute_on = timestep_end
[../]
[]
[BCs]
[./symmy]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0
[../]
[./symmx]
type = DirichletBC
variable = disp_x
boundary = left
value = 0
[../]
[./symmz]
type = DirichletBC
variable = disp_z
boundary = back
value = 0
[../]
[./axial_load]
type = NeumannBC
variable = disp_x
boundary = right
value = 10e6
[../]
[]
[Materials]
[./stress]
type = ComputeLinearViscoelasticStress
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
poissons_ratio = 0.2
youngs_modulus = 10e9
recoverable_youngs_modulus = 10e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
humidity = h
drying_creep_viscosity = 0.25
[../]
[./strain]
type = ComputeSmallStrain
displacements = 'disp_x disp_y disp_z'
[../]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
viscoelastic_model = logcreep
[../]
[]
[Postprocessors]
[./stress_xx]
type = ElementAverageValue
variable = stress_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./creep_strain_xx]
type = ElementAverageValue
variable = creep_strain_xx
block = 'ANY_BLOCK_ID 0'
[../]
[./humidity]
type = ElementAverageValue
variable = h
block = 'ANY_BLOCK_ID 0'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
l_max_its = 50
l_tol = 1e-8
nl_max_its = 20
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
dtmin = 0.01
end_time = 100
[./TimeStepper]
type = LogConstantDT
log_dt = 0.1
first_dt = 0.1
[../]
[]
[Outputs]
file_base = concrete_logcreep_drying_out
exodus = true
csv = true
[]
(assessment/asr_validation/wald2017b/analysis/A1-uniaxial.i)
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
penalty = 1e12
[]
[Problem]
type = ReferenceResidualProblem
reference_vector = 'ref'
extra_tag_vectors = 'ref'
group_variables = 'disp_x disp_y disp_z'
[]
[Mesh]
file = A1-uniaxial.e
[]
[Variables]
[./T]
order = FIRST
family = LAGRANGE
initial_condition = 10.6
[../]
[./rh]
order = FIRST
family = LAGRANGE
initial_condition = 0.8
block = 1
[../]
[]
[AuxVariables]
[./resid_x]
[../]
[./resid_y]
[../]
[./resid_z]
[../]
[./ASR_ex]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_vstrain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_xy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_yz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./ASR_strain_zx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./volumetric_strain]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_xx]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_yy]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_strain_zz]
order = CONSTANT
family = MONOMIAL
block = 1
[../]
[./thermal_conductivity]
order = CONSTANT
family = Monomial
[../]
[./thermal_capacity]
order = CONSTANT
family = Monomial
[../]
[./moisture_capacity]
order = CONSTANT
family = Monomial
[../]
[./humidity_diffusivity]
order = CONSTANT
family = Monomial
[../]
[./water_content]
order = CONSTANT
family = Monomial
[../]
[./water_hydrated]
order = CONSTANT
family = Monomial
[../]
[damage_index]
order = CONSTANT
family = MONOMIAL
[]
[./area]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_stress]
order = CONSTANT
family = MONOMIAL
[../]
[./axial_strain]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Modules/TensorMechanics/Master]
[./concrete]
block = 1
strain = FINITE
add_variables = true
eigenstrain_names = 'asr_expansion thermal_expansion'
generate_output = 'stress_xx stress_yy stress_zz stress_xy stress_yz stress_zx vonmises_stress hydrostatic_stress elastic_strain_xx elastic_strain_yy elastic_strain_zz strain_xx strain_yy strain_zz'
extra_vector_tags = 'ref'
[../]
[]
[Modules/TensorMechanics/LineElementMaster]
[./Reinforcement_block]
block = '2 '
truss = true
area = area
displacements = 'disp_x disp_y disp_z'
#Note: Intentially not including this here to have it give a nonzero
# displacement reference residual since it's an unrestrained problem
#extra_vector_tags = 'ref'
[../]
[]
[Constraints]
[./rebar_x2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_x'
primary_variable = 'disp_x'
formulation = penalty
[../]
[./rebar_y2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_y'
primary_variable = 'disp_y'
formulation = penalty
[../]
[./rebar_z2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'disp_z'
primary_variable = 'disp_z'
formulation = penalty
[../]
[./rebar_T2]
type = EqualValueEmbeddedConstraint
secondary = 2
primary = 1
variable = 'T'
primary_variable = 'T'
formulation = penalty
penalty = 1e6
[../]
[]
[Kernels]
[./T_td]
type = ConcreteThermalTimeIntegration
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_diff]
type = ConcreteThermalConduction
variable = T
block = 1
extra_vector_tags = 'ref'
[../]
[./T_conv]
type = ConcreteThermalConvection
variable = T
relative_humidity = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./T_adsorption]
type = ConcreteLatentHeat
variable = T
H = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_td]
type = ConcreteMoistureTimeIntegration
variable = rh
block = 1
extra_vector_tags = 'ref'
[../]
[./rh_diff]
type = ConcreteMoistureDiffusion
variable = rh
temperature = T
block = 1
extra_vector_tags = 'ref'
[../]
[./heat_dt]
type = TimeDerivative
variable = T
block = 2
extra_vector_tags = 'ref'
[../]
[./heat_conduction]
type = HeatConduction
variable = T
diffusion_coefficient = 53.0
block = 2
extra_vector_tags = 'ref'
[../]
[]
[AuxKernels]
[./ASR_ex]
type = MaterialRealAux
variable = ASR_ex
block = 1
property = ASR_extent
execute_on = 'timestep_end'
[../]
[./ASR_vstrain]
type = MaterialRealAux
block = 1
variable = ASR_vstrain
property = ASR_volumetric_strain
execute_on = 'timestep_end'
[../]
[./ASR_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./ASR_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_xy]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_xy
index_i = 0
index_j = 1
execute_on = 'timestep_end'
[../]
[./ASR_strain_yz]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_yz
index_i = 1
index_j = 2
execute_on = 'timestep_end'
[../]
[./ASR_strain_zx]
type = RankTwoAux
block = 1
rank_two_tensor = asr_expansion
variable = ASR_strain_zx
index_i = 0
index_j = 2
execute_on = 'timestep_end'
[../]
[./thermal_strain_xx]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_xx
index_i = 0
index_j = 0
execute_on = 'timestep_end'
[../]
[./thermal_strain_yy]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_yy
index_i = 1
index_j = 1
execute_on = 'timestep_end'
[../]
[./thermal_strain_zz]
type = RankTwoAux
block = 1
rank_two_tensor = thermal_expansion
variable = thermal_strain_zz
index_i = 2
index_j = 2
execute_on = 'timestep_end'
[../]
[./volumetric_strain]
type = RankTwoScalarAux
scalar_type = VolumetricStrain
rank_two_tensor = total_strain
variable = volumetric_strain
block = 1
[../]
[./k]
type = MaterialRealAux
variable = thermal_conductivity
property = thermal_conductivity
execute_on = 'timestep_end'
block = 1
[../]
[./capacity]
type = MaterialRealAux
variable = thermal_capacity
property = thermal_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_capacity]
type = MaterialRealAux
variable = moisture_capacity
property = moisture_capacity
execute_on = 'timestep_end'
block = 1
[../]
[./rh_duff]
type = MaterialRealAux
variable = humidity_diffusivity
property = humidity_diffusivity
execute_on = 'timestep_end'
block = 1
[../]
[./wc_duff]
type = MaterialRealAux
variable = water_content
property = moisture_content
execute_on = 'timestep_end'
block = 1
[../]
[./hydrw_duff]
type = MaterialRealAux
variable = water_hydrated
property = hydrated_water
execute_on = 'timestep_end'
block = 1
[../]
[damage_index]
type = MaterialRealAux
block = 1
variable = damage_index
property = damage_index
execute_on = timestep_end
[]
[./area]
type = ConstantAux
block = '2'
variable = area
value = 1.33e-4
execute_on = 'initial timestep_begin'
[../]
[./axial_stress]
type = MaterialRealAux
block = '2'
variable = axial_stress
property = axial_stress
[../]
[]
[Functions]
[./ramp_temp]
type = PiecewiseLinear
data_file = temperature_history.csv
format = columns
[../]
[./ramp_humidity]
type = PiecewiseLinear
data_file = humidity_history.csv
format = columns
[../]
[]
[Materials]
[./concrete]
type = ConcreteThermalMoisture
block = 1
# setup thermal property models and parameters
# options available: CONSTANT ASCE-1992 KODUR-2004 EUROCODE-2004 KIM-2003
thermal_conductivity_model = KODUR-2004
thermal_capacity_model = KODUR-2004
aggregate_type = Siliceous #options: Siliceous Carbonate
ref_density_of_concrete = 2231.0 # in kg/m^3
ref_specific_heat_of_concrete = 1100.0 # in J/(Kg.0C)
ref_thermal_conductivity_of_concrete = 3 # in W/(m.0C)
# setup moisture capacity and humidity diffusivity models
aggregate_pore_type = dense #options: dense porous
aggregate_mass = 1877.0 #mass of aggregate (kg) per m^3 of concrete
cement_type = 1 #options: 1 2 3 4
cement_mass = 354.0 #mass of cement (kg) per m^3 of concrete
water_to_cement_ratio = 0.5
concrete_cure_time = 28.0 #curing time in (days)
# options available for humidity diffusivity:
moisture_diffusivity_model = Bazant #options: Bazant Mensi
D1 = 3.0e-8
coupled_moisture_diffusivity_factor = 1.0e-2 # factor for mositure diffusivity due to heat
# coupled nonlinear variables
relative_humidity = rh
temperature = T
[../]
[./creep]
type = LinearViscoelasticStressUpdate
block = 1
[../]
[./logcreep]
type = ConcreteLogarithmicCreepModel
block = 1
poissons_ratio = 0.22
youngs_modulus = 37.3e9
recoverable_youngs_modulus = 37.3e9
recoverable_viscosity = 1
long_term_viscosity = 1
long_term_characteristic_time = 1
humidity = rh
temperature = T
activation_temperature = 23.0
[../]
[ASR_expansion]
type = ConcreteASREigenstrain
block = 1
expansion_type = Anisotropic
reference_temperature = 35.0
temperature_unit = Celsius
max_volumetric_expansion = 2.2e-2
characteristic_time = 18.9
latency_time = 18.0
characteristic_activation_energy = 5400.0
latency_activation_energy = 9400.0
stress_latency_factor = 1.0
compressive_strength = 31.0e6
compressive_stress_exponent = 0.0
expansion_stress_limit = 8.0e6
tensile_strength = 3.2e6
tensile_retention_factor = 1.0
tensile_absorption_factor = 1.0
ASR_dependent_tensile_strength = false
residual_tensile_strength_fraction = 1.0
temperature = T
relative_humidity = rh
rh_exponent = 1.0
eigenstrain_name = asr_expansion
absolute_tolerance = 1e-10
output_iteration_info_on_error = true
[]
[thermal_strain_concrete]
type = ComputeThermalExpansionEigenstrain
block = 1
temperature = T
thermal_expansion_coeff = 8.0e-6
stress_free_temperature = 10.6
eigenstrain_name = thermal_expansion
[]
[ASR_damage_concrete]
type = ConcreteASRMicrocrackingDamage
residual_youngs_modulus_fraction = 0.1
block = 1
[]
[./stress]
type = ComputeMultipleInelasticStress
block = 1
inelastic_models = 'creep'
damage_model = ASR_damage_concrete
[../]
[truss]
type = LinearElasticTruss
block = '2 '
youngs_modulus = 2e11
temperature = T
thermal_expansion_coeff = 11.3e-6
temperature_ref = 10.6
[]
[]
[UserObjects]
[./visco_update]
type = LinearViscoelasticityManager
block = 1
viscoelastic_model = logcreep
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = disp_x
boundary = '2000 2005'
value = 0.0
[../]
[./bottom]
type = DirichletBC
variable = disp_y
boundary = '2000 2001'
value = 0.0
[../]
[./back]
type = DirichletBC
variable = disp_z
boundary = '2000 2005'
value = 0.0
[../]
[./T]
type = FunctionDirichletBC
variable = T
boundary = '101 102 103 104 105 106'
function = ramp_temp
[../]
[./rh]
type = FunctionDirichletBC
variable = rh
boundary = '101 102 103 104 105 106'
function = ramp_humidity
[../]
[]
[Postprocessors]
[./nelem]
type = NumElems
[../]
[./ndof]
type = NumDOFs
[../]
[./ASR_strain]
type = ElementAverageValue
variable = ASR_vstrain
block = 1
[../]
[./ASR_strain_xx]
type = ElementAverageValue
variable = ASR_strain_xx
block = 1
[../]
[./ASR_strain_yy]
type = ElementAverageValue
variable = ASR_strain_yy
block = 1
[../]
[./ASR_strain_zz]
type = ElementAverageValue
variable = ASR_strain_zz
block = 1
[../]
[ASR_ext]
type = ElementAverageValue
variable = ASR_ex
block = 1
[]
[./vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 1
[../]
[./vstrain]
type = ElementAverageValue
variable = volumetric_strain
block = 1
[../]
[./strain_xx]
type = ElementAverageValue
variable = strain_xx
block = 1
[../]
[./strain_yy]
type = ElementAverageValue
variable = strain_yy
block = 1
[../]
[./strain_zz]
type = ElementAverageValue
variable = strain_zz
block = 1
[../]
[./temp]
type = ElementAverageValue
variable = T
block = 1
[../]
[./humidity]
type = ElementAverageValue
variable = rh
block = 1
[../]
[./thermal_strain_xx]
type = ElementAverageValue
variable = thermal_strain_xx
block = 1
[../]
[./thermal_strain_yy]
type = ElementAverageValue
variable = thermal_strain_yy
block = 1
[../]
[./thermal_strain_zz]
type = ElementAverageValue
variable = thermal_strain_zz
block = 1
[../]
[./disp_x_101]
type = SideAverageValue
variable = disp_x
boundary = 101
[../]
[./disp_x_102]
type = SideAverageValue
variable = disp_x
boundary = 102
[../]
[./disp_x_103]
type = SideAverageValue
variable = disp_x
boundary = 103
[../]
[./disp_x_104]
type = SideAverageValue
variable = disp_x
boundary = 104
[../]
[./disp_x_105]
type = SideAverageValue
variable = disp_x
boundary = 105
[../]
[./disp_x_106]
type = SideAverageValue
variable = disp_x
boundary = 106
[../]
[./disp_y_101]
type = SideAverageValue
variable = disp_y
boundary = 101
[../]
[./disp_y_102]
type = SideAverageValue
variable = disp_y
boundary = 102
[../]
[./disp_y_103]
type = SideAverageValue
variable = disp_y
boundary = 103
[../]
[./disp_y_104]
type = SideAverageValue
variable = disp_y
boundary = 104
[../]
[./disp_y_105]
type = SideAverageValue
variable = disp_y
boundary = 105
[../]
[./disp_y_106]
type = SideAverageValue
variable = disp_y
boundary = 106
[../]
[./disp_z_101]
type = SideAverageValue
variable = disp_z
boundary = 101
[../]
[./disp_z_102]
type = SideAverageValue
variable = disp_z
boundary = 102
[../]
[./disp_z_103]
type = SideAverageValue
variable = disp_z
boundary = 103
[../]
[./disp_z_104]
type = SideAverageValue
variable = disp_z
boundary = 104
[../]
[./disp_z_105]
type = SideAverageValue
variable = disp_z
boundary = 105
[../]
[./disp_z_106]
type = SideAverageValue
variable = disp_z
boundary = 106
[../]
[disp_x_p1_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 -0.08'
[../]
[disp_x_p1_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 -0.08'
[../]
[disp_x_p2_pos]
type = PointValue
variable = disp_x
point = '0.24 -0.08 0.08'
[../]
[disp_x_p2_neg]
type = PointValue
variable = disp_x
point = '-0.24 -0.08 0.08'
[../]
[disp_x_p3_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.08'
[../]
[disp_x_p3_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.08'
[../]
[disp_x_p4_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.08'
[../]
[disp_x_p4_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.08'
[../]
[disp_x_p5_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 -0.235'
[../]
[disp_x_p5_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 -0.235'
[../]
[disp_x_p6_pos]
type = PointValue
variable = disp_x
point = '0.24 0.08 0.235'
[../]
[disp_x_p6_neg]
type = PointValue
variable = disp_x
point = '-0.24 0.08 0.235'
[../]
[disp_y_p1_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 -0.08'
[../]
[disp_y_p1_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 -0.08'
[../]
[disp_y_p2_pos]
type = PointValue
variable = disp_y
point = '-0.08 0.24 0.08'
[../]
[disp_y_p2_neg]
type = PointValue
variable = disp_y
point = '-0.08 -0.24 0.08'
[../]
[disp_y_p3_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.08'
[../]
[disp_y_p3_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.08'
[../]
[disp_y_p4_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.08'
[../]
[disp_y_p4_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.08'
[../]
[disp_y_p5_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 -0.235'
[../]
[disp_y_p5_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 -0.235'
[../]
[disp_y_p6_pos]
type = PointValue
variable = disp_y
point = '0.08 0.24 0.235'
[../]
[disp_y_p6_neg]
type = PointValue
variable = disp_y
point = '0.08 -0.24 0.235'
[../]
[disp_y_p7_pos]
type = PointValue
variable = disp_y
point = '-0.235 0.24 0.08'
[../]
[disp_y_p7_neg]
type = PointValue
variable = disp_y
point = '-0.235 -0.24 0.08'
[../]
[disp_y_p8_pos]
type = PointValue
variable = disp_y
point = '0.235 0.24 0.08'
[../]
[disp_y_p8_neg]
type = PointValue
variable = disp_y
point = '0.235 -0.24 0.08'
[../]
[disp_z_p1_pos]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 0.24'
[../]
[disp_z_p1_neg]
type = PointValue
variable = disp_z
point = '-0.08 -0.08 -0.24'
[../]
[disp_z_p2_pos]
type = PointValue
variable = disp_z
point = '-0.08 0.08 0.24'
[../]
[disp_z_p2_neg]
type = PointValue
variable = disp_z
point = '-0.08 0.08 -0.24'
[../]
[disp_z_p3_pos]
type = PointValue
variable = disp_z
point = '0.08 -0.08 0.24'
[../]
[disp_z_p3_neg]
type = PointValue
variable = disp_z
point = '0.08 -0.08 -0.24'
[../]
[disp_z_p4_pos]
type = PointValue
variable = disp_z
point = '0.08 0.08 0.24'
[../]
[disp_z_p4_neg]
type = PointValue
variable = disp_z
point = '0.08 0.08 -0.24'
[../]
[disp_z_p5_pos]
type = PointValue
variable = disp_z
point = '0.235 0.08 0.24'
[../]
[disp_z_p5_neg]
type = PointValue
variable = disp_z
point = '0.235 0.08 -0.24'
[../]
[disp_z_p6_pos]
type = PointValue
variable = disp_z
point = '-0.235 0.08 0.24'
[../]
[disp_z_p6_neg]
type = PointValue
variable = disp_z
point = '-0.235 0.08 -0.24'
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
line_search = none
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
start_time = 2419200
dt = 1000000
automatic_scaling = true
end_time = 38880000
l_max_its = 20
nl_max_its = 10
nl_rel_tol = 1e-6
nl_abs_tol = 1e-10
[]
[Outputs]
perf_graph = true
csv = true
#exodus = true #Turned off to save space
[]
[Debug]
show_var_residual_norms = true
[]