www.mooseframework.org
AdvectiveFluxCalculatorBase.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "Conversion.h" // for stringify
12 
13 #include "libmesh/string_to_enum.h"
14 #include "libmesh/mesh_tools.h"
15 #include "libmesh/parallel_sync.h"
16 
17 template <>
18 InputParameters
20 {
21  InputParameters params = validParams<ElementUserObject>();
22  params.addClassDescription(
23  "Base class to compute K_ij (a measure of advective flux from node i to node j) "
24  "and R+ and R- (which quantify amount of antidiffusion to add) in the "
25  "Kuzmin-Turek FEM-TVD multidimensional scheme");
26  MooseEnum flux_limiter_type("MinMod VanLeer MC superbee None", "VanLeer");
27  params.addParam<MooseEnum>("flux_limiter_type",
28  flux_limiter_type,
29  "Type of flux limiter to use. 'None' means that no antidiffusion "
30  "will be added in the Kuzmin-Turek scheme");
31  params.addRangeCheckedParam<Real>(
32  "allowable_MB_wastage",
33  5.0,
34  "allowable_MB_wastage > 0.0",
35  "This object will issue a memory warning if the internal node-numbering data structure "
36  "wastes more than allowable_MB_wastage megabytes. This data structure uses sequential "
37  "node-numbering which is optimized for speed rather than memory efficiency");
38 
39  params.addRelationshipManager("ElementPointNeighborLayers",
40  Moose::RelationshipManagerType::GEOMETRIC |
41  Moose::RelationshipManagerType::ALGEBRAIC |
42  Moose::RelationshipManagerType::COUPLING,
43  [](const InputParameters &, InputParameters & rm_params) {
44  rm_params.set<unsigned short>("layers") = 2;
45  });
46 
47  params.set<ExecFlagEnum>("execute_on", true) = {EXEC_LINEAR};
48  return params;
49 }
50 
52  : ElementUserObject(parameters),
53  _resizing_needed(true),
54  _flux_limiter_type(getParam<MooseEnum>("flux_limiter_type").getEnum<FluxLimiterTypeEnum>()),
55  _kij(),
56  _flux_out(),
57  _dflux_out_du(),
58  _dflux_out_dKjk(),
59  _valence(),
60  _u_nodal(),
61  _u_nodal_computed_by_thread(),
62  _connections(),
63  _number_of_nodes(0),
64  _my_pid(processor_id()),
65  _nodes_to_receive(),
66  _nodes_to_send(),
67  _pairs_to_receive(),
68  _pairs_to_send(),
69  _allowable_MB_wastage(getParam<Real>("allowable_MB_wastage"))
70 {
71  if (!_execute_enum.contains(EXEC_LINEAR))
72  paramError(
73  "execute_on",
74  "The AdvectiveFluxCalculator UserObject " + name() +
75  " execute_on parameter must include, at least, 'linear'. This is to ensure that "
76  "this UserObject computes all necessary quantities just before the Kernels evaluate "
77  "their Residuals");
78 }
79 
80 void
82 {
83  // If needed, size and initialize quantities appropriately, and compute _valence
84  if (_resizing_needed)
85  {
86  /*
87  * Populate _connections for all nodes that can be seen by this processor and on relevant
88  * blocks
89  *
90  * MULTIPROC NOTE: this must loop over local elements and 2 layers of ghosted elements.
91  * The Kernel will only loop over local elements, so will only use _kij, etc, for
92  * linked node-node pairs that appear in the local elements. Nevertheless, we
93  * need to build _kij, etc, for the nodes in the ghosted elements in order to simplify
94  * Jacobian computations
95  */
97  for (const auto & elem : _fe_problem.getEvaluableElementRange())
98  if (this->hasBlocks(elem->subdomain_id()))
99  for (unsigned i = 0; i < elem->n_nodes(); ++i)
100  _connections.addGlobalNode(elem->node_id(i));
102  for (const auto & elem : _fe_problem.getEvaluableElementRange())
103  if (this->hasBlocks(elem->subdomain_id()))
104  for (unsigned i = 0; i < elem->n_nodes(); ++i)
105  for (unsigned j = 0; j < elem->n_nodes(); ++j)
106  _connections.addConnection(elem->node_id(i), elem->node_id(j));
108 
110 
111  Real mb_wasted = (_connections.sizeSequential() - _number_of_nodes) * 4.0 / 1048576;
112  gatherMax(mb_wasted);
113  if (mb_wasted > _allowable_MB_wastage)
114  mooseWarning(
115  "In at least one processor, the sequential node-numbering internal data structure used "
116  "in " +
117  name() + " is using memory inefficiently.\nThe memory wasted is " +
118  Moose::stringify(mb_wasted) +
119  " megabytes.\n The node-numbering data structure has been optimized for speed at the "
120  "expense of memory, but that may not be an appropriate optimization for your case, "
121  "because the node numbering is not close to sequential in your case.\n");
122 
123  // initialize _kij
124  _kij.resize(_number_of_nodes);
125  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
126  _kij[sequential_i].assign(
127  _connections.sequentialConnectionsToSequentialID(sequential_i).size(), 0.0);
128 
129  /*
130  * Build _valence[i], which is the number of times the sequential node i is encountered when
131  * looping over the entire mesh (and on relevant blocks)
132  *
133  * MULTIPROC NOTE: this must loop over local elements and >=1 layer of ghosted elements.
134  * The Kernel will only loop over local elements, so will only use _valence for
135  * linked node-node pairs that appear in the local elements. But other processors will
136  * loop over neighboring elements, so avoid multiple counting of the residual and Jacobian
137  * this processor must record how many times each node-node link of its local elements
138  * appears in the local elements and >=1 layer, and pass that info to the Kernel
139  */
140  _valence.assign(_number_of_nodes, 0);
141  for (const auto & elem : _fe_problem.getEvaluableElementRange())
142  if (this->hasBlocks(elem->subdomain_id()))
143  for (unsigned i = 0; i < elem->n_nodes(); ++i)
144  {
145  const dof_id_type node_i = elem->node_id(i);
146  const dof_id_type sequential_i = _connections.sequentialID(node_i);
147  _valence[sequential_i] += 1;
148  }
149 
150  _u_nodal.assign(_number_of_nodes, 0.0);
152  _flux_out.assign(_number_of_nodes, 0.0);
153  _dflux_out_du.assign(_number_of_nodes, std::map<dof_id_type, Real>());
154  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
155  zeroedConnection(_dflux_out_du[sequential_i], _connections.globalID(sequential_i));
157  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
158  {
159  const std::vector<dof_id_type> con_i =
161  const std::size_t num_con_i = con_i.size();
162  _dflux_out_dKjk[sequential_i].resize(num_con_i);
163  for (std::size_t j = 0; j < con_i.size(); ++j)
164  {
165  const dof_id_type sequential_j = con_i[j];
166  const std::size_t num_con_j =
168  _dflux_out_dKjk[sequential_i][j].assign(num_con_j, 0.0);
169  }
170  }
171 
172  if (_app.n_processors() > 1)
173  buildCommLists();
174 
175  _resizing_needed = false;
176 
177  // Clear and size all member vectors
178  _dij.assign(_number_of_nodes, {});
179  _dDij_dKij.assign(_number_of_nodes, {});
180  _dDij_dKji.assign(_number_of_nodes, {});
181  _dDii_dKij.assign(_number_of_nodes, {});
182  _dDii_dKji.assign(_number_of_nodes, {});
183  _lij.assign(_number_of_nodes, {});
184  _rP.assign(_number_of_nodes, 0.0);
185  _rM.assign(_number_of_nodes, 0.0);
186  _drP.assign(_number_of_nodes, {});
187  _drM.assign(_number_of_nodes, {});
188  _drP_dk.assign(_number_of_nodes, {});
189  _drM_dk.assign(_number_of_nodes, {});
190  _fa.assign(_number_of_nodes, {});
191  _dfa.assign(_number_of_nodes, {});
192  _dFij_dKik.assign(_number_of_nodes, {});
193  _dFij_dKjk.assign(_number_of_nodes, {});
194  }
195 }
196 
197 void
199 {
200  ElementUserObject::meshChanged();
201 
202  // Signal that _kij, _valence, etc need to be rebuilt
203  _resizing_needed = true;
204 }
205 
206 void
208 {
209  // Zero _kij and falsify _u_nodal_computed_by_thread ready for building in execute() and
210  // finalize()
212  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
213  _kij[sequential_i].assign(_connections.sequentialConnectionsToSequentialID(sequential_i).size(),
214  0.0);
215 }
216 
217 void
219 {
220  // compute _kij contributions from this element that is local to this processor
221  // and record _u_nodal
222  for (unsigned i = 0; i < _current_elem->n_nodes(); ++i)
223  {
224  const dof_id_type node_i = _current_elem->node_id(i);
225  const dof_id_type sequential_i = _connections.sequentialID(node_i);
226  if (!_u_nodal_computed_by_thread[sequential_i])
227  {
228  _u_nodal[sequential_i] = computeU(i);
229  _u_nodal_computed_by_thread[sequential_i] = true;
230  }
231  for (unsigned j = 0; j < _current_elem->n_nodes(); ++j)
232  {
233  const dof_id_type node_j = _current_elem->node_id(j);
234  for (unsigned qp = 0; qp < _qrule->n_points(); ++qp)
235  executeOnElement(node_i, node_j, i, j, qp);
236  }
237  }
238 }
239 
240 void
242  dof_id_type global_i, dof_id_type global_j, unsigned local_i, unsigned local_j, unsigned qp)
243 {
244  // KT Eqn (20)
245  const dof_id_type sequential_i = _connections.sequentialID(global_i);
246  const unsigned index_i_to_j = _connections.indexOfGlobalConnection(global_i, global_j);
247  _kij[sequential_i][index_i_to_j] += _JxW[qp] * _coord[qp] * computeVelocity(local_i, local_j, qp);
248 }
249 
250 void
252 {
253  const AdvectiveFluxCalculatorBase & afc = static_cast<const AdvectiveFluxCalculatorBase &>(uo);
254  // add the values of _kij computed by different threads
255  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
256  {
257  const std::size_t num_con_i =
259  for (std::size_t j = 0; j < num_con_i; ++j)
260  _kij[sequential_i][j] += afc._kij[sequential_i][j];
261  }
262 
263  // gather the values of _u_nodal computed by different threads
264  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
265  if (!_u_nodal_computed_by_thread[sequential_i] && afc._u_nodal_computed_by_thread[sequential_i])
266  _u_nodal[sequential_i] = afc._u_nodal[sequential_i];
267 }
268 
269 void
271 {
272  // Overall: ensure _kij is fully built, then compute Kuzmin-Turek D, L, f^a and
273  // relevant Jacobian information, and then the relevant quantities into _flux_out and
274  // _dflux_out_du, _dflux_out_dKjk
275 
276  if (_app.n_processors() > 1)
278 
279  // Calculate KuzminTurek D matrix
280  // See Eqn (32)
281  // This adds artificial diffusion, which eliminates any spurious oscillations
282  // The idea is that D will remove all negative off-diagonal elements when it is added to K
283  // This is identical to full upwinding
284  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
285  {
286  const std::vector<dof_id_type> & con_i =
288  const std::size_t num_con_i = con_i.size();
289  _dij[sequential_i].assign(num_con_i, 0.0);
290  _dDij_dKij[sequential_i].assign(num_con_i, 0.0);
291  _dDij_dKji[sequential_i].assign(num_con_i, 0.0);
292  _dDii_dKij[sequential_i].assign(num_con_i, 0.0);
293  _dDii_dKji[sequential_i].assign(num_con_i, 0.0);
294  const unsigned index_i_to_i =
295  _connections.indexOfSequentialConnection(sequential_i, sequential_i);
296  for (std::size_t j = 0; j < num_con_i; ++j)
297  {
298  const dof_id_type sequential_j = con_i[j];
299  if (sequential_i == sequential_j)
300  continue;
301  const unsigned index_j_to_i =
302  _connections.indexOfSequentialConnection(sequential_j, sequential_i);
303  const Real kij = _kij[sequential_i][j];
304  const Real kji = _kij[sequential_j][index_j_to_i];
305  if ((kij <= kji) && (kij < 0.0))
306  {
307  _dij[sequential_i][j] = -kij;
308  _dDij_dKij[sequential_i][j] = -1.0;
309  _dDii_dKij[sequential_i][j] += 1.0;
310  }
311  else if ((kji <= kij) && (kji < 0.0))
312  {
313  _dij[sequential_i][j] = -kji;
314  _dDij_dKji[sequential_i][j] = -1.0;
315  _dDii_dKji[sequential_i][j] += 1.0;
316  }
317  else
318  _dij[sequential_i][j] = 0.0;
319  _dij[sequential_i][index_i_to_i] -= _dij[sequential_i][j];
320  }
321  }
322 
323  // Calculate KuzminTurek L matrix
324  // See Fig 2: L = K + D
325  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
326  {
327  const std::vector<dof_id_type> & con_i =
329  const std::size_t num_con_i = con_i.size();
330  _lij[sequential_i].assign(num_con_i, 0.0);
331  for (std::size_t j = 0; j < num_con_i; ++j)
332  _lij[sequential_i][j] = _kij[sequential_i][j] + _dij[sequential_i][j];
333  }
334 
335  // Compute KuzminTurek R matrices
336  // See Eqns (49) and (12)
337  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
338  {
339  _rP[sequential_i] = rPlus(sequential_i, _drP[sequential_i], _drP_dk[sequential_i]);
340  _rM[sequential_i] = rMinus(sequential_i, _drM[sequential_i], _drM_dk[sequential_i]);
341  }
342 
343  // Calculate KuzminTurek f^{a} matrix
344  // This is the antidiffusive flux
345  // See Eqn (50)
346  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
347  {
348  const std::vector<dof_id_type> & con_i =
350  const unsigned num_con_i = con_i.size();
351  _fa[sequential_i].assign(num_con_i, 0.0);
352  _dfa[sequential_i].resize(num_con_i);
353  _dFij_dKik[sequential_i].resize(num_con_i);
354  _dFij_dKjk[sequential_i].resize(num_con_i);
355  for (std::size_t j = 0; j < num_con_i; ++j)
356  {
357  for (const auto & global_k : con_i)
358  _dfa[sequential_i][j][global_k] = 0;
359  const dof_id_type global_j = con_i[j];
360  const std::vector<dof_id_type> & con_j = _connections.globalConnectionsToGlobalID(global_j);
361  const unsigned num_con_j = con_j.size();
362  for (const auto & global_k : con_j)
363  _dfa[sequential_i][j][global_k] = 0;
364  _dFij_dKik[sequential_i][j].assign(num_con_i, 0.0);
365  _dFij_dKjk[sequential_i][j].assign(num_con_j, 0.0);
366  }
367  }
368 
369  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
370  {
371  const dof_id_type global_i = _connections.globalID(sequential_i);
372  const Real u_i = _u_nodal[sequential_i];
373  const std::vector<dof_id_type> & con_i =
375  const std::size_t num_con_i = con_i.size();
376  for (std::size_t j = 0; j < num_con_i; ++j)
377  {
378  const dof_id_type sequential_j = con_i[j];
379  if (sequential_i == sequential_j)
380  continue;
381  const unsigned index_j_to_i =
382  _connections.indexOfSequentialConnection(sequential_j, sequential_i);
383  const Real Lij = _lij[sequential_i][j];
384  const Real Lji = _lij[sequential_j][index_j_to_i];
385  if (Lji >= Lij) // node i is upwind of node j.
386  {
387  const Real Dij = _dij[sequential_i][j];
388  const dof_id_type global_j = _connections.globalID(sequential_j);
389  const Real u_j = _u_nodal[sequential_j];
390  Real prefactor = 0.0;
391  std::vector<Real> dprefactor_du(num_con_i,
392  0.0); // dprefactor_du[j] = d(prefactor)/du[sequential_j]
393  std::vector<Real> dprefactor_dKij(
394  num_con_i, 0.0); // dprefactor_dKij[j] = d(prefactor)/dKij[sequential_i][j].
395  Real dprefactor_dKji = 0; // dprefactor_dKji = d(prefactor)/dKij[sequential_j][index_j_to_i]
396  if (u_i >= u_j)
397  {
398  if (Lji <= _rP[sequential_i] * Dij)
399  {
400  prefactor = Lji;
401  dprefactor_dKij[j] += _dDij_dKji[sequential_j][index_j_to_i];
402  dprefactor_dKji += 1.0 + _dDij_dKij[sequential_j][index_j_to_i];
403  }
404  else
405  {
406  prefactor = _rP[sequential_i] * Dij;
407  for (std::size_t ind_j = 0; ind_j < num_con_i; ++ind_j)
408  dprefactor_du[ind_j] = _drP[sequential_i][ind_j] * Dij;
409  dprefactor_dKij[j] += _rP[sequential_i] * _dDij_dKij[sequential_i][j];
410  dprefactor_dKji += _rP[sequential_i] * _dDij_dKji[sequential_i][j];
411  for (std::size_t ind_j = 0; ind_j < num_con_i; ++ind_j)
412  dprefactor_dKij[ind_j] += _drP_dk[sequential_i][ind_j] * Dij;
413  }
414  }
415  else
416  {
417  if (Lji <= _rM[sequential_i] * Dij)
418  {
419  prefactor = Lji;
420  dprefactor_dKij[j] += _dDij_dKji[sequential_j][index_j_to_i];
421  dprefactor_dKji += 1.0 + _dDij_dKij[sequential_j][index_j_to_i];
422  }
423  else
424  {
425  prefactor = _rM[sequential_i] * Dij;
426  for (std::size_t ind_j = 0; ind_j < num_con_i; ++ind_j)
427  dprefactor_du[ind_j] = _drM[sequential_i][ind_j] * Dij;
428  dprefactor_dKij[j] += _rM[sequential_i] * _dDij_dKij[sequential_i][j];
429  dprefactor_dKji += _rM[sequential_i] * _dDij_dKji[sequential_i][j];
430  for (std::size_t ind_j = 0; ind_j < num_con_i; ++ind_j)
431  dprefactor_dKij[ind_j] += _drM_dk[sequential_i][ind_j] * Dij;
432  }
433  }
434  _fa[sequential_i][j] = prefactor * (u_i - u_j);
435  _dfa[sequential_i][j][global_i] = prefactor;
436  _dfa[sequential_i][j][global_j] = -prefactor;
437  for (std::size_t ind_j = 0; ind_j < num_con_i; ++ind_j)
438  {
439  _dfa[sequential_i][j][_connections.globalID(con_i[ind_j])] +=
440  dprefactor_du[ind_j] * (u_i - u_j);
441  _dFij_dKik[sequential_i][j][ind_j] += dprefactor_dKij[ind_j] * (u_i - u_j);
442  }
443  _dFij_dKjk[sequential_i][j][index_j_to_i] += dprefactor_dKji * (u_i - u_j);
444  }
445  }
446  }
447 
448  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
449  {
450  const std::vector<dof_id_type> & con_i =
452  const std::size_t num_con_i = con_i.size();
453  for (std::size_t j = 0; j < num_con_i; ++j)
454  {
455  const dof_id_type sequential_j = con_i[j];
456  if (sequential_i == sequential_j)
457  continue;
458  const unsigned index_j_to_i =
459  _connections.indexOfSequentialConnection(sequential_j, sequential_i);
460  if (_lij[sequential_j][index_j_to_i] < _lij[sequential_i][j]) // node_i is downwind of node_j.
461  {
462  _fa[sequential_i][j] = -_fa[sequential_j][index_j_to_i];
463  for (const auto & dof_deriv : _dfa[sequential_j][index_j_to_i])
464  _dfa[sequential_i][j][dof_deriv.first] = -dof_deriv.second;
465  for (std::size_t k = 0; k < num_con_i; ++k)
466  _dFij_dKik[sequential_i][j][k] = -_dFij_dKjk[sequential_j][index_j_to_i][k];
467  const std::size_t num_con_j =
469  for (std::size_t k = 0; k < num_con_j; ++k)
470  _dFij_dKjk[sequential_i][j][k] = -_dFij_dKik[sequential_j][index_j_to_i][k];
471  }
472  }
473  }
474 
475  // zero _flux_out and its derivatives
476  _flux_out.assign(_number_of_nodes, 0.0);
477  // The derivatives are a bit complicated.
478  // If i is upwind of a node "j" then _flux_out[i] depends on all nodes connected to i.
479  // But if i is downwind of a node "j" then _flux_out depends on all nodes connected with node
480  // j.
481  _dflux_out_du.assign(_number_of_nodes, std::map<dof_id_type, Real>());
482  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
483  for (const auto & j : _connections.globalConnectionsToSequentialID(sequential_i))
484  {
485  _dflux_out_du[sequential_i][j] = 0.0;
486  for (const auto & neighbors_j : _connections.globalConnectionsToGlobalID(j))
487  _dflux_out_du[sequential_i][neighbors_j] = 0.0;
488  }
490  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
491  {
492  const std::vector<dof_id_type> & con_i =
494  const std::size_t num_con_i = con_i.size();
495  _dflux_out_dKjk[sequential_i].resize(num_con_i);
496  for (std::size_t j = 0; j < num_con_i; ++j)
497  {
498  const dof_id_type sequential_j = con_i[j];
499  const std::vector<dof_id_type> & con_j =
501  _dflux_out_dKjk[sequential_i][j].assign(con_j.size(), 0.0);
502  }
503  }
504 
505  // Add everything together
506  // See step 3 in Fig 2, noting Eqn (36)
507  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
508  {
509  const std::vector<dof_id_type> & con_i =
511  const size_t num_con_i = con_i.size();
512  const dof_id_type index_i_to_i =
513  _connections.indexOfSequentialConnection(sequential_i, sequential_i);
514  for (std::size_t j = 0; j < num_con_i; ++j)
515  {
516  const dof_id_type sequential_j = con_i[j];
517  const dof_id_type global_j = _connections.globalID(sequential_j);
518  const Real u_j = _u_nodal[sequential_j];
519 
520  // negative sign because residual = -Lu (KT equation (19))
521  _flux_out[sequential_i] -= _lij[sequential_i][j] * u_j + _fa[sequential_i][j];
522 
523  _dflux_out_du[sequential_i][global_j] -= _lij[sequential_i][j];
524  for (const auto & dof_deriv : _dfa[sequential_i][j])
525  _dflux_out_du[sequential_i][dof_deriv.first] -= dof_deriv.second;
526 
527  _dflux_out_dKjk[sequential_i][index_i_to_i][j] -= 1.0 * u_j; // from the K in L = K + D
528 
529  if (sequential_j == sequential_i)
530  for (dof_id_type k = 0; k < num_con_i; ++k)
531  _dflux_out_dKjk[sequential_i][index_i_to_i][k] -= _dDii_dKij[sequential_i][k] * u_j;
532  else
533  _dflux_out_dKjk[sequential_i][index_i_to_i][j] -= _dDij_dKij[sequential_i][j] * u_j;
534  for (dof_id_type k = 0; k < num_con_i; ++k)
535  _dflux_out_dKjk[sequential_i][index_i_to_i][k] -= _dFij_dKik[sequential_i][j][k];
536 
537  if (sequential_j == sequential_i)
538  for (unsigned k = 0; k < con_i.size(); ++k)
539  {
540  const unsigned index_k_to_i =
541  _connections.indexOfSequentialConnection(con_i[k], sequential_i);
542  _dflux_out_dKjk[sequential_i][k][index_k_to_i] -= _dDii_dKji[sequential_i][k] * u_j;
543  }
544  else
545  {
546  const unsigned index_j_to_i =
547  _connections.indexOfSequentialConnection(sequential_j, sequential_i);
548  _dflux_out_dKjk[sequential_i][j][index_j_to_i] -= _dDij_dKji[sequential_i][j] * u_j;
549  }
550  for (unsigned k = 0;
551  k < _connections.sequentialConnectionsToSequentialID(sequential_j).size();
552  ++k)
553  _dflux_out_dKjk[sequential_i][j][k] -= _dFij_dKjk[sequential_i][j][k];
554  }
555  }
556 }
557 
558 Real
559 AdvectiveFluxCalculatorBase::rPlus(dof_id_type sequential_i,
560  std::vector<Real> & dlimited_du,
561  std::vector<Real> & dlimited_dk) const
562 {
563  const std::size_t num_con = _connections.sequentialConnectionsToSequentialID(sequential_i).size();
564  dlimited_du.assign(num_con, 0.0);
565  dlimited_dk.assign(num_con, 0.0);
567  return 0.0;
568  std::vector<Real> dp_du;
569  std::vector<Real> dp_dk;
570  const Real p = PQPlusMinus(sequential_i, PQPlusMinusEnum::PPlus, dp_du, dp_dk);
571  if (p == 0.0)
572  // Comment after Eqn (49): if P=0 then there's no antidiffusion, so no need to remove it
573  return 1.0;
574  std::vector<Real> dq_du;
575  std::vector<Real> dq_dk;
576  const Real q = PQPlusMinus(sequential_i, PQPlusMinusEnum::QPlus, dq_du, dq_dk);
577 
578  const Real r = q / p;
579  Real limited;
580  Real dlimited_dr;
581  limitFlux(1.0, r, limited, dlimited_dr);
582 
583  const Real p2 = std::pow(p, 2);
584  for (std::size_t j = 0; j < num_con; ++j)
585  {
586  const Real dr_du = dq_du[j] / p - q * dp_du[j] / p2;
587  const Real dr_dk = dq_dk[j] / p - q * dp_dk[j] / p2;
588  dlimited_du[j] = dlimited_dr * dr_du;
589  dlimited_dk[j] = dlimited_dr * dr_dk;
590  }
591  return limited;
592 }
593 
594 Real
595 AdvectiveFluxCalculatorBase::rMinus(dof_id_type sequential_i,
596  std::vector<Real> & dlimited_du,
597  std::vector<Real> & dlimited_dk) const
598 {
599  const std::size_t num_con = _connections.sequentialConnectionsToSequentialID(sequential_i).size();
600  dlimited_du.assign(num_con, 0.0);
601  dlimited_dk.assign(num_con, 0.0);
603  return 0.0;
604  std::vector<Real> dp_du;
605  std::vector<Real> dp_dk;
606  const Real p = PQPlusMinus(sequential_i, PQPlusMinusEnum::PMinus, dp_du, dp_dk);
607  if (p == 0.0)
608  // Comment after Eqn (49): if P=0 then there's no antidiffusion, so no need to remove it
609  return 1.0;
610  std::vector<Real> dq_du;
611  std::vector<Real> dq_dk;
612  const Real q = PQPlusMinus(sequential_i, PQPlusMinusEnum::QMinus, dq_du, dq_dk);
613 
614  const Real r = q / p;
615  Real limited;
616  Real dlimited_dr;
617  limitFlux(1.0, r, limited, dlimited_dr);
618 
619  const Real p2 = std::pow(p, 2);
620  for (std::size_t j = 0; j < num_con; ++j)
621  {
622  const Real dr_du = dq_du[j] / p - q * dp_du[j] / p2;
623  const Real dr_dk = dq_dk[j] / p - q * dp_dk[j] / p2;
624  dlimited_du[j] = dlimited_dr * dr_du;
625  dlimited_dk[j] = dlimited_dr * dr_dk;
626  }
627  return limited;
628 }
629 
630 void
631 AdvectiveFluxCalculatorBase::limitFlux(Real a, Real b, Real & limited, Real & dlimited_db) const
632 {
633  limited = 0.0;
634  dlimited_db = 0.0;
636  return;
637 
638  if ((a >= 0.0 && b <= 0.0) || (a <= 0.0 && b >= 0.0))
639  return;
640  const Real s = (a > 0.0 ? 1.0 : -1.0);
641 
642  const Real lal = std::abs(a);
643  const Real lbl = std::abs(b);
644  const Real dlbl = (b >= 0.0 ? 1.0 : -1.0); // d(lbl)/db
645  switch (_flux_limiter_type)
646  {
648  {
649  if (lal <= lbl)
650  {
651  limited = s * lal;
652  dlimited_db = 0.0;
653  }
654  else
655  {
656  limited = s * lbl;
657  dlimited_db = s * dlbl;
658  }
659  return;
660  }
662  {
663  limited = s * 2 * lal * lbl / (lal + lbl);
664  dlimited_db = s * 2 * lal * (dlbl / (lal + lbl) - lbl * dlbl / std::pow(lal + lbl, 2));
665  return;
666  }
668  {
669  const Real av = 0.5 * std::abs(a + b);
670  if (2 * lal <= av && lal <= lbl)
671  {
672  // 2 * lal is the smallest
673  limited = s * 2.0 * lal;
674  dlimited_db = 0.0;
675  }
676  else if (2 * lbl <= av && lbl <= lal)
677  {
678  // 2 * lbl is the smallest
679  limited = s * 2.0 * lbl;
680  dlimited_db = s * 2.0 * dlbl;
681  }
682  else
683  {
684  // av is the smallest
685  limited = s * av;
686  // if (a>0 and b>0) then d(av)/db = 0.5 = 0.5 * dlbl
687  // if (a<0 and b<0) then d(av)/db = -0.5 = 0.5 * dlbl
688  // if a and b have different sign then limited=0, above
689  dlimited_db = s * 0.5 * dlbl;
690  }
691  return;
692  }
694  {
695  const Real term1 = std::min(2.0 * lal, lbl);
696  const Real term2 = std::min(lal, 2.0 * lbl);
697  if (term1 >= term2)
698  {
699  if (2.0 * lal <= lbl)
700  {
701  limited = s * 2 * lal;
702  dlimited_db = 0.0;
703  }
704  else
705  {
706  limited = s * lbl;
707  dlimited_db = s * dlbl;
708  }
709  }
710  else
711  {
712  if (lal <= 2.0 * lbl)
713  {
714  limited = s * lal;
715  dlimited_db = 0.0;
716  }
717  else
718  {
719  limited = s * 2.0 * lbl;
720  dlimited_db = s * 2.0 * dlbl;
721  }
722  }
723  return;
724  }
725  default:
726  return;
727  }
728 }
729 
730 const std::map<dof_id_type, Real> &
732 {
733  return _dflux_out_du[_connections.sequentialID(node_i)];
734 }
735 
736 const std::vector<std::vector<Real>> &
738 {
739  return _dflux_out_dKjk[_connections.sequentialID(node_i)];
740 }
741 
742 Real
744 {
745  return _flux_out[_connections.sequentialID(node_i)];
746 }
747 
748 unsigned
750 {
751  return _valence[_connections.sequentialID(node_i)];
752 }
753 
754 void
755 AdvectiveFluxCalculatorBase::zeroedConnection(std::map<dof_id_type, Real> & the_map,
756  dof_id_type node_i) const
757 {
758  the_map.clear();
759  for (const auto & node_j : _connections.globalConnectionsToGlobalID(node_i))
760  the_map[node_j] = 0.0;
761 }
762 
763 Real
765  const PQPlusMinusEnum pq_plus_minus,
766  std::vector<Real> & derivs,
767  std::vector<Real> & dpqdk) const
768 {
769  // Find the value of u at sequential_i
770  const Real u_i = _u_nodal[sequential_i];
771 
772  // Connections to sequential_i
773  const std::vector<dof_id_type> con_i =
775  const std::size_t num_con = con_i.size();
776  // The neighbor number of sequential_i to sequential_i
777  const unsigned i_index_i = _connections.indexOfSequentialConnection(sequential_i, sequential_i);
778 
779  // Initialize the results
780  Real result = 0.0;
781  derivs.assign(num_con, 0.0);
782  dpqdk.assign(num_con, 0.0);
783 
784  // Sum over all nodes connected with node_i.
785  for (std::size_t j = 0; j < num_con; ++j)
786  {
787  const dof_id_type sequential_j = con_i[j];
788  if (sequential_j == sequential_i)
789  continue;
790  const Real kentry = _kij[sequential_i][j];
791 
792  // Find the value of u at node_j
793  const Real u_j = _u_nodal[sequential_j];
794  const Real ujminusi = u_j - u_i;
795 
796  // Evaluate the i-j contribution to the result
797  switch (pq_plus_minus)
798  {
800  {
801  if (ujminusi < 0.0 && kentry < 0.0)
802  {
803  result += kentry * ujminusi;
804  derivs[j] += kentry;
805  derivs[i_index_i] -= kentry;
806  dpqdk[j] += ujminusi;
807  }
808  break;
809  }
811  {
812  if (ujminusi > 0.0 && kentry < 0.0)
813  {
814  result += kentry * ujminusi;
815  derivs[j] += kentry;
816  derivs[i_index_i] -= kentry;
817  dpqdk[j] += ujminusi;
818  }
819  break;
820  }
822  {
823  if (ujminusi > 0.0 && kentry > 0.0)
824  {
825  result += kentry * ujminusi;
826  derivs[j] += kentry;
827  derivs[i_index_i] -= kentry;
828  dpqdk[j] += ujminusi;
829  }
830  break;
831  }
833  {
834  if (ujminusi < 0.0 && kentry > 0.0)
835  {
836  result += kentry * ujminusi;
837  derivs[j] += kentry;
838  derivs[i_index_i] -= kentry;
839  dpqdk[j] += ujminusi;
840  }
841  break;
842  }
843  }
844  }
845 
846  return result;
847 }
848 
849 void
851 {
865  _nodes_to_receive.clear();
866  for (const auto & elem : _fe_problem.getEvaluableElementRange())
867  if (this->hasBlocks(elem->subdomain_id()))
868  {
869  const processor_id_type elem_pid = elem->processor_id();
870  if (elem_pid != _my_pid)
871  {
872  if (_nodes_to_receive.find(elem_pid) == _nodes_to_receive.end())
873  _nodes_to_receive[elem_pid] = std::vector<dof_id_type>();
874  for (unsigned i = 0; i < elem->n_nodes(); ++i)
875  if (std::find(_nodes_to_receive[elem_pid].begin(),
876  _nodes_to_receive[elem_pid].end(),
877  elem->node_id(i)) == _nodes_to_receive[elem_pid].end())
878  _nodes_to_receive[elem_pid].push_back(elem->node_id(i));
879  }
880  }
881 
882  // exchange this info with other processors, building _nodes_to_send at the same time
883  _nodes_to_send.clear();
884  auto nodes_action_functor = [this](processor_id_type pid, const std::vector<dof_id_type> & nts) {
885  _nodes_to_send[pid] = nts;
886  };
887  Parallel::push_parallel_vector_data(this->comm(), _nodes_to_receive, nodes_action_functor);
888 
889  // At the moment, _nodes_to_send and _nodes_to_receive contain global node numbers
890  // It is slightly more efficient to convert to sequential node numbers
891  // so that we don't have to keep doing things like
892  // _u_nodal[_connections.sequentialID(_nodes_to_send[pid][i])
893  // every time we send/receive
894  for (auto & kv : _nodes_to_send)
895  {
896  const processor_id_type pid = kv.first;
897  const std::size_t num_nodes = kv.second.size();
898  for (unsigned i = 0; i < num_nodes; ++i)
900  }
901  for (auto & kv : _nodes_to_receive)
902  {
903  const processor_id_type pid = kv.first;
904  const std::size_t num_nodes = kv.second.size();
905  for (unsigned i = 0; i < num_nodes; ++i)
907  }
908 
909  // Build pairs_to_receive
910  _pairs_to_receive.clear();
911  for (const auto & elem : _fe_problem.getEvaluableElementRange())
912  if (this->hasBlocks(elem->subdomain_id()))
913  {
914  const processor_id_type elem_pid = elem->processor_id();
915  if (elem_pid != _my_pid)
916  {
917  if (_pairs_to_receive.find(elem_pid) == _pairs_to_receive.end())
918  _pairs_to_receive[elem_pid] = std::vector<std::pair<dof_id_type, dof_id_type>>();
919  for (unsigned i = 0; i < elem->n_nodes(); ++i)
920  for (unsigned j = 0; j < elem->n_nodes(); ++j)
921  {
922  std::pair<dof_id_type, dof_id_type> the_pair(elem->node_id(i), elem->node_id(j));
923  if (std::find(_pairs_to_receive[elem_pid].begin(),
924  _pairs_to_receive[elem_pid].end(),
925  the_pair) == _pairs_to_receive[elem_pid].end())
926  _pairs_to_receive[elem_pid].push_back(the_pair);
927  }
928  }
929  }
930 
931  _pairs_to_send.clear();
932  auto pairs_action_functor = [this](processor_id_type pid,
933  const std::vector<std::pair<dof_id_type, dof_id_type>> & pts) {
934  _pairs_to_send[pid] = pts;
935  };
936  Parallel::push_parallel_vector_data(this->comm(), _pairs_to_receive, pairs_action_functor);
937 
938  // _pairs_to_send and _pairs_to_receive have been built using global node IDs
939  // since all processors know about that. However, using global IDs means
940  // that every time we send/receive, we keep having to do things like
941  // _kij[_connections.sequentialID(_pairs_to_send[pid][i].first)][_connections.indexOfGlobalConnection(_pairs_to_send[pid][i].first,
942  // _pairs_to_send[pid][i].second)] which is quite inefficient. So:
943  for (auto & kv : _pairs_to_send)
944  {
945  const processor_id_type pid = kv.first;
946  const std::size_t num_pairs = kv.second.size();
947  for (unsigned i = 0; i < num_pairs; ++i)
948  {
950  _pairs_to_send[pid][i].first, _pairs_to_send[pid][i].second);
951  _pairs_to_send[pid][i].first = _connections.sequentialID(_pairs_to_send[pid][i].first);
952  }
953  }
954  for (auto & kv : _pairs_to_receive)
955  {
956  const processor_id_type pid = kv.first;
957  const std::size_t num_pairs = kv.second.size();
958  for (unsigned i = 0; i < num_pairs; ++i)
959  {
961  _pairs_to_receive[pid][i].first, _pairs_to_receive[pid][i].second);
962  _pairs_to_receive[pid][i].first = _connections.sequentialID(_pairs_to_receive[pid][i].first);
963  }
964  }
965 }
966 
967 void
969 {
970  // Exchange ghosted _u_nodal information with other processors
971  std::map<processor_id_type, std::vector<Real>> unodal_to_send;
972  for (const auto & kv : _nodes_to_send)
973  {
974  const processor_id_type pid = kv.first;
975  unodal_to_send[pid] = std::vector<Real>();
976  for (const auto & nd : kv.second)
977  unodal_to_send[pid].push_back(_u_nodal[nd]);
978  }
979 
980  auto unodal_action_functor = [this](processor_id_type pid,
981  const std::vector<Real> & unodal_received) {
982  const std::size_t msg_size = unodal_received.size();
983  mooseAssert(msg_size == _nodes_to_receive[pid].size(),
984  "Message size, " << msg_size
985  << ", incompatible with nodes_to_receive, which has size "
986  << _nodes_to_receive[pid].size());
987  for (unsigned i = 0; i < msg_size; ++i)
988  _u_nodal[_nodes_to_receive[pid][i]] = unodal_received[i];
989  };
990  Parallel::push_parallel_vector_data(this->comm(), unodal_to_send, unodal_action_functor);
991 
992  // Exchange _kij information with other processors
993  std::map<processor_id_type, std::vector<Real>> kij_to_send;
994  for (const auto & kv : _pairs_to_send)
995  {
996  const processor_id_type pid = kv.first;
997  kij_to_send[pid] = std::vector<Real>();
998  for (const auto & pr : kv.second)
999  kij_to_send[pid].push_back(_kij[pr.first][pr.second]);
1000  }
1001 
1002  auto kij_action_functor = [this](processor_id_type pid, const std::vector<Real> & kij_received) {
1003  const std::size_t msg_size = kij_received.size();
1004  mooseAssert(msg_size == _pairs_to_receive[pid].size(),
1005  "Message size, " << msg_size
1006  << ", incompatible with pairs_to_receive, which has size "
1007  << _pairs_to_receive[pid].size());
1008  for (unsigned i = 0; i < msg_size; ++i)
1009  _kij[_pairs_to_receive[pid][i].first][_pairs_to_receive[pid][i].second] += kij_received[i];
1010  };
1011  Parallel::push_parallel_vector_data(this->comm(), kij_to_send, kij_action_functor);
1012 }
AdvectiveFluxCalculatorBase::_dDii_dKji
std::vector< std::vector< Real > > _dDii_dKji
dDii_dKji[i][j] = d(D[i][i])/d(K[j][i])
Definition: AdvectiveFluxCalculatorBase.h:296
AdvectiveFluxCalculatorBase::_nodes_to_send
std::map< processor_id_type, std::vector< dof_id_type > > _nodes_to_send
_nodes_to_send[proc_id] = list of sequential nodal IDs.
Definition: AdvectiveFluxCalculatorBase.h:226
AdvectiveFluxCalculatorBase::_drP
std::vector< std::vector< Real > > _drP
drP[i][j] = d(rP[i])/d(u[j]). Here j indexes the j^th node connected to i
Definition: AdvectiveFluxCalculatorBase.h:303
AdvectiveFluxCalculatorBase::threadJoin
virtual void threadJoin(const UserObject &uo) override
Definition: AdvectiveFluxCalculatorBase.C:251
AdvectiveFluxCalculatorBase.h
AdvectiveFluxCalculatorBase::execute
virtual void execute() override
Definition: AdvectiveFluxCalculatorBase.C:218
AdvectiveFluxCalculatorBase::_u_nodal_computed_by_thread
std::vector< bool > _u_nodal_computed_by_thread
_u_nodal_computed_by_thread(i) = true if _u_nodal[i] has been computed in execute() by the thread on ...
Definition: AdvectiveFluxCalculatorBase.h:201
AdvectiveFluxCalculatorBase::PQPlusMinusEnum::QPlus
AdvectiveFluxCalculatorBase::_nodes_to_receive
std::map< processor_id_type, std::vector< dof_id_type > > _nodes_to_receive
_nodes_to_receive[proc_id] = list of sequential nodal IDs.
Definition: AdvectiveFluxCalculatorBase.h:218
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum::MinMod
PorousFlowConnectedNodes::globalConnectionsToSequentialID
const std::vector< dof_id_type > & globalConnectionsToSequentialID(dof_id_type sequential_node_ID) const
Return all the nodes (global node IDs) connected to the given sequential node ID.
Definition: PorousFlowConnectedNodes.C:153
AdvectiveFluxCalculatorBase::_dFij_dKik
std::vector< std::vector< std::vector< Real > > > _dFij_dKik
dFij_dKik[sequential_i][j][k] = d(fa[sequential_i][j])/d(K[sequential_i][k]).
Definition: AdvectiveFluxCalculatorBase.h:321
AdvectiveFluxCalculatorBase::_flux_limiter_type
enum AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum _flux_limiter_type
AdvectiveFluxCalculatorBase::zeroedConnection
void zeroedConnection(std::map< dof_id_type, Real > &the_map, dof_id_type node_i) const
Clears the_map, then, using _kij, constructs the_map so that the_map[node_id] = 0....
Definition: AdvectiveFluxCalculatorBase.C:755
AdvectiveFluxCalculatorBase::computeVelocity
virtual Real computeVelocity(unsigned i, unsigned j, unsigned qp) const =0
Computes the transfer velocity between current node i and current node j at the current qp in the cur...
AdvectiveFluxCalculatorBase::getdFluxOutdKjk
const std::vector< std::vector< Real > > & getdFluxOutdKjk(dof_id_type node_i) const
Returns r where r[j][k] = d(flux out of global node i)/dK[connected node j][connected node k] used in...
Definition: AdvectiveFluxCalculatorBase.C:737
PorousFlowConnectedNodes::numNodes
std::size_t numNodes() const
number of nodes known by this class
Definition: PorousFlowConnectedNodes.C:31
AdvectiveFluxCalculatorBase::_drM_dk
std::vector< std::vector< Real > > _drM_dk
drM_dk[i][j] = d(rM[i])/d(K[i][j]). Here j indexes the j^th node connected to i
Definition: AdvectiveFluxCalculatorBase.h:309
AdvectiveFluxCalculatorBase::getdFluxOutdu
const std::map< dof_id_type, Real > & getdFluxOutdu(dof_id_type node_i) const
Returns r where r[j] = d(flux out of global node i)/du(global node j) used in Jacobian computations.
Definition: AdvectiveFluxCalculatorBase.C:731
AdvectiveFluxCalculatorBase::_dflux_out_du
std::vector< std::map< dof_id_type, Real > > _dflux_out_du
_dflux_out_du[i][j] = d(flux_out[i])/d(u[j]).
Definition: AdvectiveFluxCalculatorBase.h:183
AdvectiveFluxCalculatorBase::getValence
unsigned getValence(dof_id_type node_i) const
Returns the valence of the global node i Valence is the number of times the node is encountered in a ...
Definition: AdvectiveFluxCalculatorBase.C:749
AdvectiveFluxCalculatorBase::PQPlusMinusEnum
PQPlusMinusEnum
Signals to the PQPlusMinus method what should be computed.
Definition: AdvectiveFluxCalculatorBase.h:255
AdvectiveFluxCalculatorBase::PQPlusMinus
Real PQPlusMinus(dof_id_type sequential_i, const PQPlusMinusEnum pq_plus_minus, std::vector< Real > &derivs, std::vector< Real > &dpq_dk) const
Returns the value of P_{i}^{+}, P_{i}^{-}, Q_{i}^{+} or Q_{i}^{-} (depending on pq_plus_minus) which ...
Definition: AdvectiveFluxCalculatorBase.C:764
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum::superbee
AdvectiveFluxCalculatorBase::_dij
std::vector< std::vector< Real > > _dij
Vectors used in finalize()
Definition: AdvectiveFluxCalculatorBase.h:288
AdvectiveFluxCalculatorBase::executeOnElement
virtual void executeOnElement(dof_id_type global_i, dof_id_type global_j, unsigned local_i, unsigned local_j, unsigned qp)
This is called by multiple times in execute() in a double loop over _current_elem's nodes (local_i an...
Definition: AdvectiveFluxCalculatorBase.C:241
PorousFlowConnectedNodes::globalID
dof_id_type globalID(dof_id_type sequential_node_ID) const
Return the global node ID (node number in the mesh) corresponding to the provided sequential node ID.
Definition: PorousFlowConnectedNodes.C:78
PorousFlowConnectedNodes::sequentialConnectionsToSequentialID
const std::vector< dof_id_type > & sequentialConnectionsToSequentialID(dof_id_type sequential_node_ID) const
Return all the nodes (sequential node IDs) connected to the given sequential node ID All elements of ...
Definition: PorousFlowConnectedNodes.C:133
AdvectiveFluxCalculatorBase::_number_of_nodes
std::size_t _number_of_nodes
Number of nodes held by the _connections object.
Definition: AdvectiveFluxCalculatorBase.h:207
AdvectiveFluxCalculatorBase::_rP
std::vector< Real > _rP
Definition: AdvectiveFluxCalculatorBase.h:299
AdvectiveFluxCalculatorBase::_dDii_dKij
std::vector< std::vector< Real > > _dDii_dKij
dDii_dKij[i][j] = d(D[i][i])/d(K[i][j])
Definition: AdvectiveFluxCalculatorBase.h:294
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum::VanLeer
AdvectiveFluxCalculatorBase::_rM
std::vector< Real > _rM
Definition: AdvectiveFluxCalculatorBase.h:300
validParams< AdvectiveFluxCalculatorBase >
InputParameters validParams< AdvectiveFluxCalculatorBase >()
Definition: AdvectiveFluxCalculatorBase.C:19
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum::MC
AdvectiveFluxCalculatorBase::buildCommLists
virtual void buildCommLists()
When using multiple processors, other processors will compute:
Definition: AdvectiveFluxCalculatorBase.C:850
PorousFlowConnectedNodes::indexOfSequentialConnection
unsigned indexOfSequentialConnection(dof_id_type sequential_node_ID_from, dof_id_type sequential_node_ID_to) const
Return the index of sequential_node_ID_to in the sequentialConnectionsToSequentialID(sequential_node_...
Definition: PorousFlowConnectedNodes.C:189
AdvectiveFluxCalculatorBase::rPlus
Real rPlus(dof_id_type sequential_i, std::vector< Real > &dlimited_du, std::vector< Real > &dlimited_dk) const
Returns the value of R_{i}^{+}, Eqn (49) of KT.
Definition: AdvectiveFluxCalculatorBase.C:559
AdvectiveFluxCalculatorBase::timestepSetup
virtual void timestepSetup() override
Definition: AdvectiveFluxCalculatorBase.C:81
AdvectiveFluxCalculatorBase::computeU
virtual Real computeU(unsigned i) const =0
Computes the value of u at the local node id of the current element (_current_elem)
AdvectiveFluxCalculatorBase::_u_nodal
std::vector< Real > _u_nodal
_u_nodal[i] = value of _u at sequential node number i
Definition: AdvectiveFluxCalculatorBase.h:198
PorousFlowConnectedNodes::clear
void clear()
clear all data in readiness for adding global nodes and connections
Definition: PorousFlowConnectedNodes.C:17
AdvectiveFluxCalculatorBase::finalize
virtual void finalize() override
Definition: AdvectiveFluxCalculatorBase.C:270
AdvectiveFluxCalculatorBase::meshChanged
virtual void meshChanged() override
Definition: AdvectiveFluxCalculatorBase.C:198
AdvectiveFluxCalculatorBase::exchangeGhostedInfo
virtual void exchangeGhostedInfo()
Sends and receives multi-processor information regarding u_nodal and k_ij.
Definition: AdvectiveFluxCalculatorBase.C:968
AdvectiveFluxCalculatorBase::_my_pid
processor_id_type _my_pid
processor ID of this object
Definition: AdvectiveFluxCalculatorBase.h:210
name
const std::string name
Definition: Setup.h:21
AdvectiveFluxCalculatorBase::_connections
PorousFlowConnectedNodes _connections
Holds the sequential and global nodal IDs, and info regarding mesh connections between them.
Definition: AdvectiveFluxCalculatorBase.h:204
AdvectiveFluxCalculatorBase::_resizing_needed
bool _resizing_needed
whether _kij, etc, need to be sized appropriately (and valence recomputed) at the start of the timest...
Definition: AdvectiveFluxCalculatorBase.h:130
AdvectiveFluxCalculatorBase::AdvectiveFluxCalculatorBase
AdvectiveFluxCalculatorBase(const InputParameters &parameters)
Definition: AdvectiveFluxCalculatorBase.C:51
AdvectiveFluxCalculatorBase::initialize
virtual void initialize() override
Definition: AdvectiveFluxCalculatorBase.C:207
PorousFlowConnectedNodes::indexOfGlobalConnection
unsigned indexOfGlobalConnection(dof_id_type global_node_ID_from, dof_id_type global_node_ID_to) const
Return the index of global_node_ID_to in the globalConnectionsToGlobalID(global_node_ID_from) vector.
Definition: PorousFlowConnectedNodes.C:172
PorousFlowConnectedNodes::globalConnectionsToGlobalID
const std::vector< dof_id_type > & globalConnectionsToGlobalID(dof_id_type global_node_ID) const
Return all the nodes (global node IDs) connected to the given global node ID.
Definition: PorousFlowConnectedNodes.C:143
AdvectiveFluxCalculatorBase::limitFlux
void limitFlux(Real a, Real b, Real &limited, Real &dlimited_db) const
flux limiter, L, on Page 135 of Kuzmin and Turek
Definition: AdvectiveFluxCalculatorBase.C:631
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum::None
AdvectiveFluxCalculatorBase
Base class to compute Advective fluxes.
Definition: AdvectiveFluxCalculatorBase.h:33
PorousFlowConnectedNodes::finalizeAddingGlobalNodes
void finalizeAddingGlobalNodes()
Signal that all global node IDs have been added to the internal data structures.
Definition: PorousFlowConnectedNodes.C:50
AdvectiveFluxCalculatorBase::_pairs_to_send
std::map< processor_id_type, std::vector< std::pair< dof_id_type, dof_id_type > > > _pairs_to_send
_pairs_to_send[proc_id] indicates the k(i, j) pairs that we will send to proc_id _pairs_to_send is fi...
Definition: AdvectiveFluxCalculatorBase.h:244
AdvectiveFluxCalculatorBase::_kij
std::vector< std::vector< Real > > _kij
Kuzmin-Turek K_ij matrix.
Definition: AdvectiveFluxCalculatorBase.h:176
AdvectiveFluxCalculatorBase::_pairs_to_receive
std::map< processor_id_type, std::vector< std::pair< dof_id_type, dof_id_type > > > _pairs_to_receive
_pairs_to_receive[proc_id] indicates the k(i, j) pairs that will be sent to us from proc_id _pairs_to...
Definition: AdvectiveFluxCalculatorBase.h:235
AdvectiveFluxCalculatorBase::getFluxOut
Real getFluxOut(dof_id_type node_i) const
Returns the flux out of lobal node id.
Definition: AdvectiveFluxCalculatorBase.C:743
PorousFlowConnectedNodes::sequentialID
dof_id_type sequentialID(dof_id_type global_node_ID) const
Return the sequential node ID corresponding to the global node ID This is guaranteed to lie in the ra...
Definition: PorousFlowConnectedNodes.C:87
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum
FluxLimiterTypeEnum
Determines Flux Limiter type (Page 135 of Kuzmin and Turek) "None" means that limitFlux=0 always,...
Definition: AdvectiveFluxCalculatorBase.h:169
AdvectiveFluxCalculatorBase::_allowable_MB_wastage
const Real _allowable_MB_wastage
A mooseWarning is issued if mb_wasted = (_connections.sizeSequential() - _connections....
Definition: AdvectiveFluxCalculatorBase.h:252
PorousFlowConnectedNodes::finalizeAddingConnections
void finalizeAddingConnections()
Signal that all global node IDs have been added to the internal data structures.
Definition: PorousFlowConnectedNodes.C:108
AdvectiveFluxCalculatorBase::_drP_dk
std::vector< std::vector< Real > > _drP_dk
drP_dk[i][j] = d(rP[i])/d(K[i][j]). Here j indexes the j^th node connected to i
Definition: AdvectiveFluxCalculatorBase.h:307
AdvectiveFluxCalculatorBase::_dDij_dKji
std::vector< std::vector< Real > > _dDij_dKji
dDij_dKji[i][j] = d(D[i][j])/d(K[j][i]) for i!=j
Definition: AdvectiveFluxCalculatorBase.h:292
AdvectiveFluxCalculatorBase::_flux_out
std::vector< Real > _flux_out
_flux_out[i] = flux of "heat" from sequential node i
Definition: AdvectiveFluxCalculatorBase.h:179
PorousFlowConnectedNodes::addGlobalNode
void addGlobalNode(dof_id_type global_node_ID)
Add the given global_node_ID to the internal data structures If the global node ID has already been a...
Definition: PorousFlowConnectedNodes.C:39
PorousFlowConnectedNodes::addConnection
void addConnection(dof_id_type global_node_from, dof_id_type global_node_to)
Specifies that global_node_to is connected to global_node_from.
Definition: PorousFlowConnectedNodes.C:96
AdvectiveFluxCalculatorBase::rMinus
Real rMinus(dof_id_type sequential_i, std::vector< Real > &dlimited_du, std::vector< Real > &dlimited_dk) const
Returns the value of R_{i}^{-}, Eqn (49) of KT.
Definition: AdvectiveFluxCalculatorBase.C:595
AdvectiveFluxCalculatorBase::_lij
std::vector< std::vector< Real > > _lij
Definition: AdvectiveFluxCalculatorBase.h:298
AdvectiveFluxCalculatorBase::PQPlusMinusEnum::PPlus
AdvectiveFluxCalculatorBase::PQPlusMinusEnum::PMinus
PorousFlowConnectedNodes::sizeSequential
std::size_t sizeSequential() const
Return the size of _sequential_id, for checking memory efficiency.
Definition: PorousFlowConnectedNodes.C:69
AdvectiveFluxCalculatorBase::_dfa
std::vector< std::vector< std::map< dof_id_type, Real > > > _dfa
dfa[sequential_i][j][global_k] = d(fa[sequential_i][j])/du[global_k].
Definition: AdvectiveFluxCalculatorBase.h:316
AdvectiveFluxCalculatorBase::_drM
std::vector< std::vector< Real > > _drM
drM[i][j] = d(rM[i])/d(u[j]). Here j indexes the j^th node connected to i
Definition: AdvectiveFluxCalculatorBase.h:305
AdvectiveFluxCalculatorBase::_dflux_out_dKjk
std::vector< std::vector< std::vector< Real > > > _dflux_out_dKjk
_dflux_out_dKjk[sequential_i][j][k] = d(flux_out[sequential_i])/d(K[j][k]).
Definition: AdvectiveFluxCalculatorBase.h:190
AdvectiveFluxCalculatorBase::PQPlusMinusEnum::QMinus
AdvectiveFluxCalculatorBase::_dFij_dKjk
std::vector< std::vector< std::vector< Real > > > _dFij_dKjk
dFij_dKjk[sequential_i][j][k] = d(fa[sequential_i][j])/d(K[sequential_j][k]).
Definition: AdvectiveFluxCalculatorBase.h:326
AdvectiveFluxCalculatorBase::_valence
std::vector< unsigned > _valence
_valence[i] = number of times, in a loop over elements seen by this processor (viz,...
Definition: AdvectiveFluxCalculatorBase.h:195
AdvectiveFluxCalculatorBase::_dDij_dKij
std::vector< std::vector< Real > > _dDij_dKij
dDij_dKij[i][j] = d(D[i][j])/d(K[i][j]) for i!=j
Definition: AdvectiveFluxCalculatorBase.h:290
AdvectiveFluxCalculatorBase::_fa
std::vector< std::vector< Real > > _fa
fa[sequential_i][j] sequential_j is the j^th connection to sequential_i
Definition: AdvectiveFluxCalculatorBase.h:312