www.mooseframework.org
NonlocalIntegratedBC.C
Go to the documentation of this file.
1 /****************************************************************/
2 /* DO NOT MODIFY THIS HEADER */
3 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
4 /* */
5 /* (c) 2010 Battelle Energy Alliance, LLC */
6 /* ALL RIGHTS RESERVED */
7 /* */
8 /* Prepared by Battelle Energy Alliance, LLC */
9 /* Under Contract No. DE-AC07-05ID14517 */
10 /* With the U. S. Department of Energy */
11 /* */
12 /* See COPYRIGHT for full restrictions */
13 /****************************************************************/
14 
15 #include "NonlocalIntegratedBC.h"
16 #include "Assembly.h"
17 #include "MooseVariable.h"
18 #include "Problem.h"
19 #include "SubProblem.h"
20 #include "SystemBase.h"
21 #include "MooseMesh.h"
22 
23 #include "libmesh/threads.h"
24 #include "libmesh/quadrature.h"
25 
26 template <>
29 {
31  return params;
32 }
33 
35  : IntegratedBC(parameters)
36 {
37  _mesh.errorIfDistributedMesh("NonlocalIntegratedBC");
38  mooseWarning("NonlocalIntegratedBC is a computationally expensive experimental capability used "
39  "only for integral terms.");
40 }
41 
42 void
44 {
45  DenseMatrix<Number> & ke = _assembly.jacobianBlock(_var.number(), _var.number());
46  _local_ke.resize(ke.m(), ke.n());
47  _local_ke.zero();
48 
49  for (_j = 0; _j < _phi.size();
50  _j++) // looping order for _i & _j are reversed for performance improvement
51  {
53  for (_i = 0; _i < _test.size(); _i++)
54  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
56  }
57 
58  ke += _local_ke;
59 
61  {
62  unsigned int rows = ke.m();
63  DenseVector<Number> diag(rows);
64  for (unsigned int i = 0; i < rows; i++)
65  diag(i) = _local_ke(i, i);
66 
67  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
68  for (const auto & var : _diag_save_in)
69  var->sys().solution().add_vector(diag, var->dofIndices());
70  }
71 }
72 
73 void
75 {
76  if (jvar == _var.number())
78  else
79  {
80  MooseVariable & jv = _sys.getVariable(_tid, jvar);
81  DenseMatrix<Number> & ke = _assembly.jacobianBlock(_var.number(), jvar);
82 
83  for (_j = 0; _j < _phi.size();
84  _j++) // looping order for _i & _j are reversed for performance improvement
85  {
86  getUserObjectJacobian(jvar, jv.dofIndices()[_j]);
87  for (_i = 0; _i < _test.size(); _i++)
88  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
89  ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpOffDiagJacobian(jvar);
90  }
91  }
92 }
93 
94 void
96 {
97  DenseMatrix<Number> & keg = _assembly.jacobianBlockNonlocal(_var.number(), _var.number());
98  // compiling set of global IDs for the local DOFs on the element
99  std::set<dof_id_type> local_dofindices(_var.dofIndices().begin(), _var.dofIndices().end());
100  // storing the global IDs for all the DOFs of the variable
101  const std::vector<dof_id_type> & var_alldofindices = _var.allDofIndices();
102  unsigned int n_total_dofs = var_alldofindices.size();
103 
104  for (_k = 0; _k < n_total_dofs;
105  _k++) // looping order for _i & _k are reversed for performance improvement
106  {
107  // eliminating the local components
108  auto it = local_dofindices.find(var_alldofindices[_k]);
109  if (it == local_dofindices.end())
110  {
111  getUserObjectJacobian(_var.number(), var_alldofindices[_k]);
112  // skip global DOFs that do not contribute to the jacobian
113  if (!globalDoFEnabled(_var, var_alldofindices[_k]))
114  continue;
115 
116  for (_i = 0; _i < _test.size(); _i++)
117  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
118  keg(_i, _k) += _JxW[_qp] * _coord[_qp] * computeQpNonlocalJacobian(var_alldofindices[_k]);
119  }
120  }
121 }
122 
123 void
125 {
126  if (jvar == _var.number())
128  else
129  {
130  MooseVariable & jv = _sys.getVariable(_tid, jvar);
131  DenseMatrix<Number> & keg = _assembly.jacobianBlockNonlocal(_var.number(), jvar);
132  // compiling set of global IDs for the local DOFs on the element
133  std::set<dof_id_type> local_dofindices(jv.dofIndices().begin(), jv.dofIndices().end());
134  // storing the global IDs for all the DOFs of the variable
135  const std::vector<dof_id_type> & jv_alldofindices = jv.allDofIndices();
136  unsigned int n_total_dofs = jv_alldofindices.size();
137 
138  for (_k = 0; _k < n_total_dofs;
139  _k++) // looping order for _i & _k are reversed for performance improvement
140  {
141  // eliminating the local components
142  auto it = local_dofindices.find(jv_alldofindices[_k]);
143  if (it == local_dofindices.end())
144  {
145  getUserObjectJacobian(jvar, jv_alldofindices[_k]);
146  // skip global DOFs that do not contribute to the jacobian
147  if (!globalDoFEnabled(jv, jv_alldofindices[_k]))
148  continue;
149 
150  for (_i = 0; _i < _test.size(); _i++)
151  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
152  keg(_i, _k) += _JxW[_qp] * _coord[_qp] *
153  computeQpNonlocalOffDiagJacobian(jvar, jv_alldofindices[_k]);
154  }
155  }
156  }
157 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:105
const std::vector< dof_id_type > & allDofIndices() const
Get all global dofindices for the variable.
void mooseWarning(Args &&...args) const
Definition: MooseObject.h:89
InputParameters validParams< NonlocalIntegratedBC >()
virtual Real computeQpNonlocalOffDiagJacobian(unsigned int, dof_id_type)
virtual bool globalDoFEnabled(MooseVariable &, dof_id_type)
optimization option for executing nonlocal jacobian calculation only for nonzero elements ...
InputParameters validParams< IntegratedBC >()
Definition: IntegratedBC.C:28
Class for stuff related to variables.
Definition: MooseVariable.h:43
virtual void computeNonlocalJacobian()
computeNonlocalJacobian and computeNonlocalOffDiagJacobian methods are introduced for providing the j...
const MooseArray< Real > & _coord
coordinate transformation
Definition: IntegratedBC.h:91
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const MooseArray< Real > & _JxW
transformed Jacobian weights
Definition: IntegratedBC.h:89
virtual Real computeQpJacobian()
Definition: IntegratedBC.C:198
void errorIfDistributedMesh(std::string name) const
Generate a unified error message if the underlying libMesh mesh is a DistributedMesh.
Definition: MooseMesh.C:2473
virtual void getUserObjectJacobian(unsigned int, dof_id_type)
Optimization option for getting jocobinas from userobject once per dof.
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:98
NonlocalIntegratedBC(const InputParameters &parameters)
virtual MooseVariable & getVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a variable of with specified name.
Definition: SystemBase.C:103
MooseVariable & _var
variable this BC works on
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: IntegratedBC.h:128
Assembly & _assembly
Reference to assembly.
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:33
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
Definition: IntegratedBC.C:204
unsigned int number() const
Get variable number coming from libMesh.
std::vector< dof_id_type > & dofIndices()
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:250
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
Definition: IntegratedBC.h:120
unsigned int _qp
quadrature point index
Definition: IntegratedBC.h:83
unsigned int _j
Definition: IntegratedBC.h:93
SystemBase & _sys
Reference to SystemBase.
std::vector< MooseVariable * > _diag_save_in
Definition: IntegratedBC.h:129
virtual Real computeQpNonlocalJacobian(dof_id_type)
Compute this IntegratedBC&#39;s contribution to the Jacobian corresponding to nolocal dof at the current ...
virtual void computeNonlocalOffDiagJacobian(unsigned int jvar)
Computes d-residual / d-jvar...
virtual void computeJacobianBlock(unsigned int jvar)
Computes d-ivar-residual / d-jvar...
virtual void computeJacobian()
computeJacobian and computeQpOffDiagJacobian methods are almost same as IntegratedBC except for few a...
DenseMatrix< Number > & jacobianBlock(unsigned int ivar, unsigned int jvar)
Definition: Assembly.C:887
unsigned int _i
i-th, j-th index for enumerating test and shape functions
Definition: IntegratedBC.h:93
THREAD_ID _tid
Thread id.
DenseMatrix< Number > & jacobianBlockNonlocal(unsigned int ivar, unsigned int jvar)
Definition: Assembly.C:894
QBase *& _qrule
active quadrature rule
Definition: IntegratedBC.h:85
MooseMesh & _mesh
Mesh this BC is defined on.