https://mooseframework.inl.gov
MooseVariableDataLinearFV.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 "MooseVariableDataFV.h"
12 #include "MooseVariableField.h"
13 #include "Assembly.h"
14 #include "MooseError.h"
15 #include "DisplacedSystem.h"
16 #include "TimeIntegrator.h"
17 #include "MooseVariableFV.h"
18 #include "MooseTypes.h"
19 #include "MooseMesh.h"
20 #include "Attributes.h"
21 #include "FVDirichletBCBase.h"
22 #include "SubProblem.h"
23 #include "FVKernel.h"
24 #include "ADUtils.h"
25 #include "MooseLinearVariableFV.h"
26 
27 #include "libmesh/quadrature.h"
28 #include "libmesh/fe_base.h"
29 #include "libmesh/system.h"
30 #include "libmesh/type_n_tensor.h"
31 
32 template <typename OutputType>
35  SystemBase & sys,
36  THREAD_ID tid,
37  Moose::ElementType element_type,
38  const Elem * const & elem)
39  : MooseVariableDataBase<OutputType>(var, sys, tid),
40  _var(var),
41  _fe_type(_var.feType()),
42  _var_num(_var.number()),
43  _assembly(_subproblem.assembly(_tid, var.kind() == Moose::VAR_SOLVER ? sys.number() : 0)),
44  _element_type(element_type),
45  _time_integrator(_sys.queryTimeIntegrator(_var_num)),
46  _elem(elem),
47  _displaced(dynamic_cast<const DisplacedSystem *>(&_sys) ? true : false),
48  _qrule(nullptr)
49 {
50 }
51 
52 template <typename OutputType>
53 void
55 {
56  switch (gm_type)
57  {
58  case Moose::Volume:
59  {
60  _qrule = _assembly.qRule();
61  // TODO: set integration multiplier to cell volume
62  break;
63  }
64  case Moose::Face:
65  {
66  _qrule = _assembly.qRuleFace();
67  // TODO: set integration multiplier to face area
68  break;
69  }
70  }
71 }
72 
73 template <typename OutputType>
74 void
76 {
77  auto && active_coupleable_matrix_tags =
78  _sys.subproblem().getActiveFEVariableCoupleableMatrixTags(_tid);
79  mooseAssert(_qrule, "We should have a non-null qrule");
80  const auto nqp = _qrule->n_points();
81 
82  for (auto tag : _required_vector_tags)
83  if (_need_vector_tag_u[tag])
84  {
85  _vector_tag_u[tag].resize(nqp);
86  assignForAllQps(0, _vector_tag_u[tag], nqp);
87  }
88 
89  for (auto tag : active_coupleable_matrix_tags)
90  if (_need_matrix_tag_u[tag])
91  {
92  _matrix_tag_u[tag].resize(nqp);
93  assignForAllQps(0, _matrix_tag_u[tag], nqp);
94  }
95 }
96 
97 template <typename OutputType>
98 const std::vector<dof_id_type> &
100 {
101  Moose::initDofIndices(*this, *_elem);
102  return _dof_indices;
103 }
104 
105 template <typename OutputType>
106 void
108 {
109  initDofIndices();
110  initializeSolnVars();
111 
112  unsigned int num_dofs = _dof_indices.size();
113 
114  if (num_dofs > 0)
115  fetchDoFValues();
116  else
117  // We don't have any dofs. There's nothing to do
118  return;
119 
120  mooseAssert(num_dofs == 1 && _vector_tags_dof_u[_solution_tag].size() == 1,
121  "There should only be one dof per elem for FV variables");
122 
123  const auto nqp = _qrule->n_points();
124  auto && active_coupleable_matrix_tags =
125  _sys.subproblem().getActiveFEVariableCoupleableMatrixTags(_tid);
126 
127  for (const auto qp : make_range(nqp))
128  {
129  for (auto tag : _required_vector_tags)
130  if (_need_vector_tag_u[tag])
131  _vector_tag_u[tag][qp] = _vector_tags_dof_u[tag][0];
132 
133  for (auto tag : active_coupleable_matrix_tags)
134  if (_need_matrix_tag_u[tag])
135  _matrix_tag_u[tag][qp] = _matrix_tags_dof_u[tag][0];
136  }
137 }
138 
139 template <typename OutputType>
140 void
142 {
143  mooseAssert(index == 0, "We only ever have one dof value locally");
144  _vector_tags_dof_u[_solution_tag][index] = value;
145  _has_dof_values = true;
146 
147  auto & u = _vector_tag_u[_solution_tag];
148  // Update the qp values as well
149  for (const auto qp : index_range(u))
150  u[qp] = value;
151 }
152 
153 template <typename OutputType>
154 void
155 MooseVariableDataLinearFV<OutputType>::setDofValues(const DenseVector<OutputData> & values)
156 {
157  auto & dof_values = _vector_tags_dof_u[_solution_tag];
158  for (unsigned int i = 0; i < values.size(); i++)
159  dof_values[i] = values(i);
160  _has_dof_values = true;
161 }
162 
163 template <typename OutputType>
164 void
166  std::vector<dof_id_type> & dof_indices) const
167 {
168  _dof_map.dof_indices(elem, dof_indices, _var_num);
169 }
170 
171 template <typename OutputType>
174 {
175  return _var;
176 }
177 
178 template class MooseVariableDataLinearFV<Real>;
179 // TODO: implement vector fv variable support. This will require some template
180 // specializations for various member functions in this and the FV variable
181 // classes. And then you will need to uncomment out the line below:
182 // template class MooseVariableDataLinearFV<RealVectorValue>;
void initDofIndices(T &data, const Elem &elem)
Moose::DOFType< OutputType >::type OutputData
MooseVariableDataLinearFV(const MooseLinearVariableFV< OutputType > &var, SystemBase &sys, THREAD_ID tid, Moose::ElementType element_type, const Elem *const &elem)
ElementType
Definition: MooseTypes.h:763
Base class for a system (of equations)
Definition: SystemBase.h:84
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
void setGeometry(Moose::GeometryType gm_type)
Set the geometry type before calculating variables values.
void computeValues()
Compute the variable values.
void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const
Get the dof indices for an element.
Class holding the data members for linear finite volume variables.
GeometryType
Definition: MooseTypes.h:248
void setDofValue(const OutputData &value, unsigned int index)
Set local DOF value at index to value .
IntRange< T > make_range(T beg, T end)
virtual const MooseLinearVariableFV< OutputType > & var() const override
Get the corresponding variable.
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
const std::vector< dof_id_type > & initDofIndices()
Fetch and return the dof indices of this variable on the current element.
This class provides variable solution interface for linear finite volume problems.
Definition: FVUtils.h:24
auto index_range(const T &sizable)
void setDofValues(const DenseVector< OutputData > &values)
Set local DOF values to the entries of values .
unsigned int THREAD_ID
Definition: MooseTypes.h:209