- markersA list of marker names to combine into a single marker.
C++ Type:std::vector<MarkerName>
Controllable:No
Description:A list of marker names to combine into a single marker.
ComboMarker
A marker that converts many markers into a single marker by considering the maximum value of the listed markers (i.e., refinement takes precedent).
Description
The ComboMarker
is used to combine multiple markers into a single marker. This is done by taking the maximum value of the marker value from the supplied markers. Therefore, refinement of an element takes precedence.
Example Input Syntax
[Adaptivity]
[./Markers]
[./box]
type = BoxMarker
bottom_left = '0.3 0.3 0'
top_right = '0.6 0.6 0'
inside = refine
outside = do_nothing
[../]
[./combo]
type = ComboMarker
markers = 'box box2'
[../]
[./box2]
type = BoxMarker
bottom_left = '0.5 0.5 0'
top_right = '0.8 0.8 0'
inside = refine
outside = coarsen
[../]
[../]
[]
Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:No
Description:Set the enabled status of the MooseObject.
- outputsVector of output names where you would like to restrict the output of variables(s) associated with this object
C++ Type:std::vector<OutputName>
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
Advanced Parameters
Input Files
- (test/tests/dirackernels/point_caching/point_caching_adaptive_refinement.i)
- (test/tests/markers/dont_mark/dont_mark_test.i)
- (modules/phase_field/test/tests/GBType/GB_Type_Phase1.i)
- (test/tests/executioners/adapt_and_modify/adapt_and_modify_heavy.i)
- (modules/navier_stokes/examples/laser-welding/2d.i)
- (test/tests/markers/combo_marker/combo_marker_block_restricted.i)
- (modules/navier_stokes/examples/laser-welding/3d.i)
- (test/tests/executioners/adapt_and_modify/adapt_and_modify.i)
- (modules/phase_field/examples/nucleation/refine.i)
- (test/tests/markers/combo_marker/combo_marker_test.i)
(test/tests/markers/combo_marker/combo_marker_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
nz = 0
zmax = 0
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = 1
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = 2
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Adaptivity]
[./Markers]
[./box]
type = BoxMarker
bottom_left = '0.3 0.3 0'
top_right = '0.6 0.6 0'
inside = refine
outside = do_nothing
[../]
[./combo]
type = ComboMarker
markers = 'box box2'
[../]
[./box2]
type = BoxMarker
bottom_left = '0.5 0.5 0'
top_right = '0.8 0.8 0'
inside = refine
outside = coarsen
[../]
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/dirackernels/point_caching/point_caching_adaptive_refinement.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = QUAD4
uniform_refine = 2
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
active = 'diff'
[./diff]
type = Diffusion
variable = u
[../]
[]
[DiracKernels]
active = 'point_source'
[./point_source]
type = CachingPointSource
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = 3
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = 1
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Adaptivity]
steps = 3
marker = 'combo'
[./Markers]
[./combo]
# In a real problem you would want to mark based on an error
# indicator, but we want this test to run consistently in
# parallel, so we just mark elements within a box for
# refinement. The boxes here are based on the 8x8
# uniformly-refined initial grid.
type = ComboMarker
markers = 'box1 box2 box3'
[../]
[./box1]
type = BoxMarker
bottom_left = '0.125 0.625 0'
top_right = '0.375 0.875 0'
inside = refine
outside = dont_mark
[../]
[./box2]
type = BoxMarker
bottom_left = '0.625 0.625 0'
top_right = '0.875 0.875 0'
inside = refine
outside = dont_mark
[../]
[./box3]
type = BoxMarker
bottom_left = '0.625 0.125 0'
top_right = '0.875 0.375 0'
inside = refine
outside = dont_mark
[../]
[../]
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(test/tests/markers/dont_mark/dont_mark_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
nz = 0
zmax = 0
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = 1
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = 2
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Adaptivity]
[./Markers]
[./box]
type = BoxMarker
bottom_left = '0.3 0.3 0'
top_right = '0.6 0.6 0'
inside = refine
outside = coarsen
[../]
[./combo]
type = ComboMarker
markers = 'box box2'
[../]
[./box2]
type = BoxMarker
bottom_left = '0.5 0.5 0'
top_right = '0.8 0.8 0'
inside = dont_mark
outside = refine
[../]
[../]
[]
[Outputs]
exodus = true
[]
(modules/phase_field/test/tests/GBType/GB_Type_Phase1.i)
# MOOSE input file
# Written by Pierre-Clement Simon - Idaho National Laboratory
#
# Project:
# TRISO fuel fission gas transport: Silver diffusion in silicon carbide
#
# Published with:
# ---
#
# Phase Field Model: Isotropic diffusion equation
# type: Transient
# Grain structure: Single grain
# BCs: Fixed value on the right, flux on the left
#
#
# Info:
# - Input file used to generate polycrystals for SiC
#
# Updates from previous file:
# -
#
# Units
# length: --
# time: --
# energy: --
# quantity: --
# This simulation predicts GB migration of a 2D copper polycrystal with 15 grains
# Mesh adaptivity (new system) and time step adaptivity are used
# An AuxVariable is used to calculate the grain boundary locations
# Postprocessors are used to record time step and the number of grains
# We are not using the GrainTracker in this example so the number
# of order paramaters must match the number of grains.
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
# Two Parallel Grains
filename = 'EBSD_ThreeGrains.txt'
[]
[]
[GlobalParams]
# Parameters used by several kernels that are defined globally to simplify input file
op_num = 6 # Number of grains
var_name_base = gr # Base name of grains
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
enable_var_coloring = true
# output_adjacency_matrix = true
[]
[grain_tracker]
type = GrainTracker
threshold = 0.001
connecting_threshold = 0.008
compute_var_to_feature_map = true
compute_halo_maps = true # For displaying HALO fields
remap_grains = true
polycrystal_ic_uo = ebsd
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[]
[Variables]
# Variable block, where all variables in the simulation are declared
[./PolycrystalVariables]
# Custom action that created all of the grain variables and sets their initial condition
[../]
[]
[AuxVariables]
# Dependent variables
[./bnds]
# Variable used to visualize the grain boundaries in the simulation
[../]
[./unique_grains]
order = CONSTANT
family = MONOMIAL
[../]
[./aphi1]
order = CONSTANT
family = MONOMIAL
[../]
[./bPhi]
order = CONSTANT
family = MONOMIAL
[../]
[./cphi2]
order = CONSTANT
family = MONOMIAL
[../]
[./ebsd_numbers]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
# Kernel block, where the kernels defining the residual equations are set up.
[./PolycrystalKernel]
# Custom action creating all necessary kernels for grain growth. All input parameters are up in GlobalParams
[../]
[]
[AuxKernels]
# AuxKernel block, defining the equations used to calculate the auxvars
[./bnds_aux]
# AuxKernel that calculates the GB term
type = BndsCalcAux
variable = bnds
execute_on = 'initial timestep_end'
[../]
# generate the unique ID from grain_tracker
[./unique_grains]
type = FeatureFloodCountAux
variable = unique_grains
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
field_display = UNIQUE_REGION
[../]
# The phi will output the Euler angle from EBSD data, and the data structure
# will change with the guide from grain_tracker
[./aphi1]
type = OutputEulerAngles
variable = aphi1
euler_angle_provider = ebsd_reader
grain_tracker = grain_tracker
output_euler_angle = 'phi1'
execute_on = 'INITIAL TIMESTEP_END'
[../]
[./bPhi]
type = OutputEulerAngles
variable = bPhi
euler_angle_provider = ebsd_reader
grain_tracker = grain_tracker
output_euler_angle = 'Phi'
execute_on = 'INITIAL TIMESTEP_END'
[../]
[./cphi2]
type = OutputEulerAngles
variable = cphi2
euler_angle_provider = ebsd_reader
grain_tracker = grain_tracker
output_euler_angle = 'phi2'
execute_on = 'INITIAL TIMESTEP_END'
[../]
# Import the unique grain ID from ebsd data, and the data structure
# will change with the guide from grain_tracker
[ebsd_numbers]
type = EBSDReaderAvgDataAux
data_name = feature_id
ebsd_reader = ebsd_reader
grain_tracker = grain_tracker
variable = ebsd_numbers
execute_on = 'initial timestep_end'
[../]
[]
[BCs]
# Boundary Condition block
[./Periodic]
[./top_bottom]
auto_direction = 'x y' # Makes problem periodic in the x and y directions
[../]
[../]
[]
[Materials]
[./CuGrGr]
# Material properties
type = GBEvolution # Quantitative material properties for copper grain growth. Dimensions are nm and ns
GBmob0 = 2.5e-6 # Mobility prefactor for Cu from schonfelder1997molecular bibtex entry
GBenergy = 0.708 # GB energy for Cu from schonfelder1997molecular bibtex entry
Q = 0.23 # Activation energy for grain growth from Schonfelder 1997
T = 450 # Constant temperature of the simulation (for mobility calculation)
wGB = 6 # Width of the diffuse GB
[../]
[./GB_type]
# The new developed Miso Bnds Aux Kernel
type = ComputeGBMisorientationType
ebsd_reader = ebsd_reader
grain_tracker = grain_tracker
output_properties = 'gb_type'
outputs = exodus
[../]
[]
[Postprocessors]
# Scalar postprocessors
[./dt]
# Outputs the current time step
type = TimestepSize
[../]
[n_elements]
type = NumElements
execute_on = 'initial timestep_end'
[]
[n_nodes]
type = NumNodes
execute_on = 'initial timestep_end'
[]
[DOFs]
type = NumDOFs
[]
[]
[Adaptivity]
initial_steps = 1
max_h_level = 1
marker = combined
[./Indicators]
[./error]
type = GradientJumpIndicator
variable = bnds
[../]
[../]
[./Markers]
[./bound_adapt]
type = ValueThresholdMarker
third_state = DO_NOTHING
coarsen = 0.999 #1.0
refine = 0.95 #0.95
variable = bnds
invert = true
[../]
[./errorfrac]
type = ErrorFractionMarker
coarsen = 0.1
indicator = error
refine = 0.7
[../]
[./combined]
type = ComboMarker
markers = 'bound_adapt errorfrac'
[../]
[../]
[]
[Executioner]
type = Transient # Type of executioner, here it is transient with an adaptive time step
scheme = bdf2 # Type of time integration (2nd order backward euler), defaults to 1st order backward euler
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type -pc_hypre_boomeramg_strong_threshold'
petsc_options_value = ' hypre boomeramg 0.7'
l_max_its = 30 # Max number of linear iterations
l_tol = 1e-4 # Relative tolerance for linear solves
nl_max_its = 40 # Max number of nonlinear iterations
nl_abs_tol = 1e-11 # Relative tolerance for nonlienar solves
nl_rel_tol = 1e-10 # Absolute tolerance for nonlienar solves
[TimeStepper]
type = IterationAdaptiveDT
cutback_factor = 0.9
dt = 1
growth_factor = 1.1
optimal_iterations = 7
[]
start_time = 0.0
num_steps = 2
[]
[Outputs]
perf_graph = true
exodus = true
[./console]
type = Console
max_rows = 10
[../]
[]
(test/tests/executioners/adapt_and_modify/adapt_and_modify_heavy.i)
[Mesh]
# This example uses Adaptivity Indicators, which are written out as
# CONSTANT MONOMIAL variables, which don't currently work correctly
# 2122 for more information.
type = GeneratedMesh
dim = 2
nx = 15
ny = 15
# parallel_type = replicated
[]
[Variables]
[./u]
[../]
[]
[AuxVariables]
[./elem]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[./td]
type = TimeDerivative
variable = u
[../]
[]
[AuxKernels]
[./elem]
type = UniqueIDAux
variable = elem
execute_on = timestep_begin
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[UserObjects]
[./rh_uo]
type = RandomHitUserObject
execute_on = timestep_begin
num_hits = 1
[../]
[./rhsm]
type = RandomHitSolutionModifier
execute_on = custom
modify = u
random_hits = rh_uo
amount = 10
[../]
[]
[Executioner]
type = AdaptAndModify
num_steps = 400
dt = 2e-4
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
adapt_cycles = 2
[]
[Adaptivity]
marker = rhm # Switch to combo to get the effect of both
[./Indicators]
[./gji]
type = GradientJumpIndicator
variable = u
[../]
[../]
[./Markers]
[./rhm]
type = RandomHitMarker
random_hits = rh_uo
[../]
[./efm]
type = ErrorFractionMarker
coarsen = 0.2
indicator = gji
refine = 0.8
[../]
[./combo]
type = ComboMarker
markers = 'efm rhm'
[../]
[../]
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/examples/laser-welding/2d.i)
endtime=5e-4 # s
timestep=${fparse endtime/100} # s
surfacetemp=300 # K
power=190 # W
R=1.8257418583505537e-4 # m
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -.45e-3 # m
xmax = 0.45e-3 # m
ymin = -.9e-4 # m
ymax = 0
nx = 25
ny = 5
displacements = 'disp_x disp_y'
[]
[GlobalParams]
temperature = T
[]
[Variables]
[vel]
family = LAGRANGE_VEC
[]
[T]
[]
[p]
[]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[vel_x_aux]
[InitialCondition]
type = ConstantIC
value = 1e-15
[]
[]
[vel_y_aux]
[InitialCondition]
type = ConstantIC
value = 1e-15
[]
[]
[]
[AuxKernels]
[vel_x_value]
type = VectorVariableComponentAux
variable = vel_x_aux
vector_variable = vel
component = x
[]
[vel_y_value]
type = VectorVariableComponentAux
variable = vel_y_aux
vector_variable = vel
component = y
[]
[]
[ICs]
[T]
type = FunctionIC
variable = T
function = '(${surfacetemp} - 300) / .7e-3 * y + ${surfacetemp}'
[]
[]
[Kernels]
[disp_x]
type = Diffusion
variable = disp_x
[]
[disp_y]
type = Diffusion
variable = disp_y
[]
[mass]
type = INSADMass
variable = p
use_displaced_mesh = true
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
use_displaced_mesh = true
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
use_displaced_mesh = true
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
use_displaced_mesh = true
[]
[momentum_mesh_advection]
type = INSADMomentumMeshAdvection
variable = vel
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
use_displaced_mesh = true
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
pressure = p
integrate_p_by_parts = true
use_displaced_mesh = true
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
material_velocity = relative_velocity
use_displaced_mesh = true
[]
[temperature_time]
type = INSADHeatConductionTimeDerivative
variable = T
use_displaced_mesh = true
[]
[temperature_advection]
type = INSADEnergyAdvection
variable = T
use_displaced_mesh = true
[]
[temperature_mesh_advection]
type = INSADEnergyMeshAdvection
variable = T
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[]
[temperature_conduction]
type = ADHeatConduction
variable = T
thermal_conductivity = 'k'
use_displaced_mesh = true
[]
[temperature_supg]
type = INSADEnergySUPG
variable = T
velocity = vel
use_displaced_mesh = true
[]
[]
[BCs]
[x_no_disp]
type = DirichletBC
variable = disp_x
boundary = 'bottom'
value = 0
[]
[y_no_disp]
type = DirichletBC
variable = disp_y
boundary = 'bottom'
value = 0
[]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'bottom right left'
[]
[T_cold]
type = DirichletBC
variable = T
boundary = 'bottom'
value = 300
[]
[radiation_flux]
type = FunctionRadiativeBC
variable = T
boundary = 'top'
emissivity_function = '1'
Tinfinity = 300
stefan_boltzmann_constant = 5.67e-8
use_displaced_mesh = true
[]
[weld_flux]
type = GaussianEnergyFluxBC
variable = T
boundary = 'top'
P0 = ${power}
R = ${R}
x_beam_coord = '-0.35e-3 +0.7e-3*t/${endtime}'
y_beam_coord = '0'
use_displaced_mesh = true
[]
[vapor_recoil]
type = INSADVaporRecoilPressureMomentumFluxBC
variable = vel
boundary = 'top'
use_displaced_mesh = true
[]
[surface_tension]
type = INSADSurfaceTensionBC
variable = vel
boundary = 'top'
use_displaced_mesh = true
include_gradient_terms = true
[]
[displace_x_top]
type = INSADDisplaceBoundaryBC
boundary = 'top'
variable = 'disp_x'
velocity = 'vel'
component = 0
associated_subdomain = 0
[]
[displace_y_top]
type = INSADDisplaceBoundaryBC
boundary = 'top'
variable = 'disp_y'
velocity = 'vel'
component = 1
associated_subdomain = 0
[]
[displace_x_top_dummy]
type = INSADDummyDisplaceBoundaryIntegratedBC
boundary = 'top'
variable = 'disp_x'
velocity = 'vel'
component = 0
[]
[displace_y_top_dummy]
type = INSADDummyDisplaceBoundaryIntegratedBC
boundary = 'top'
variable = 'disp_y'
velocity = 'vel'
component = 1
[]
[]
[Materials]
[ins_mat]
type = INSADStabilized3Eqn
velocity = vel
pressure = p
temperature = T
use_displaced_mesh = true
[]
[steel]
type = AriaLaserWeld304LStainlessSteel
temperature = T
beta = 1e7
use_displaced_mesh = true
[]
[steel_boundary]
type = AriaLaserWeld304LStainlessSteelBoundary
boundary = 'top'
temperature = T
use_displaced_mesh = true
[]
[const]
type = GenericConstantMaterial
prop_names = 'abs sb_constant'
prop_values = '1 5.67e-8'
use_displaced_mesh = true
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu NONZERO strumpack'
[]
[]
[Executioner]
type = Transient
end_time = ${endtime}
dtmin = 1e-8
dtmax = ${timestep}
petsc_options = '-snes_converged_reason -ksp_converged_reason -options_left'
solve_type = 'NEWTON'
line_search = 'none'
nl_max_its = 12
l_max_its = 100
[TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 7
dt = ${timestep}
linear_iteration_ratio = 1e6
growth_factor = 1.5
[]
[]
[Outputs]
[exodus]
type = Exodus
output_material_properties = true
show_material_properties = 'mu'
[]
checkpoint = true
perf_graph = true
[]
[Debug]
show_var_residual_norms = true
[]
[Adaptivity]
marker = combo
max_h_level = 4
[Indicators]
[error_T]
type = GradientJumpIndicator
variable = T
[]
[error_dispz]
type = GradientJumpIndicator
variable = disp_y
[]
[]
[Markers]
[errorfrac_T]
type = ErrorFractionMarker
refine = 0.4
coarsen = 0.2
indicator = error_T
[]
[errorfrac_dispz]
type = ErrorFractionMarker
refine = 0.4
coarsen = 0.2
indicator = error_dispz
[]
[combo]
type = ComboMarker
markers = 'errorfrac_T errorfrac_dispz'
[]
[]
[]
[Postprocessors]
[num_dofs]
type = NumDOFs
system = 'NL'
[]
[nl]
type = NumNonlinearIterations
[]
[tot_nl]
type = CumulativeValuePostprocessor
postprocessor = 'nl'
[]
[]
(test/tests/markers/combo_marker/combo_marker_block_restricted.i)
[Mesh]
[base]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
elem_type = QUAD4
[]
[subdomain_1]
type = SubdomainBoundingBoxGenerator
input = base
bottom_left = '0.1 0.1 0'
block_id = 1
top_right = '0.4 0.9 0'
[]
[subdomain_2]
type = SubdomainBoundingBoxGenerator
input = subdomain_1
bottom_left = '0.1 0.1 0'
block_id = 2
top_right = '0.4 0.9 0'
location = OUTSIDE
[]
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[]
[BCs]
[left]
type = DirichletBC
variable = u
boundary = 1
value = 0
[]
[right]
type = DirichletBC
variable = u
boundary = 2
value = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Adaptivity]
initial_marker = combo
initial_steps = 1
[Markers]
[box]
type = BoxMarker
bottom_left = '0.0 0. 0'
top_right = '1 0.5 0'
inside = refine
outside = do_nothing
block = '1'
[]
[box2]
type = BoxMarker
bottom_left = '0.5 0.5 0'
top_right = '0.8 0.8 0'
inside = refine
outside = coarsen
block = '2'
[]
[combo]
type = ComboMarker
markers = 'box box2'
[]
[]
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/examples/laser-welding/3d.i)
period=1.25e-3
endtime=${period}
timestep=1.25e-5
surfacetemp=300
sb=5.67e-8
[GlobalParams]
temperature = T
[]
[Mesh]
type = GeneratedMesh
dim = 3
xmin = -.35e-3
xmax = 0.35e-3
ymin = -.35e-3
ymax = .35e-3
zmin = -.7e-3
zmax = 0
nx = 2
ny = 2
nz = 2
displacements = 'disp_x disp_y disp_z'
uniform_refine = 2
[]
[Variables]
[vel]
family = LAGRANGE_VEC
[]
[T]
[]
[p]
[]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[]
[ICs]
[T]
type = FunctionIC
variable = T
function = '(${surfacetemp} - 300) / .7e-3 * z + ${surfacetemp}'
[]
[]
[Kernels]
[disp_x]
type = Diffusion
variable = disp_x
[]
[disp_y]
type = Diffusion
variable = disp_y
[]
[disp_z]
type = Diffusion
variable = disp_z
[]
[mass]
type = INSADMass
variable = p
use_displaced_mesh = true
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
use_displaced_mesh = true
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
use_displaced_mesh = true
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
use_displaced_mesh = true
[]
[momentum_mesh_advection]
type = INSADMomentumMeshAdvection
variable = vel
disp_x = disp_x
disp_y = disp_y
disp_z = disp_z
use_displaced_mesh = true
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
use_displaced_mesh = true
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
pressure = p
integrate_p_by_parts = true
use_displaced_mesh = true
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
material_velocity = relative_velocity
use_displaced_mesh = true
[]
[temperature_time]
type = INSADHeatConductionTimeDerivative
variable = T
use_displaced_mesh = true
[]
[temperature_advection]
type = INSADEnergyAdvection
variable = T
use_displaced_mesh = true
[]
[temperature_mesh_advection]
type = INSADEnergyMeshAdvection
variable = T
disp_x = disp_x
disp_y = disp_y
disp_z = disp_z
use_displaced_mesh = true
[]
[temperature_conduction]
type = ADHeatConduction
variable = T
thermal_conductivity = 'k'
use_displaced_mesh = true
[]
[temperature_supg]
type = INSADEnergySUPG
variable = T
velocity = vel
use_displaced_mesh = true
[]
[]
[BCs]
[x_no_disp]
type = DirichletBC
variable = disp_x
boundary = 'back'
value = 0
[]
[y_no_disp]
type = DirichletBC
variable = disp_y
boundary = 'back'
value = 0
[]
[z_no_disp]
type = DirichletBC
variable = disp_z
boundary = 'back'
value = 0
[]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'bottom right left top back'
[]
[T_cold]
type = DirichletBC
variable = T
boundary = 'back'
value = 300
[]
[radiation_flux]
type = FunctionRadiativeBC
variable = T
boundary = 'front'
emissivity_function = '1'
Tinfinity = 300
stefan_boltzmann_constant = ${sb}
use_displaced_mesh = true
[]
[weld_flux]
type = GaussianEnergyFluxBC
variable = T
boundary = 'front'
P0 = 159.96989792079225
R = 1.8257418583505537e-4
x_beam_coord = '2e-4 * cos(t * 2 * pi / ${period})'
y_beam_coord = '2e-4 * sin(t * 2 * pi / ${period})'
z_beam_coord = 0
use_displaced_mesh = true
[]
[vapor_recoil]
type = INSADVaporRecoilPressureMomentumFluxBC
variable = vel
boundary = 'front'
use_displaced_mesh = true
[]
[surface_tension]
type = INSADSurfaceTensionBC
variable = vel
boundary = 'front'
use_displaced_mesh = true
[]
[displace_x_top]
type = INSADDisplaceBoundaryBC
boundary = 'front'
variable = 'disp_x'
velocity = 'vel'
component = 0
associated_subdomain = 0
[]
[displace_y_top]
type = INSADDisplaceBoundaryBC
boundary = 'front'
variable = 'disp_y'
velocity = 'vel'
component = 1
associated_subdomain = 0
[]
[displace_z_top]
type = INSADDisplaceBoundaryBC
boundary = 'front'
variable = 'disp_z'
velocity = 'vel'
component = 2
associated_subdomain = 0
[]
[displace_x_top_dummy]
type = INSADDummyDisplaceBoundaryIntegratedBC
boundary = 'front'
variable = 'disp_x'
velocity = 'vel'
component = 0
[]
[displace_y_top_dummy]
type = INSADDummyDisplaceBoundaryIntegratedBC
boundary = 'front'
variable = 'disp_y'
velocity = 'vel'
component = 1
[]
[displace_z_top_dummy]
type = INSADDummyDisplaceBoundaryIntegratedBC
boundary = 'front'
variable = 'disp_z'
velocity = 'vel'
component = 2
[]
[]
[Materials]
[ins_mat]
type = INSADStabilized3Eqn
velocity = vel
pressure = p
temperature = T
use_displaced_mesh = true
[]
[steel]
type = AriaLaserWeld304LStainlessSteel
temperature = T
beta = 1e7
use_displaced_mesh = true
[]
[steel_boundary]
type = AriaLaserWeld304LStainlessSteelBoundary
boundary = 'front'
temperature = T
use_displaced_mesh = true
[]
[const]
type = GenericConstantMaterial
prop_names = 'abs sb_constant'
prop_values = '1 ${sb}'
use_displaced_mesh = true
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu NONZERO strumpack'
[]
[]
[Executioner]
type = Transient
end_time = ${endtime}
dtmin = 1e-8
dtmax = ${timestep}
petsc_options = '-snes_converged_reason -ksp_converged_reason -options_left'
solve_type = 'NEWTON'
line_search = 'none'
nl_max_its = 12
l_max_its = 100
[TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 7
dt = ${timestep}
linear_iteration_ratio = 1e6
growth_factor = 1.5
[]
[]
[Outputs]
[exodus]
type = Exodus
output_material_properties = true
show_material_properties = 'mu'
[]
checkpoint = true
perf_graph = true
[]
[Debug]
show_var_residual_norms = true
[]
[Adaptivity]
marker = combo
max_h_level = 4
[Indicators]
[error_T]
type = GradientJumpIndicator
variable = T
[]
[error_dispz]
type = GradientJumpIndicator
variable = disp_z
[]
[]
[Markers]
[errorfrac_T]
type = ErrorFractionMarker
refine = 0.4
coarsen = 0.2
indicator = error_T
[]
[errorfrac_dispz]
type = ErrorFractionMarker
refine = 0.4
coarsen = 0.2
indicator = error_dispz
[]
[combo]
type = ComboMarker
markers = 'errorfrac_T errorfrac_dispz'
[]
[]
[]
[Postprocessors]
[num_dofs]
type = NumDOFs
system = 'NL'
[]
[nl]
type = NumNonlinearIterations
[]
[tot_nl]
type = CumulativeValuePostprocessor
postprocessor = 'nl'
[]
[]
(test/tests/executioners/adapt_and_modify/adapt_and_modify.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 15
ny = 15
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[./td]
type = TimeDerivative
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[UserObjects]
[./rh_uo]
type = RandomHitUserObject
execute_on = 'initial timestep_begin'
num_hits = 1
[../]
[./rhsm]
type = RandomHitSolutionModifier
execute_on = 'custom'
modify = u
random_hits = rh_uo
amount = 1000
[../]
[]
[Executioner]
type = AdaptAndModify
num_steps = 4
dt = 1e-3
solve_type = 'PJFNK'
nl_rel_tol = 1e-15
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
adapt_cycles = 2
[]
[Adaptivity]
marker = rhm # Switch to combo to get the effect of both
[./Indicators]
[./gji]
type = GradientJumpIndicator
variable = u
[../]
[../]
[./Markers]
[./rhm]
type = RandomHitMarker
random_hits = rh_uo
[../]
[./efm]
type = ErrorFractionMarker
coarsen = 0.001
indicator = gji
refine = 0.8
[../]
[./combo]
type = ComboMarker
markers = 'efm rhm'
[../]
[../]
[]
[Outputs]
exodus = true
[]
(modules/phase_field/examples/nucleation/refine.i)
#
# Example derived from cahn_hilliard.i demonstrating the use of Adaptivity
# with the DiscreteNucleation system. The DiscreteNucleationMarker triggers
# mesh refinement for the nucleus geometry. It is up to the user to specify
# refinement for the physics. In this example this is done using a GradientJumpIndicator
# with a ValueThresholdMarker. The nucleation system marker and the physics marker
# must be combined using a ComboMarker to combine their effect.
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
xmax = 500
ymax = 500
elem_type = QUAD
[]
[Modules]
[./PhaseField]
[./Conserved]
[./c]
free_energy = F
mobility = M
kappa = kappa_c
solve_type = REVERSE_SPLIT
[../]
[../]
[../]
[]
[ICs]
[./c_IC]
type = ConstantIC
variable = c
value = 0.2
[../]
[]
[Materials]
[./pfmobility]
type = GenericConstantMaterial
prop_names = 'M kappa_c'
prop_values = '1 25'
[../]
[./chemical_free_energy]
# simple double well free energy
type = DerivativeParsedMaterial
property_name = Fc
coupled_variables = 'c'
constant_names = 'barr_height cv_eq'
constant_expressions = '0.1 0'
expression = 16*barr_height*c^2*(1-c)^2 # +0.01*(c*plog(c,0.005)+(1-c)*plog(1-c,0.005))
derivative_order = 2
outputs = exodus
[../]
[./probability]
# This is a made up toy nucleation rate it should be replaced by
# classical nucleation theory in a real simulation.
type = ParsedMaterial
property_name = P
coupled_variables = c
expression = 'if(c<0.21,c*1e-8,0)'
outputs = exodus
[../]
[./nucleation]
# The nucleation material is configured to insert nuclei into the free energy
# tht force the concentration to go to 0.95, and holds this enforcement for 500
# time units.
type = DiscreteNucleation
property_name = Fn
op_names = c
op_values = 0.90
penalty = 5
penalty_mode = MIN
map = map
outputs = exodus
[../]
[./free_energy]
# add the chemical and nucleation free energy contributions together
type = DerivativeSumMaterial
derivative_order = 2
coupled_variables = c
sum_materials = 'Fc Fn'
[../]
[]
[UserObjects]
[./inserter]
# The inserter runs at the end of each time step to add nucleation events
# that happened during the timestep (if it converged) to the list of nuclei
type = DiscreteNucleationInserter
hold_time = 50
probability = P
radius = 10
[../]
[./map]
# The map UO runs at the beginning of a timestep and generates a per-element/qp
# map of nucleus locations. The map is only regenerated if the mesh changed or
# the list of nuclei was modified.
# The map converts the nucleation points into finite area objects with a given radius.
type = DiscreteNucleationMap
periodic = c
inserter = inserter
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[BCs]
[./Periodic]
[./all]
auto_direction = 'x y'
[../]
[../]
[]
[Postprocessors]
[./dt]
type = TimestepSize
[../]
[./ndof]
type = NumDOFs
[../]
[./rate]
type = DiscreteNucleationData
value = RATE
inserter = inserter
[../]
[./dtnuc]
type = DiscreteNucleationTimeStep
inserter = inserter
p2nucleus = 0.0005
dt_max = 10
[../]
[./update]
type = DiscreteNucleationData
value = UPDATE
inserter = inserter
[../]
[./count]
type = DiscreteNucleationData
value = COUNT
inserter = inserter
[../]
[]
[Adaptivity]
[./Indicators]
[./jump]
type = GradientJumpIndicator
variable = c
[../]
[../]
[./Markers]
[./nuc]
type = DiscreteNucleationMarker
map = map
[../]
[./grad]
type = ValueThresholdMarker
variable = jump
coarsen = 0.1
refine = 0.2
[../]
[./combo]
type = ComboMarker
markers = 'nuc grad'
[../]
[../]
marker = combo
cycles_per_step = 3
recompute_markers_during_cycles = true
max_h_level = 3
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'asm lu '
nl_max_its = 20
l_tol = 1.0e-4
nl_rel_tol = 1.0e-10
nl_abs_tol = 1.0e-10
start_time = 0.0
num_steps = 120
[./TimeStepper]
type = IterationAdaptiveDT
dt = 10
growth_factor = 1.5
cutback_factor = 0.5
optimal_iterations = 8
iteration_window = 2
timestep_limiting_postprocessor = dtnuc
[../]
[]
[Outputs]
exodus = true
csv = true
print_linear_residuals = false
[]
(test/tests/markers/combo_marker/combo_marker_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
nz = 0
zmax = 0
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = 1
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = 2
value = 1
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Adaptivity]
[./Markers]
[./box]
type = BoxMarker
bottom_left = '0.3 0.3 0'
top_right = '0.6 0.6 0'
inside = refine
outside = do_nothing
[../]
[./combo]
type = ComboMarker
markers = 'box box2'
[../]
[./box2]
type = BoxMarker
bottom_left = '0.5 0.5 0'
top_right = '0.8 0.8 0'
inside = refine
outside = coarsen
[../]
[../]
[]
[Outputs]
exodus = true
[]