- diffusion_coefficientThe name of the diffusivity, can be scalar, vector, or matrix material property.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The name of the diffusivity, can be scalar, vector, or matrix material property.
- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
ArrayDiffusion
Description
This array kernel implements the following piece of a weak form:
which is expanded as where is the test function, is the finite element solution and is the diffusion coefficients. is an array variable that has number of components. is the scalings of all components of the array variable. The size of the vector test function is the same as the size of . The kernel can be further spelled out as (1) where denotes elements of the mesh; is the number of shape functions on the local element; is the number of quadrature points for doing the spatial integration over an element. and are the expansion coefficients for the test function and the solution respectively, often referred as degrees of freedom. Subscript indicates that the associated quantities are evaluated on a quadrature point. is the determinant of Jacobian, that transform a physical element to a reference element where the shape functions are defined, times local quadrature weighting. It is noted that the test functions of all components are identical because all components of an array variable share the same finite element family and order.
We can rearrange the fully expanded Eq. (1) into where the underlined term is the array provided by ArrayDiffusion::computeQpResidual. The element, shape function and quadrature point summations are taken care of by MOOSE. It is noted that since test functions are arbitrary, can be viewed as an index indicator with which the local residual goes into the global residual vector. Note that represents element-wise multiplication, i.e. is equal to vector whose ith element is , where and are two generic vectors.
In general, the diffusion coefficient is a square matrix with the size of the number of components. When it is a diagonal matrix, it can be represented by an array. In such a case, the components are not coupled with this array diffusion kernel. If all elements of the diffusion coefficient vector are the same, we can use a scalar diffusion coefficient. Thus this kernel gives users an option to set the type of diffusion coefficient with a parameter named as diffusion_coefficient_type. Users can set it to scalar, array or full corresponding to scalar, diagonal matrix and full matrix respectively. Its default value is array.
With some further transformation, the kernel becomes where the underlined part is the local Jacobian evaluated by ArrayDiffusion::computeQpJacobian and ArrayDiffusion::computeQpOffDiagJacobian.
Example Input Syntax
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[]
(test/tests/kernels/array_kernels/array_diffusion_test.i)Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Unit:(no unit assumed)
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacements
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Unit:(no unit assumed)
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Unit:(no unit assumed)
Controllable:No
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Unit:(no unit assumed)
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
Input Files
- (test/tests/kernels/array_kernels/array_diffusion_reaction_other_coupling.i)
- (test/tests/transfers/multiapp_copy_transfer/array_variable_transfer/sub.i)
- (test/tests/multiapps/relaxation/picard_relaxed_array_parent.i)
- (test/tests/transfers/multiapp_variable_value_sample_transfer/parent_array.i)
- (test/tests/kernels/array_kernels/array_body_force.i)
- (test/tests/kernels/array_kernels/array_diffusion_test.i)
- (test/tests/kernels/hfem/robin_adapt.i)
- (test/tests/kernels/hfem/array_neumann.i)
- (test/tests/kernels/array_kernels/array_diffusion_reaction_dg.i)
- (test/tests/kernels/hfem/array_robin.i)
- (test/tests/kernels/hfem/variable_robin.i)
- (test/tests/kernels/array_kernels/array_custom_coupling_test.i)
- (test/tests/problems/eigen_problem/arraykernels/ne_array_kernels.i)
- (test/tests/kernels/hfem/array_dirichlet_transform.i)
- (test/tests/kernels/array_kernels/array_diffusion_reaction_transient.i)
- (test/tests/kernels/hfem/array_dirichlet_transform_bc.i)
- (test/tests/kernels/hfem/variable_dirichlet.i)
- (test/tests/bcs/array_vacuum/array_vacuum.i)
- (test/tests/tag/coupled_array_grad_dot.i)
- (test/tests/kernels/hfem/array_dirichlet_pjfnk.i)
- (test/tests/variables/array_variable/array_variable_size_one_test.i)
- (test/tests/preconditioners/fsp/array-test.i)
- (test/tests/tag/tag-array-grad.i)
- (test/tests/kernels/hfem/robin_displaced.i)
- (modules/optimization/test/tests/executioners/steady_and_adjoint/array_variable.i)
- (test/tests/kernels/array_coupled_time_derivative/test_jacobian.i)
- (test/tests/scaling/array-variables/test.i)
- (test/tests/kernels/hfem/robin_dist.i)
- (test/tests/problems/eigen_problem/jfnk_mo/ne_mo_with_linear_aux.i)
- (test/tests/tag/tag-array-var.i)
- (test/tests/problems/eigen_problem/jfnk_mo/ne_array_mo.i)
- (test/tests/auxkernels/array_var_component/array_var_component.i)
- (test/tests/kernels/array_kernels/array_save_in.i)
- (test/tests/kernels/array_kernels/array_diffusion_reaction.i)
- (test/tests/kernels/array_kernels/array_diffusion_reaction_coupling.i)
- (test/tests/kernels/hfem/array_dirichlet.i)
- (test/tests/multiapps/relaxation/picard_relaxed_array_sub.i)
- (test/tests/interfaces/coupleable/array_coupling_by_name.i)
(framework/src/kernels/ArrayDiffusion.C)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#include "ArrayDiffusion.h"
registerMooseObject("MooseApp", ArrayDiffusion);
InputParameters
ArrayDiffusion::validParams()
{
InputParameters params = ArrayKernel::validParams();
params.addRequiredParam<MaterialPropertyName>(
"diffusion_coefficient",
"The name of the diffusivity, can be scalar, vector, or matrix material property.");
params.addClassDescription(
"The array Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak "
"form of $(\\nabla \\phi_i, \\nabla u_h)$.");
return params;
}
ArrayDiffusion::ArrayDiffusion(const InputParameters & parameters)
: ArrayKernel(parameters),
_d(hasMaterialProperty<Real>("diffusion_coefficient")
? &getMaterialProperty<Real>("diffusion_coefficient")
: nullptr),
_d_array(hasMaterialProperty<RealEigenVector>("diffusion_coefficient")
? &getMaterialProperty<RealEigenVector>("diffusion_coefficient")
: nullptr),
_d_2d_array(hasMaterialProperty<RealEigenMatrix>("diffusion_coefficient")
? &getMaterialProperty<RealEigenMatrix>("diffusion_coefficient")
: nullptr)
{
if (!_d && !_d_array && !_d_2d_array)
{
MaterialPropertyName mat = getParam<MaterialPropertyName>("diffusion_coefficient");
mooseError("Property " + mat + " is of unsupported type for ArrayDiffusion");
}
}
void
ArrayDiffusion::initQpResidual()
{
if (_d_array)
{
mooseAssert((*_d_array)[_qp].size() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
}
else if (_d_2d_array)
{
mooseAssert((*_d_2d_array)[_qp].cols() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
mooseAssert((*_d_2d_array)[_qp].rows() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
}
}
void
ArrayDiffusion::computeQpResidual(RealEigenVector & residual)
{
// WARNING: the noalias() syntax is an Eigen optimization tactic, it avoids creating
// a temporary object for the matrix multiplication on the right-hand-side. However,
// it should be used with caution because it could cause unintended results,
// developers should NOT use it if the vector on the left-hand-side appears on the
// right-hand-side, for instance:
// vector = matrix * vector;
// See http://eigen.tuxfamily.org/dox/group__TopicAliasing.html for more details.
if (_d)
residual.noalias() = (*_d)[_qp] * _grad_u[_qp] * _array_grad_test[_i][_qp];
else if (_d_array)
residual.noalias() = (*_d_array)[_qp].asDiagonal() * _grad_u[_qp] * _array_grad_test[_i][_qp];
else
residual.noalias() = (*_d_2d_array)[_qp] * _grad_u[_qp] * _array_grad_test[_i][_qp];
}
RealEigenVector
ArrayDiffusion::computeQpJacobian()
{
if (_d)
return RealEigenVector::Constant(_var.count(),
_grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d)[_qp]);
else if (_d_array)
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_array)[_qp];
else
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_2d_array)[_qp].diagonal();
}
RealEigenMatrix
ArrayDiffusion::computeQpOffDiagJacobian(const MooseVariableFEBase & jvar)
{
if (jvar.number() == _var.number() && _d_2d_array)
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_2d_array)[_qp];
else
return ArrayKernel::computeQpOffDiagJacobian(jvar);
}
(framework/src/kernels/ArrayDiffusion.C)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#include "ArrayDiffusion.h"
registerMooseObject("MooseApp", ArrayDiffusion);
InputParameters
ArrayDiffusion::validParams()
{
InputParameters params = ArrayKernel::validParams();
params.addRequiredParam<MaterialPropertyName>(
"diffusion_coefficient",
"The name of the diffusivity, can be scalar, vector, or matrix material property.");
params.addClassDescription(
"The array Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak "
"form of $(\\nabla \\phi_i, \\nabla u_h)$.");
return params;
}
ArrayDiffusion::ArrayDiffusion(const InputParameters & parameters)
: ArrayKernel(parameters),
_d(hasMaterialProperty<Real>("diffusion_coefficient")
? &getMaterialProperty<Real>("diffusion_coefficient")
: nullptr),
_d_array(hasMaterialProperty<RealEigenVector>("diffusion_coefficient")
? &getMaterialProperty<RealEigenVector>("diffusion_coefficient")
: nullptr),
_d_2d_array(hasMaterialProperty<RealEigenMatrix>("diffusion_coefficient")
? &getMaterialProperty<RealEigenMatrix>("diffusion_coefficient")
: nullptr)
{
if (!_d && !_d_array && !_d_2d_array)
{
MaterialPropertyName mat = getParam<MaterialPropertyName>("diffusion_coefficient");
mooseError("Property " + mat + " is of unsupported type for ArrayDiffusion");
}
}
void
ArrayDiffusion::initQpResidual()
{
if (_d_array)
{
mooseAssert((*_d_array)[_qp].size() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
}
else if (_d_2d_array)
{
mooseAssert((*_d_2d_array)[_qp].cols() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
mooseAssert((*_d_2d_array)[_qp].rows() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
}
}
void
ArrayDiffusion::computeQpResidual(RealEigenVector & residual)
{
// WARNING: the noalias() syntax is an Eigen optimization tactic, it avoids creating
// a temporary object for the matrix multiplication on the right-hand-side. However,
// it should be used with caution because it could cause unintended results,
// developers should NOT use it if the vector on the left-hand-side appears on the
// right-hand-side, for instance:
// vector = matrix * vector;
// See http://eigen.tuxfamily.org/dox/group__TopicAliasing.html for more details.
if (_d)
residual.noalias() = (*_d)[_qp] * _grad_u[_qp] * _array_grad_test[_i][_qp];
else if (_d_array)
residual.noalias() = (*_d_array)[_qp].asDiagonal() * _grad_u[_qp] * _array_grad_test[_i][_qp];
else
residual.noalias() = (*_d_2d_array)[_qp] * _grad_u[_qp] * _array_grad_test[_i][_qp];
}
RealEigenVector
ArrayDiffusion::computeQpJacobian()
{
if (_d)
return RealEigenVector::Constant(_var.count(),
_grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d)[_qp]);
else if (_d_array)
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_array)[_qp];
else
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_2d_array)[_qp].diagonal();
}
RealEigenMatrix
ArrayDiffusion::computeQpOffDiagJacobian(const MooseVariableFEBase & jvar)
{
if (jvar.number() == _var.number() && _d_2d_array)
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_2d_array)[_qp];
else
return ArrayKernel::computeQpOffDiagJacobian(jvar);
}
(framework/src/kernels/ArrayDiffusion.C)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#include "ArrayDiffusion.h"
registerMooseObject("MooseApp", ArrayDiffusion);
InputParameters
ArrayDiffusion::validParams()
{
InputParameters params = ArrayKernel::validParams();
params.addRequiredParam<MaterialPropertyName>(
"diffusion_coefficient",
"The name of the diffusivity, can be scalar, vector, or matrix material property.");
params.addClassDescription(
"The array Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak "
"form of $(\\nabla \\phi_i, \\nabla u_h)$.");
return params;
}
ArrayDiffusion::ArrayDiffusion(const InputParameters & parameters)
: ArrayKernel(parameters),
_d(hasMaterialProperty<Real>("diffusion_coefficient")
? &getMaterialProperty<Real>("diffusion_coefficient")
: nullptr),
_d_array(hasMaterialProperty<RealEigenVector>("diffusion_coefficient")
? &getMaterialProperty<RealEigenVector>("diffusion_coefficient")
: nullptr),
_d_2d_array(hasMaterialProperty<RealEigenMatrix>("diffusion_coefficient")
? &getMaterialProperty<RealEigenMatrix>("diffusion_coefficient")
: nullptr)
{
if (!_d && !_d_array && !_d_2d_array)
{
MaterialPropertyName mat = getParam<MaterialPropertyName>("diffusion_coefficient");
mooseError("Property " + mat + " is of unsupported type for ArrayDiffusion");
}
}
void
ArrayDiffusion::initQpResidual()
{
if (_d_array)
{
mooseAssert((*_d_array)[_qp].size() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
}
else if (_d_2d_array)
{
mooseAssert((*_d_2d_array)[_qp].cols() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
mooseAssert((*_d_2d_array)[_qp].rows() == _var.count(),
"diffusion_coefficient size is inconsistent with the number of components of array "
"variable");
}
}
void
ArrayDiffusion::computeQpResidual(RealEigenVector & residual)
{
// WARNING: the noalias() syntax is an Eigen optimization tactic, it avoids creating
// a temporary object for the matrix multiplication on the right-hand-side. However,
// it should be used with caution because it could cause unintended results,
// developers should NOT use it if the vector on the left-hand-side appears on the
// right-hand-side, for instance:
// vector = matrix * vector;
// See http://eigen.tuxfamily.org/dox/group__TopicAliasing.html for more details.
if (_d)
residual.noalias() = (*_d)[_qp] * _grad_u[_qp] * _array_grad_test[_i][_qp];
else if (_d_array)
residual.noalias() = (*_d_array)[_qp].asDiagonal() * _grad_u[_qp] * _array_grad_test[_i][_qp];
else
residual.noalias() = (*_d_2d_array)[_qp] * _grad_u[_qp] * _array_grad_test[_i][_qp];
}
RealEigenVector
ArrayDiffusion::computeQpJacobian()
{
if (_d)
return RealEigenVector::Constant(_var.count(),
_grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d)[_qp]);
else if (_d_array)
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_array)[_qp];
else
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_2d_array)[_qp].diagonal();
}
RealEigenMatrix
ArrayDiffusion::computeQpOffDiagJacobian(const MooseVariableFEBase & jvar)
{
if (jvar.number() == _var.number() && _d_2d_array)
return _grad_phi[_j][_qp] * _grad_test[_i][_qp] * (*_d_2d_array)[_qp];
else
return ArrayKernel::computeQpOffDiagJacobian(jvar);
}
(test/tests/kernels/array_kernels/array_diffusion_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
scaling = '.9 .9'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/array_kernels/array_diffusion_reaction_other_coupling.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[v]
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[diffv]
type = Diffusion
variable = v
[]
[vu]
type = ArrayCoupledForce
variable = u
v = v
coef = '0 0.5'
[]
[uv]
type = CoupledArrayForce
variable = v
v = u
coef = '0.05 0'
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[leftv]
type = DirichletBC
variable = v
boundary = 1
value = 0
[]
[rightv]
type = DirichletBC
variable = v
boundary = 2
value = 2
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Postprocessors]
[intu0]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 0
[]
[intu1]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 1
[]
[intv]
type = ElementIntegralVariablePostprocessor
variable = v
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/transfers/multiapp_copy_transfer/array_variable_transfer/sub.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Executioner]
type = Transient
num_steps = 1
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/multiapps/relaxation/picard_relaxed_array_parent.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
[]
[]
[Variables]
[u]
components = 2
[]
[]
[AuxVariables]
[v]
components = 2
initial_condition = '1 1'
[]
[inverse_v]
components = 2
initial_condition = '1 1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[time]
type = ArrayTimeDerivative
variable = u
time_derivative_coefficient = tc
[]
[force_u]
type = ArrayCoupledForce
variable = u
v = inverse_v
is_v_array = true
coef = '1 1'
[]
[]
[AuxKernels]
[invert_v]
type = ArrayQuotientAux
variable = inverse_v
denominator = v
numerator = '20 20'
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = left
values = '0 0'
[]
[Neumann_right]
type = ArrayNeumannBC
variable = u
boundary = right
value = '1 1'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '0.1 0.1'
[]
[tc]
type = GenericConstantArray
prop_name = tc
prop_value = '1 1'
[]
[]
[Postprocessors]
[picard_its]
type = NumFixedPointIterations
execute_on = 'initial timestep_end'
[]
[]
[Executioner]
type = Transient
num_steps = 4
dt = 0.5
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
fixed_point_max_its = 30
nl_abs_tol = 1e-14
relaxation_factor = 0.8
transformed_variables = u
[]
[Outputs]
exodus = true
execute_on = 'INITIAL TIMESTEP_END'
[]
[MultiApps]
[sub]
type = TransientMultiApp
execute_on = timestep_begin
input_files = picard_relaxed_array_sub.i
[]
[]
[Transfers]
[v_from_sub]
type = MultiAppCopyTransfer
from_multi_app = sub
source_variable = v
variable = v
[]
[u_to_sub]
type = MultiAppCopyTransfer
to_multi_app = sub
source_variable = u
variable = u
[]
[]
(test/tests/transfers/multiapp_variable_value_sample_transfer/parent_array.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[source]
type = ArrayBodyForce
variable = u
function = '1 x'
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 'left right bottom top'
values = '0 0'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 1
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
[MultiApps]
[pp_sub0]
type = TransientMultiApp
input_files = pp_sub.i
execute_on = timestep_end
positions = '0.5 0.5 0 0.7 0.7 0'
[]
[pp_sub1]
type = TransientMultiApp
input_files = pp_sub.i
execute_on = timestep_end
positions = '0.5 0.5 0 0.7 0.7 0'
[]
[]
[Transfers]
[sample_pp_transfer0]
type = MultiAppVariableValueSamplePostprocessorTransfer
to_multi_app = pp_sub0
postprocessor = from_parent
source_variable = u
source_variable_component = 0
[]
[sample_pp_transfer1]
type = MultiAppVariableValueSamplePostprocessorTransfer
to_multi_app = pp_sub1
postprocessor = from_parent
source_variable = u
source_variable_component = 1
[]
[]
(test/tests/kernels/array_kernels/array_body_force.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
scaling = '.9 .9'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[source]
type = ArrayBodyForce
variable = u
function = '1 x'
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 'left right bottom top'
values = '0 0'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/array_kernels/array_diffusion_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
scaling = '.9 .9'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/hfem/robin_adapt.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayVacuumBC
boundary = 'left right top bottom'
variable = u
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
exodus = true
[]
[Adaptivity]
steps = 1
marker = box
max_h_level = 2
initial_steps = 2
[Markers]
[box]
bottom_left = '0 0 0'
inside = refine
top_right = '0.5 0.5 0'
outside = do_nothing
type = BoxMarker
[]
[]
[]
(test/tests/kernels/hfem/array_neumann.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
block = 0
reaction_coefficient = rc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayNeumannBC
boundary = 'left right top bottom'
variable = u
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstantArray
prop_name = rc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/kernels/array_kernels/array_diffusion_reaction_dg.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 4
ny = 4
[]
[subdomain1]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0.5 0.5 0'
top_right = '1 1 0'
block_id = 1
[]
[]
[Variables]
[u]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[DGKernels]
[dgdiff]
type = ArrayDGDiffusion
variable = u
diff = dc
[]
[]
[BCs]
[left]
type = ArrayVacuumBC
variable = u
boundary = 1
[]
[right]
type = ArrayPenaltyDirichletBC
variable = u
boundary = 2
value = '1 2'
penalty = 4
[]
[]
[Materials]
[dc0]
type = GenericConstantArray
block = 0
prop_name = dc
prop_value = '1 1'
[]
[dc1]
type = GenericConstantArray
block = 1
prop_name = dc
prop_value = '2 1'
[]
[rc]
type = GenericConstant2DArray
block = '0 1'
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Postprocessors]
[intu0]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 0
[]
[intu1]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/hfem/array_robin.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayVacuumBC
boundary = 'left right top bottom'
variable = u
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/kernels/hfem/variable_robin.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[uhat]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[lambdab]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[uhat_reaction]
type = ArrayReaction
variable = uhat
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
reaction_coefficient = rc
[]
[uhat_coupled]
type = ArrayCoupledForce
variable = uhat
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
v = lambdab
is_v_array = true
coef = '1 1'
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayHFEMDirichletBC
boundary = 'left right top bottom'
variable = u
lowerd_variable = lambdab
uhat = uhat
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstantArray
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
prop_name = rc
prop_value = '0.5 0.5'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/kernels/array_kernels/array_custom_coupling_test.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 4
ny = 4
dim = 2
[]
[]
[Variables]
[u]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[DGKernels]
[dgdiff]
type = ArrayDGDiffusion
variable = u
diff = dc
[]
[]
[BCs]
[left]
type = ArrayVacuumBC
variable = u
boundary = 1
[]
[right]
type = ArrayPenaltyDirichletBC
variable = u
boundary = 2
value = '1 2'
penalty = 4
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[pbp]
type = PBP
solve_order = 'u_0 u_1'
preconditioner = 'AMG AMG'
off_diag_row = 'u_0 u_1'
off_diag_column = 'u_0 u_1'
[]
[]
[Executioner]
type = Steady
solve_type = JFNK
petsc_options = '-mat_view'
[]
[Outputs]
exodus = true
[]
(test/tests/problems/eigen_problem/arraykernels/ne_array_kernels.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 10
elem_type = QUAD4
nx = 8
ny = 8
[]
# the minimum eigenvalue of this problem is 2*(PI/a)^2;
# Its inverse is 0.5*(a/PI)^2 = 5.0660591821169. Here a is equal to 10.
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[v]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[diffu]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reactionu]
type = ArrayReaction
variable = u
reaction_coefficient = rc
extra_vector_tags = 'eigen'
[]
[diffv]
type = ArrayDiffusion
variable = v
diffusion_coefficient = dc
[]
[reactionv]
type = ArrayReaction
variable = v
reaction_coefficient = rc
extra_vector_tags = 'eigen'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1. 1.'
[]
[rc]
type = GenericConstantArray
prop_name = rc
prop_value = '-1 -1'
[]
[]
[BCs]
[hom_u]
type = ArrayDirichletBC
variable = u
values = '0 0'
boundary = '0 1 2 3'
[]
[eigenhom_u]
type = EigenArrayDirichletBC
variable = u
boundary = '0 1 2 3'
[]
[hom_v]
type = ArrayDirichletBC
variable = v
values = '0 0'
boundary = '0 1 2 3'
[]
[eigenhom_v]
type = EigenArrayDirichletBC
variable = v
boundary = '0 1 2 3'
[]
[]
[Executioner]
type = Eigenvalue
[]
[VectorPostprocessors]
[./eigenvalues]
type = Eigenvalues
execute_on = 'timestep_end'
[../]
[]
[Outputs]
csv = true
execute_on = 'timestep_end'
[]
(test/tests/kernels/hfem/array_dirichlet_transform.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[lambdab]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusionTest
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayHFEMDirichletBC
boundary = 'left right top bottom'
variable = u
lowerd_variable = lambdab
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/kernels/array_kernels/array_diffusion_reaction_transient.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[dudt]
type = ArrayTimeDerivative
variable = u
time_derivative_coefficient = tc
[]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[tc]
type = GenericConstantArray
prop_name = tc
prop_value = '1 1'
[]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Postprocessors]
[intu0]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 0
[]
[intu1]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 1
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
dt = 0.1
num_steps = 10
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/hfem/array_dirichlet_transform_bc.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[lambdab]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayHFEMDirichletTestBC
boundary = 'left right top bottom'
variable = u
lowerd_variable = lambdab
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/kernels/hfem/variable_dirichlet.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[lambdab]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[uhat]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayHFEMDirichletBC
boundary = 'left right top bottom'
variable = u
lowerd_variable = lambdab
uhat = uhat
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/bcs/array_vacuum/array_vacuum.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[BCs]
[left]
type = ArrayVacuumBC
variable = u
boundary = 1
alpha = '1 1.2'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/tag/coupled_array_grad_dot.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 3
[]
[Variables]
[u]
order = FIRST
components = 2
family = L2_LAGRANGE
[]
[v]
order = FIRST
components = 2
family = L2_LAGRANGE
[]
[]
[Kernels]
[u_coupled_time_derivative]
type = ArrayCoupledTimeDerivative
variable = u
v = v
[]
[u_time_derivative]
type = ArrayTimeDerivative
variable = u
[]
[u_diffusion]
type = ArrayDiffusion
variable = u
diffusion_coefficient = u_dc
[]
[v_time_derivative]
type = ArrayTimeDerivative
variable = v
[]
[v_diffusion]
type = ArrayDiffusion
variable = v
diffusion_coefficient = v_dc
[]
[]
[ICs]
[u]
type = ArrayFunctionIC
variable = u
function = '2*(x+1) 3*(x+1)'
[]
[v]
type = ArrayFunctionIC
variable = v
function = '0.1*(x+1) 0.2*(x+1)'
[]
[]
[Materials]
[u_dc]
type = GenericConstantArray
prop_name = u_dc
prop_value = '1 1'
[]
[v_dc]
type = GenericConstantArray
prop_name = v_dc
prop_value = '2 2'
[]
[]
[AuxVariables]
[u_grad_dot_x]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[]
[AuxKernels]
[u_grad_dot_x]
type = CoupledArrayGradDotAux
variable = u_grad_dot_x
v = u
grad_component = x
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
dt = 0.1
dtmin = 0.1
num_steps = 3
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/hfem/array_dirichlet_pjfnk.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[lambdab]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
block = 0
reaction_coefficient = re
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusionTest
variable = u
lowerd_variable = lambda
for_pjfnk = true
[]
[]
[BCs]
[all]
type = ArrayHFEMDirichletTestBC
boundary = 'left right top bottom'
variable = u
lowerd_variable = lambdab
for_pjfnk = true
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[re]
type = GenericConstantArray
prop_name = re
prop_value = '0.1 0.1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/variables/array_variable/array_variable_size_one_test.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
[]
[Variables]
[u]
order = FIRST
family = MONOMIAL
components = 1
array = true
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[source]
type = ArrayBodyForce
variable = u
function = '1'
[]
[]
[BCs]
[all]
type = ArrayVacuumBC
boundary = 'left right top bottom'
variable = u
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1'
[]
[rc]
type = GenericConstantArray
prop_name = rc
prop_value = '2'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
[out]
type = Exodus
[]
[]
(test/tests/preconditioners/fsp/array-test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[v]
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[diffv]
type = Diffusion
variable = v
[]
[vu]
type = ArrayCoupledForce
variable = u
v = v
coef = '0 0.5'
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[leftv]
type = DirichletBC
variable = v
boundary = 1
value = 0
[]
[rightv]
type = DirichletBC
variable = v
boundary = 2
value = 2
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[FSP]
type = FSP
topsplit = 'uv'
[uv]
splitting = 'u v'
# Generally speaking, there are four types of splitting we could choose
# <additive,multiplicative,symmetric_multiplicative,schur>
splitting_type = symmetric_multiplicative
[]
[u]
vars = 'u'
petsc_options_iname = '-pc_type -ksp_type'
petsc_options_value = ' hypre preonly'
[]
[v]
vars = 'v'
petsc_options_iname = '-pc_type -ksp_type'
petsc_options_value = ' hypre preonly'
[]
[]
[]
[Postprocessors]
[intu0]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 0
[]
[intu1]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 1
[]
[intv]
type = ElementIntegralVariablePostprocessor
variable = v
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/tag/tag-array-grad.i)
[Mesh]
[gen]
type = CartesianMeshGenerator
dim = 2
dx = '1 1'
ix = '2 2'
dy = '1 1'
iy = '2 2'
subdomain_id = '0 0 0 1'
[]
[]
[Variables]
[u]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[DGKernels]
[dgdiff]
type = ArrayDGDiffusion
variable = u
diff = dc
[]
[]
[BCs]
[left]
type = ArrayVacuumBC
variable = u
boundary = 1
[]
[right]
type = ArrayPenaltyDirichletBC
variable = u
boundary = 2
value = '1 2'
penalty = 4
[]
[]
[Materials]
[dc0]
type = GenericConstantArray
block = 0
prop_name = dc
prop_value = '1 1'
[]
[dc1]
type = GenericConstantArray
block = 1
prop_name = dc
prop_value = '2 1'
[]
[rc]
type = GenericConstant2DArray
block = '0 1'
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[AuxVariables]
[u_tag_x]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[u_tag_y]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[]
[AuxKernels]
[u_tag_x]
type = TagVectorArrayVariableGradientAux
variable = u_tag_x
v = u
grad_component = x
vector_tag = 'SOLUTION'
[]
[u_tag_y]
type = TagVectorArrayVariableGradientAux
variable = u_tag_y
v = u
grad_component = y
vector_tag = 'NONTIME'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/hfem/robin_displaced.i)
[Mesh]
[square]
type = CartesianMeshGenerator
dx = '0.3125 0.3125 0.3125'
dy = '0.3333333333333 0.3333333333333 0.3333333333333'
dim = 2
[]
displacements = 'x_disp y_disp'
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[x_disp]
block = 0
[]
[y_disp]
block = 0
initial_condition = 0
[]
[]
[AuxKernels]
[x_disp]
type = ParsedAux
variable = x_disp
use_xyzt = true
expression = 'x/15'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
use_displaced_mesh = true
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
use_displaced_mesh = true
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
use_displaced_mesh = true
[]
[]
[BCs]
[all]
type = ArrayVacuumBC
boundary = 'left right top bottom'
variable = u
use_displaced_mesh = true
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
use_displaced_mesh = true
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
use_displaced_mesh = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/array_variable.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
xmax = 1
ymax = 1
nx = 10
ny = 10
[]
[]
[Problem]
nl_sys_names = 'nl0 adjoint'
[]
[Variables]
[u]
components = 2
[]
[u_adjoint]
components = 2
solver_sys = adjoint
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = D
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = A
[]
[src_u]
type = ArrayBodyForce
variable = u
function = '1 0'
[]
[src_adjoint]
type = ArrayBodyForce
variable = u_adjoint
function = '0 1'
[]
[]
[Materials]
[diffc]
type = GenericConstantArray
prop_name = 'D'
prop_value = '1 1'
[]
[coeff]
type = GenericConstant2DArray
prop_name = 'A'
prop_value = '0 -10; -1 0'
[]
[]
[BCs]
[dirichlet]
type = ArrayDirichletBC
variable = u
boundary = 'top right'
values = '0 0'
[]
[]
[Executioner]
type = SteadyAndAdjoint
forward_system = nl0
adjoint_system = adjoint
nl_rel_tol = 1e-12
l_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/array_coupled_time_derivative/test_jacobian.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 3
[]
[Variables]
[u]
components = 2
[]
[v]
components = 2
[]
[]
[Kernels]
[u_coupled_time_derivative]
type = ArrayCoupledTimeDerivative
variable = u
v = v
[]
[u_time_derivative]
type = ArrayTimeDerivative
variable = u
[]
[u_diffusion]
type = ArrayDiffusion
variable = u
diffusion_coefficient = u_dc
[]
[v_time_derivative]
type = ArrayTimeDerivative
variable = v
[]
[v_diffusion]
type = ArrayDiffusion
variable = v
diffusion_coefficient = v_dc
[]
[]
[ICs]
[u]
type = ArrayFunctionIC
variable = u
function = '2*(x+1) 3*(x+1)'
[]
[v]
type = ArrayFunctionIC
variable = v
function = '0.1*(x+1) 0.2*(x+1)'
[]
[]
[Materials]
[u_dc]
type = GenericConstantArray
prop_name = u_dc
prop_value = '1 1'
[]
[v_dc]
type = GenericConstantArray
prop_name = v_dc
prop_value = '2 2'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
dt = 0.1
dtmin = 0.1
num_steps = 3
solve_type = 'NEWTON'
petsc_options = '-snes_test_jacobian -snes_test_jacobian_view'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[]
[Outputs]
exodus = true
[]
(test/tests/scaling/array-variables/test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -200 1000'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
automatic_scaling = true
verbose = true
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/hfem/robin_dist.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
parallel_type = DISTRIBUTED
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayVacuumBC
boundary = 'left right top bottom'
variable = u
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/problems/eigen_problem/jfnk_mo/ne_mo_with_linear_aux.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 10
elem_type = QUAD4
nx = 8
ny = 8
[]
# the minimum eigenvalue of this problem is 2*(PI/a)^2;
# Its inverse is 0.5*(a/PI)^2 = 5.0660591821169. Here a is equal to 10.
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[v]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[AuxVariables]
[w][]
[]
[AuxKernels]
[w]
type = ArrayVariableComponent
variable = w
array_variable = u
component = 0
execute_on = linear
[]
[]
[Kernels]
[diffu]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reactionu]
type = ArrayReaction
variable = u
reaction_coefficient = rc
extra_vector_tags = 'eigen'
[]
[diffv]
type = ArrayDiffusion
variable = v
diffusion_coefficient = dc
[]
[reactionv]
type = ArrayReaction
variable = v
reaction_coefficient = rc
extra_vector_tags = 'eigen'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1. 1.'
[]
[rc]
type = GenericConstantArray
prop_name = rc
prop_value = '-1 -1'
[]
[]
[BCs]
[hom_u]
type = ArrayDirichletBC
variable = u
values = '0 0'
boundary = '0 1 2 3'
[]
[eigenhom_u]
type = EigenArrayDirichletBC
variable = u
boundary = '0 1 2 3'
[]
[hom_v]
type = ArrayDirichletBC
variable = v
values = '0 0'
boundary = '0 1 2 3'
[]
[eigenhom_v]
type = EigenArrayDirichletBC
variable = v
boundary = '0 1 2 3'
[]
[]
[Executioner]
type = Eigenvalue
solve_type = PJFNKMO
constant_matrices = true
[]
[VectorPostprocessors]
[eigenvalues]
type = Eigenvalues
execute_on = 'timestep_end'
[]
[]
[Postprocessors]
[unorm]
type = ElementArrayL2Norm
variable = u
component = 0
[]
[wnorm]
type = ElementL2Norm
variable = w
[]
[]
[Outputs]
csv = true
execute_on = 'timestep_end'
[]
(test/tests/tag/tag-array-var.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[AuxVariables]
[u_tag]
components = 2
[]
[]
[AuxKernels]
[u_tag]
type = TagVectorArrayVariableAux
variable = u_tag
v = u
vector_tag = 'nontime'
[]
[]
[Variables]
[u]
components = 2
[]
[]
[Kernels]
[time]
type = ArrayTimeDerivative
variable = u
[]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
num_steps = 1
[]
[Outputs]
exodus = true
[]
(test/tests/problems/eigen_problem/jfnk_mo/ne_array_mo.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 10
elem_type = QUAD4
nx = 8
ny = 8
[]
[left_block]
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0 0'
top_right = '5 10 0'
block_id = 1
input = gen
[]
[refine]
type = RefineBlockGenerator
input = left_block
refinement = '0 0'
block = '1 0'
[]
[]
# the minimum eigenvalue of this problem is 2*(PI/a)^2;
# Its inverse is 0.5*(a/PI)^2 = 5.0660591821169. Here a is equal to 10.
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[v]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[diffu]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reactionu]
type = ArrayReaction
variable = u
reaction_coefficient = rc
extra_vector_tags = 'eigen'
[]
[diffv]
type = ArrayDiffusion
variable = v
diffusion_coefficient = dc
[]
[reactionv]
type = ArrayReaction
variable = v
reaction_coefficient = rc
extra_vector_tags = 'eigen'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1. 1.'
[]
[rc]
type = GenericConstantArray
prop_name = rc
prop_value = '-1 -1'
[]
[]
[BCs]
[hom_u]
type = ArrayDirichletBC
variable = u
values = '0 0'
boundary = '0 1 2 3'
[]
[eigenhom_u]
type = EigenArrayDirichletBC
variable = u
boundary = '0 1 2 3'
[]
[hom_v]
type = ArrayDirichletBC
variable = v
values = '0 0'
boundary = '0 1 2 3'
[]
[eigenhom_v]
type = EigenArrayDirichletBC
variable = v
boundary = '0 1 2 3'
[]
[]
[Executioner]
type = Eigenvalue
solve_type = PJFNKMO
constant_matrices = true
[]
[VectorPostprocessors]
[./eigenvalues]
type = Eigenvalues
execute_on = 'timestep_end'
[../]
[]
[Outputs]
csv = true
execute_on = 'timestep_end'
[]
(test/tests/auxkernels/array_var_component/array_var_component.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[AuxVariables]
[u0][]
[]
[AuxKernels]
[u0]
type = ArrayVariableComponent
variable = u0
array_variable = u
component = 0
[]
[]
[Postprocessors]
[intu0]
type = ElementIntegralVariablePostprocessor
variable = u0
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/array_kernels/array_save_in.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 4
ny = 4
[]
[subdomain1]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0.5 0.5 0'
top_right = '1 1 0'
block_id = 1
[]
[]
[Variables]
[u]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[]
[AuxVariables]
[u_diff_save_in]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[u_vacuum_save_in]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[u_dg_save_in]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[u_diff_diag_save_in]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[u_vacuum_diag_save_in]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[u_dg_diag_save_in]
order = FIRST
family = L2_LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
save_in = u_diff_save_in
diag_save_in = u_diff_diag_save_in
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[DGKernels]
[dgdiff]
type = ArrayDGDiffusion
variable = u
diff = dc
save_in = u_dg_save_in
diag_save_in = u_dg_diag_save_in
[]
[]
[BCs]
[left]
type = ArrayVacuumBC
variable = u
boundary = 1
save_in = u_vacuum_save_in
diag_save_in = u_vacuum_diag_save_in
[]
[right]
type = ArrayPenaltyDirichletBC
variable = u
boundary = 2
value = '1 2'
penalty = 4
[]
[]
[Materials]
[dc0]
type = GenericConstantArray
block = 0
prop_name = dc
prop_value = '1 1'
[]
[dc1]
type = GenericConstantArray
block = 1
prop_name = dc
prop_value = '2 1'
[]
[rc]
type = GenericConstant2DArray
block = '0 1'
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Postprocessors]
[intu0]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 0
[]
[intu1]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/array_kernels/array_diffusion_reaction.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/array_kernels/array_diffusion_reaction_coupling.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
components = 2
[]
[v]
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[reaction]
type = ArrayReaction
variable = u
reaction_coefficient = rc
[]
[diffv]
type = Diffusion
variable = v
[]
[vu]
type = ArrayCoupledForce
variable = u
v = v
coef = '0 0.5'
[]
[]
[BCs]
[left]
type = ArrayDirichletBC
variable = u
boundary = 1
values = '0 0'
[]
[right]
type = ArrayDirichletBC
variable = u
boundary = 2
values = '1 2'
[]
[leftv]
type = DirichletBC
variable = v
boundary = 1
value = 0
[]
[rightv]
type = DirichletBC
variable = v
boundary = 2
value = 2
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[rc]
type = GenericConstant2DArray
prop_name = rc
prop_value = '1 0; -0.1 1'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Postprocessors]
[intu0]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 0
[]
[intu1]
type = ElementIntegralArrayVariablePostprocessor
variable = u
component = 1
[]
[intv]
type = ElementIntegralVariablePostprocessor
variable = v
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/hfem/array_dirichlet.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 3
ny = 3
dim = 2
[]
build_all_side_lowerd_mesh = true
[]
[Variables]
[u]
order = THIRD
family = MONOMIAL
block = 0
components = 2
[]
[lambda]
order = CONSTANT
family = MONOMIAL
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[lambdab]
order = CONSTANT
family = MONOMIAL
block = BOUNDARY_SIDE_LOWERD_SUBDOMAIN
components = 2
[]
[]
[AuxVariables]
[v]
order = CONSTANT
family = MONOMIAL
block = 0
initial_condition = '1'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
block = 0
diffusion_coefficient = dc
[]
[source]
type = ArrayCoupledForce
variable = u
v = v
coef = '1 2'
block = 0
[]
[]
[DGKernels]
[surface]
type = ArrayHFEMDiffusion
variable = u
lowerd_variable = lambda
[]
[]
[BCs]
[all]
type = ArrayHFEMDirichletBC
boundary = 'left right top bottom'
variable = u
lowerd_variable = lambdab
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Postprocessors]
[intu]
type = ElementIntegralArrayVariablePostprocessor
variable = u
block = 0
[]
[lambdanorm]
type = ElementArrayL2Norm
variable = lambda
block = INTERNAL_SIDE_LOWERD_SUBDOMAIN
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -snes_linesearch_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu basic mumps'
[]
[Outputs]
[out]
# we hide lambda because it may flip sign due to element
# renumbering with distributed mesh
type = Exodus
hide = lambda
[]
[]
(test/tests/multiapps/relaxation/picard_relaxed_array_sub.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
[]
[]
[Variables]
[v]
components = 2
[]
[]
[AuxVariables]
[u]
components = 2
[]
[]
[Kernels]
[diff_v]
type = ArrayDiffusion
variable = v
diffusion_coefficient = dc
[]
[force_v]
type = ArrayCoupledForce
variable = v
v = u
is_v_array = true
coef = '1 1'
[]
[time_v]
type = ArrayTimeDerivative
variable = v
time_derivative_coefficient = tc
[]
[]
[BCs]
[left_v]
type = ArrayDirichletBC
variable = v
boundary = left
values = '2 2'
[]
[right_v]
type = ArrayDirichletBC
variable = v
boundary = right
values = '1 1'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[tc]
type = GenericConstantArray
prop_name = tc
prop_value = '1 1'
[]
[]
[Executioner]
type = Transient
num_steps = 20
dt = 0.1
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
nl_abs_tol = 1e-10
[]
[Outputs]
exodus = true
[]
(test/tests/interfaces/coupleable/array_coupling_by_name.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
[u]
components = 2
[]
[]
[AuxVariables]
[v]
components = 4
initial_condition = '1 2 3 4'
[]
[]
[Kernels]
[diff]
type = ArrayDiffusion
variable = u
diffusion_coefficient = dc
[]
[force]
type = ArrayCoupledForceVar
variable = u
v_coef = 'v 1.5,2.5,3.5,4.5'
[]
[]
[BCs]
[vac]
type = ArrayVacuumBC
variable = u
boundary = '0 1 2 3'
[]
[]
[Materials]
[dc]
type = GenericConstantArray
prop_name = dc
prop_value = '1 1'
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]