https://mooseframework.inl.gov
ADVolumeJunctionBaseUserObject.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "MooseVariableScalar.h"
12 #include "FEProblemBase.h"
13 #include "metaphysicl/parallel_numberarray.h"
14 #include "metaphysicl/parallel_dualnumber.h"
15 #include "metaphysicl/parallel_semidynamicsparsenumberarray.h"
16 #include "libmesh/parallel_algebra.h"
17 
20 {
22 
23  params.addParam<bool>(
24  "use_scalar_variables", true, "True if the junction variables are scalar variables");
26  "junction_subdomain_id",
28  "Junction subdomain ID (required if 'use_scalar_variables' is 'false')");
29 
30  params.addRequiredParam<Real>("volume", "Volume of the junction");
31  params.addRequiredParam<std::vector<UserObjectName>>(
32  "numerical_flux_names",
33  "The names of the user objects that compute the numerical flux at each flow channel.");
34  params.addClassDescription("User object to compute fluxes and residuals for a volume junction");
35  return params;
36 }
37 
39  : ADFlowJunctionUserObject(params),
40  _use_scalar_variables(getParam<bool>("use_scalar_variables")),
41  _junction_subdomain_id(getParam<subdomain_id_type>("junction_subdomain_id")),
42  _volume(getParam<Real>("volume")),
43  _numerical_flux_names(getParam<std::vector<UserObjectName>>("numerical_flux_names"))
44 {
46  mooseError(name(),
47  ": The number of supplied numerical flux objects '",
48  _numerical_flux_names.size(),
49  "' does not match the number of connections '",
51  "'.");
52 
53  if (!_use_scalar_variables && !isParamSetByUser("junction_subdomain_id"))
54  mooseError("If 'use_scalar_variables' is set to false, 'junction_subdomain_id' is required.");
55 }
56 
57 void
59 {
62 
63  _scalar_dofs.resize(_n_scalar_eq);
66  _residual.resize(_n_scalar_eq);
67 }
68 
69 void
71 {
72  _flux.assign(_n_connections, std::vector<ADReal>(_n_flux_eq, 0.0));
73  for (auto & i : _residual)
74  {
75  i.value() = 0;
76  i.derivatives() = DNDerivativeType();
77  }
78  _connection_indices.clear();
79 
80  auto junction_vars = getJunctionVariables();
81 
82  // Cache the junction variable values and get their dof indices
84  {
85  for (unsigned int i = 0; i < _n_scalar_eq; i++)
86  {
88  _scalar_dofs[i] = 0;
89  }
90 
91  for (const Elem * elem : *_mesh.getActiveLocalElementRange())
92  {
93  if (elem->subdomain_id() == _junction_subdomain_id)
94  {
95  // Reinitialize the element
97  _fe_problem.prepare(elem, _tid);
99 
100  // Cache the junction variable values and Dof indices
101  for (unsigned int i = 0; i < _n_scalar_eq; i++)
102  {
104 
105  auto && dofs = junction_vars[i]->dofIndices();
106  mooseAssert(dofs.size() == 1,
107  "There should be exactly 1 coupled DoF index for the variable '" +
108  junction_vars[i]->name() + "'.");
109  _scalar_dofs[i] = dofs[0];
110  }
111  }
112  }
113 
115  comm().sum(_scalar_dofs);
116  }
117  else
118  {
119  for (unsigned int i = 0; i < _junction_var_values.size(); i++)
121 
122  for (unsigned int i = 0; i < _n_scalar_eq; i++)
123  {
124  auto && dofs = junction_vars[i]->dofIndices();
125  mooseAssert(dofs.size() == 1,
126  "There should be exactly 1 coupled DoF index for the variable '" +
127  junction_vars[i]->name() + "'.");
128  _scalar_dofs[i] = dofs[0];
129  }
130  }
131 }
132 
133 void
135 {
136  // Get the connection index
137  const unsigned int c = getBoundaryIDIndex();
138  _connection_indices.push_back(c);
139 
140  // Get flow channel Dofs and basic function values
141  _flow_channel_dofs[c].clear();
142  for (unsigned int j = 0; j < _n_flux_eq; j++)
143  {
145 
146  auto && dofs = var->dofIndices();
147  for (unsigned int k = 0; k < dofs.size(); k++)
148  _flow_channel_dofs[c].push_back(dofs[k]);
149  }
150 }
151 
152 void
154 {
156 
157  const unsigned int c = getBoundaryIDIndex();
159 }
160 
161 void
163 {
164  const auto & volume_junction_uo = static_cast<const ADVolumeJunctionBaseUserObject &>(uo);
165 
166  // Store the data computed/retrieved in the other threads
167  for (unsigned int i = 0; i < volume_junction_uo._connection_indices.size(); i++)
168  {
169  const unsigned int c = volume_junction_uo._connection_indices[i];
170 
171  _flux[c] = volume_junction_uo._flux[c];
172  _flow_channel_dofs[c] = volume_junction_uo._flow_channel_dofs[c];
173  }
174 
175  // Add the scalar residuals from the other threads
176  for (unsigned int i = 0; i < _n_scalar_eq; i++)
177  _residual[i] += volume_junction_uo._residual[i];
178 }
179 
180 const std::vector<ADReal> &
182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184 
185  return _residual;
186 }
187 
188 const std::vector<ADReal> &
189 ADVolumeJunctionBaseUserObject::getFlux(const unsigned int & connection_index) const
190 {
191  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
192 
193  checkValidConnectionIndex(connection_index);
194  return _flux[connection_index];
195 }
196 
197 std::vector<const MooseVariableBase *>
199 {
200  std::vector<const MooseVariableBase *> vars(_scalar_variable_names.size());
201  for (unsigned int i = 0; i < _scalar_variable_names.size(); i++)
203  return vars;
204 }
205 
206 const MooseVariableBase *
207 ADVolumeJunctionBaseUserObject::getJunctionVar(const std::string & var_name, unsigned int i) const
208 {
209  const MooseVariableBase * var;
211  var = getScalarVar(var_name, i);
212  else
213  var = getVar(var_name, i);
214  return var;
215 }
216 
217 const ADVariableValue &
219  unsigned int i) const
220 {
221  return _use_scalar_variables ? adCoupledScalarValue(var_name, i) : adCoupledValue(var_name, i);
222 }
std::vector< std::vector< ADReal > > _flux
Cached flux vector for each connection.
std::vector< std::vector< dof_id_type > > _flow_channel_dofs
Degrees of freedom for flow channel variables, for each connection.
libMesh::ConstElemRange * getActiveLocalElementRange()
const unsigned int invalid_uint
virtual void computeFluxesAndResiduals(const unsigned int &c)=0
Computes and stores the fluxes, the scalar residuals, and their Jacobians.
virtual void prepare(const Elem *elem, const THREAD_ID tid) override
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< MOOSE_AD_MAX_DOFS_PER_ELEM > > DNDerivativeType
unsigned int getBoundaryIDIndex()
Gets the index of the currently executing boundary within the vector of boundary IDs given to this Si...
std::vector< std::string > _scalar_variable_names
Vector of coupled variable names for each scalar variable.
char ** vars
const ADVariableValue & adCoupledValue(const std::string &var_name, unsigned int comp=0) const
const Parallel::Communicator & comm() const
unsigned int _n_scalar_eq
Number of scalar residual components.
unsigned int _n_flux_eq
Number of flow channel flux components.
MooseVariable * getVar(const std::string &var_name, unsigned int comp)
const unsigned int _n_connections
Number of connected flow channels.
ADVolumeJunctionBaseUserObject(const InputParameters &params)
Constructor.
virtual const std::string & name() const
const std::vector< UserObjectName > & _numerical_flux_names
Names of numerical flux user objects for each connected flow channel.
void addRequiredParam(const std::string &name, const std::string &doc_string)
std::vector< ADReal > _residual
Cached scalar residual vector.
const subdomain_id_type _junction_subdomain_id
Junction subdomain ID.
virtual std::vector< const MooseVariableBase * > getJunctionVariables() const
Gets the junction variables.
std::vector< const ADVariableValue * > _junction_var_values
std::vector< unsigned int > _connection_indices
Connection indices for this thread.
std::vector< std::string > _flow_variable_names
Vector of coupled variable names for each flow variable.
const std::vector< ADReal > & getResidual() const
Returns the residual vector for the scalar variables.
virtual void threadJoin(const UserObject &uo) override
const std::vector< dof_id_type > & dofIndices() const final
virtual void reinitElem(const Elem *elem, const THREAD_ID tid) override
const std::vector< ADReal > & getFlux(const unsigned int &connection_index) const override
Gets the flux vector for a connection.
virtual void storeConnectionData()
Stores data (connection index, face shape functions, DoFs associated with flow channel variables) rel...
const MooseVariableScalar * getScalarVar(const std::string &var_name, unsigned int comp) const
Base class for computing and caching flux and residual vectors for a volume junction.
Provides common interfaces for flow junction user objects.
virtual void setCurrentSubdomainID(const Elem *elem, const THREAD_ID tid) override
std::vector< dof_id_type > _scalar_dofs
Degrees of freedom for scalar variables.
const bool _use_scalar_variables
True if the junction variables are scalar variables.
bool isParamSetByUser(const std::string &nm) const
const ADVariableValue & adCoupledScalarValue(const std::string &var_name, unsigned int comp=0) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MooseMesh & _mesh
FEProblemBase & _fe_problem
const THREAD_ID _tid
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
void checkValidConnectionIndex(const unsigned int &connection_index) const
Checks that a connection index is valid.
static const std::string k
Definition: NS.h:130
static InputParameters validParams()
const ADVariableValue & coupledJunctionValue(const std::string &var_name, unsigned int i=0) const
Gets an AD junction variable value.
const MooseVariableBase * getJunctionVar(const std::string &var_name, unsigned int i=0) const
Gets a junction variable.