www.mooseframework.org
PorousFlowAdvectiveFluxCalculatorBase.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 "Assembly.h"
12 #include "libmesh/string_to_enum.h"
13 #include "libmesh/parallel_sync.h"
14 
15 template <>
16 InputParameters
18 {
19  InputParameters params = validParams<AdvectiveFluxCalculatorBase>();
20  params.addClassDescription(
21  "Base class to compute the advective flux of fluid in PorousFlow situations. The velocity "
22  "is U * (-permeability * (grad(P) - density * gravity)), while derived classes define U. "
23  "The Kuzmin-Turek FEM-TVD multidimensional stabilization scheme is used");
24  params.addRequiredParam<RealVectorValue>("gravity",
25  "Gravitational acceleration vector downwards (m/s^2)");
26  params.addRequiredParam<UserObjectName>(
27  "PorousFlowDictator", "The UserObject that holds the list of PorousFlow variable names");
28  params.addParam<unsigned int>(
29  "phase", 0, "The index corresponding to the phase for this UserObject");
30  MooseEnum families("LAGRANGE MONOMIAL HERMITE SCALAR HIERARCHIC CLOUGH XYZ SZABAB BERNSTEIN");
31  params.addParam<MooseEnum>(
32  "fe_family",
33  families,
34  "FE Family to use (eg Lagrange). You only need to specify this is your porous_flow_vars in "
35  "your PorousFlowDictator have different families or orders");
36  MooseEnum orders("CONSTANT FIRST SECOND THIRD FOURTH");
37  params.addParam<MooseEnum>(
38  "fe_order",
39  orders,
40  "FE Order to use (eg First). You only need to specify this is your porous_flow_vars in your "
41  "PorousFlowDictator have different families or orders");
42  return params;
43 }
44 
46  const InputParameters & parameters)
47  : AdvectiveFluxCalculatorBase(parameters),
48  _dictator(getUserObject<PorousFlowDictator>("PorousFlowDictator")),
49  _num_vars(_dictator.numVariables()),
50  _gravity(getParam<RealVectorValue>("gravity")),
51  _phase(getParam<unsigned int>("phase")),
52  _permeability(getMaterialProperty<RealTensorValue>("PorousFlow_permeability_qp")),
53  _dpermeability_dvar(
54  getMaterialProperty<std::vector<RealTensorValue>>("dPorousFlow_permeability_qp_dvar")),
55  _dpermeability_dgradvar(getMaterialProperty<std::vector<std::vector<RealTensorValue>>>(
56  "dPorousFlow_permeability_qp_dgradvar")),
57  _fluid_density_qp(getMaterialProperty<std::vector<Real>>("PorousFlow_fluid_phase_density_qp")),
58  _dfluid_density_qp_dvar(getMaterialProperty<std::vector<std::vector<Real>>>(
59  "dPorousFlow_fluid_phase_density_qp_dvar")),
60  _grad_p(getMaterialProperty<std::vector<RealGradient>>("PorousFlow_grad_porepressure_qp")),
61  _dgrad_p_dgrad_var(getMaterialProperty<std::vector<std::vector<Real>>>(
62  "dPorousFlow_grad_porepressure_qp_dgradvar")),
63  _dgrad_p_dvar(getMaterialProperty<std::vector<std::vector<RealGradient>>>(
64  "dPorousFlow_grad_porepressure_qp_dvar")),
65  _fe_type(isParamValid("fe_family") && isParamValid("fe_order")
66  ? FEType(Utility::string_to_enum<Order>(getParam<MooseEnum>("fe_order")),
67  Utility::string_to_enum<FEFamily>(getParam<MooseEnum>("fe_family")))
68  : _dictator.feType()),
69  _phi(_assembly.fePhi<Real>(_fe_type)),
70  _grad_phi(_assembly.feGradPhi<Real>(_fe_type)),
71  _du_dvar(),
72  _du_dvar_computed_by_thread(),
73  _dkij_dvar(),
74  _dflux_out_dvars(),
75  _triples_to_receive(),
76  _triples_to_send(),
77  _perm_derivs(_dictator.usePermDerivs())
78 {
79  if (_phase >= _dictator.numPhases())
80  paramError("phase",
81  "Phase number entered is greater than the number of phases specified in the "
82  "Dictator. Remember that indexing starts at 0");
83 
84  if (isParamValid("fe_family") && !isParamValid("fe_order"))
85  paramError("fe_order", "If you specify fe_family you must also specify fe_order");
86  if (isParamValid("fe_order") && !isParamValid("fe_family"))
87  paramError("fe_family", "If you specify fe_order you must also specify fe_family");
88  if (!_dictator.consistentFEType() && !isParamValid("fe_family"))
89  paramError("fe_family",
90  "The PorousFlowDictator cannot determine the appropriate FE type to use because "
91  "your porous_flow_vars are of different types. You must specify the appropriate "
92  "fe_family and fe_order to use.");
93 }
94 
95 Real
96 PorousFlowAdvectiveFluxCalculatorBase::computeVelocity(unsigned i, unsigned j, unsigned qp) const
97 {
98  // The following is but one choice for PorousFlow situations
99  // If you change this, you will probably have to change
100  // - the derivative in executeOnElement
101  // - computeU
102  // - computedU_dvar
103  return -_grad_phi[i][qp] *
104  (_permeability[qp] * (_grad_p[qp][_phase] - _fluid_density_qp[qp][_phase] * _gravity)) *
105  _phi[j][qp];
106 }
107 
108 void
110  dof_id_type global_i, dof_id_type global_j, unsigned local_i, unsigned local_j, unsigned qp)
111 {
112  AdvectiveFluxCalculatorBase::executeOnElement(global_i, global_j, local_i, local_j, qp);
113  const dof_id_type sequential_i = _connections.sequentialID(global_i);
114  const unsigned j = _connections.indexOfGlobalConnection(global_i, global_j);
115 
116  // compute d(Kij)/d(porous_flow_variables)
117  for (unsigned local_k = 0; local_k < _current_elem->n_nodes(); ++local_k)
118  {
119  const dof_id_type global_k = _current_elem->node_id(local_k);
120  for (unsigned pvar = 0; pvar < _num_vars; ++pvar)
121  {
122  RealVectorValue deriv =
123  _permeability[qp] *
124  (_grad_phi[local_k][qp] * _dgrad_p_dgrad_var[qp][_phase][pvar] -
125  _phi[local_k][qp] * _dfluid_density_qp_dvar[qp][_phase][pvar] * _gravity);
126  deriv += _permeability[qp] * (_dgrad_p_dvar[qp][_phase][pvar] * _phi[local_k][qp]);
127 
128  if (_perm_derivs)
129  {
130  deriv += _dpermeability_dvar[qp][pvar] * _phi[local_k][qp] *
132  for (unsigned i = 0; i < LIBMESH_DIM; ++i)
133  deriv += _dpermeability_dgradvar[qp][i][pvar] * _grad_phi[local_k][qp](i) *
135  }
136 
137  _dkij_dvar[sequential_i][j][global_k][pvar] +=
138  _JxW[qp] * _coord[qp] * (-_grad_phi[local_i][qp] * deriv * _phi[local_j][qp]);
139  }
140  }
141 }
142 
143 void
145 {
146  const bool resizing_was_needed =
147  _resizing_needed; // _resizing_needed gets set to false at the end of
148  // AdvectiveFluxCalculatorBase::timestepSetup()
150 
151  // clear and appropriately size all the derivative info
152  // d(U)/d(porous_flow_variables) and
153  // d(Kij)/d(porous_flow_variables) and
154  // d(flux_out)/d(porous_flow_variables)
155  if (resizing_was_needed)
156  {
157  const std::size_t num_nodes = _connections.numNodes();
158  _du_dvar.assign(num_nodes, std::vector<Real>(_num_vars, 0.0));
159  _du_dvar_computed_by_thread.assign(num_nodes, false);
160  _dflux_out_dvars.assign(num_nodes, std::map<dof_id_type, std::vector<Real>>());
161  _dkij_dvar.resize(num_nodes);
162  for (dof_id_type sequential_i = 0; sequential_i < num_nodes; ++sequential_i)
163  {
164  const std::vector<dof_id_type> con_i =
166  const std::size_t num_con_i = con_i.size();
167  _dkij_dvar[sequential_i].assign(num_con_i, std::map<dof_id_type, std::vector<Real>>());
168  for (unsigned j = 0; j < num_con_i; ++j)
169  for (const auto & global_neighbor_to_i : con_i)
170  _dkij_dvar[sequential_i][j][global_neighbor_to_i] = std::vector<Real>(_num_vars, 0.0);
171  }
172  }
173 }
174 
175 void
177 {
179  const std::size_t num_nodes = _connections.numNodes();
180  _du_dvar_computed_by_thread.assign(num_nodes, false);
181  for (dof_id_type sequential_i = 0; sequential_i < num_nodes; ++sequential_i)
182  {
183  const std::vector<dof_id_type> & con_i =
185  const std::size_t num_con_i = con_i.size();
186  for (unsigned j = 0; j < num_con_i; ++j)
187  for (const auto & global_neighbor_to_i : con_i)
188  _dkij_dvar[sequential_i][j][global_neighbor_to_i] = std::vector<Real>(_num_vars, 0.0);
189  }
190 }
191 
192 void
194 {
196 
197  // compute d(U)/d(porous_flow_variables) for nodes in _current_elem and for this
198  // execution thread. In threadJoin all these computations get gathered
199  // using _du_dvar_computed_by_thread
200  for (unsigned i = 0; i < _current_elem->n_nodes(); ++i)
201  {
202  const dof_id_type global_i = _current_elem->node_id(i);
203  const dof_id_type sequential_i = _connections.sequentialID(global_i);
204  if (_du_dvar_computed_by_thread[sequential_i])
205  continue;
206  for (unsigned pvar = 0; pvar < _num_vars; ++pvar)
207  _du_dvar[sequential_i][pvar] = computedU_dvar(i, pvar);
208  _du_dvar_computed_by_thread[sequential_i] = true;
209  }
210 }
211 
212 void
214 {
217  static_cast<const PorousFlowAdvectiveFluxCalculatorBase &>(uo);
218  // add the values of _dkij_dvar computed by different threads
219  const std::size_t num_nodes = _connections.numNodes();
220  for (dof_id_type sequential_i = 0; sequential_i < num_nodes; ++sequential_i)
221  {
222  const std::vector<dof_id_type> & con_i =
224  const std::size_t num_con_i = con_i.size();
225  for (unsigned j = 0; j < num_con_i; ++j)
226  for (const auto & global_derivs : pfafc._dkij_dvar[sequential_i][j])
227  for (unsigned pvar = 0; pvar < _num_vars; ++pvar)
228  _dkij_dvar[sequential_i][j][global_derivs.first][pvar] += global_derivs.second[pvar];
229  }
230 
231  // gather the values of _du_dvar computed by different threads
232  for (dof_id_type sequential_i = 0; sequential_i < _number_of_nodes; ++sequential_i)
233  if (!_du_dvar_computed_by_thread[sequential_i] &&
234  pfafc._du_dvar_computed_by_thread[sequential_i])
235  _du_dvar[sequential_i] = pfafc._du_dvar[sequential_i];
236 }
237 
238 const std::map<dof_id_type, std::vector<Real>> &
239 PorousFlowAdvectiveFluxCalculatorBase::getdK_dvar(dof_id_type node_i, dof_id_type node_j) const
240 {
241  const dof_id_type sequential_i = _connections.sequentialID(node_i);
242  const unsigned j = _connections.indexOfGlobalConnection(node_i, node_j);
243  return _dkij_dvar[sequential_i][j];
244 }
245 
246 const std::map<dof_id_type, std::vector<Real>> &
248 {
249  return _dflux_out_dvars[_connections.sequentialID(node_id)];
250 }
251 
252 void
254 {
256 
257  // compute d(flux_out)/d(porous flow variable)
259  for (const auto & node_i : _connections.globalIDs())
260  {
261  const dof_id_type sequential_i = _connections.sequentialID(node_i);
262  _dflux_out_dvars[sequential_i].clear();
263 
264  const std::map<dof_id_type, Real> & dflux_out_du =
266  for (const auto & node_du : dflux_out_du)
267  {
268  const dof_id_type j = node_du.first;
269  const Real dflux_out_du_j = node_du.second;
270  _dflux_out_dvars[sequential_i][j] = _du_dvar[_connections.sequentialID(j)];
271  for (unsigned pvar = 0; pvar < _num_vars; ++pvar)
272  _dflux_out_dvars[sequential_i][j][pvar] *= dflux_out_du_j;
273  }
274 
275  // _dflux_out_dvars is now sized correctly, because getdFluxOutdu(i) contains all nodes
276  // connected to i and all nodes connected to nodes connected to i. The
277  // getdFluxOutdKij contains no extra nodes, so just += the dflux/dK terms
278  const std::vector<std::vector<Real>> & dflux_out_dKjk =
280  const std::vector<dof_id_type> & con_i = _connections.globalConnectionsToGlobalID(node_i);
281  for (std::size_t index_j = 0; index_j < con_i.size(); ++index_j)
282  {
283  const dof_id_type node_j = con_i[index_j];
284  const std::vector<dof_id_type> & con_j = _connections.globalConnectionsToGlobalID(node_j);
285  for (std::size_t index_k = 0; index_k < con_j.size(); ++index_k)
286  {
287  const dof_id_type node_k = con_j[index_k];
288  const Real dflux_out_dK_jk = dflux_out_dKjk[index_j][index_k];
289  const std::map<dof_id_type, std::vector<Real>> & dkj_dvarl = getdK_dvar(node_j, node_k);
290  for (const auto & nodel_deriv : dkj_dvarl)
291  {
292  const dof_id_type l = nodel_deriv.first;
293  for (unsigned pvar = 0; pvar < _num_vars; ++pvar)
294  _dflux_out_dvars[sequential_i][l][pvar] += dflux_out_dK_jk * nodel_deriv.second[pvar];
295  }
296  }
297  }
298  }
299 }
300 
301 void
303 {
304  // build nodes and pairs to exchange
306 
307  // Build _triples_to_receive
308  // tpr_global is essentially _triples_to_receive, but its key-pairs have global nodal IDs: later
309  // we build _triples_to_receive by flattening the data structure and making these key-pairs into
310  // sequential nodal IDs
311  std::map<processor_id_type,
312  std::map<std::pair<dof_id_type, dof_id_type>, std::vector<dof_id_type>>>
313  tpr_global;
314  for (const auto & elem : _fe_problem.getEvaluableElementRange())
315  if (this->hasBlocks(elem->subdomain_id()))
316  {
317  const processor_id_type elem_pid = elem->processor_id();
318  if (elem_pid != _my_pid)
319  {
320  if (tpr_global.find(elem_pid) == tpr_global.end())
321  tpr_global[elem_pid] =
322  std::map<std::pair<dof_id_type, dof_id_type>, std::vector<dof_id_type>>();
323  for (unsigned i = 0; i < elem->n_nodes(); ++i)
324  for (unsigned j = 0; j < elem->n_nodes(); ++j)
325  {
326  std::pair<dof_id_type, dof_id_type> the_pair(elem->node_id(i), elem->node_id(j));
327  if (tpr_global[elem_pid].find(the_pair) == tpr_global[elem_pid].end())
328  tpr_global[elem_pid][the_pair] = std::vector<dof_id_type>();
329 
330  for (const auto & global_neighbor_to_i :
331  _connections.globalConnectionsToGlobalID(elem->node_id(i)))
332  if (std::find(tpr_global[elem_pid][the_pair].begin(),
333  tpr_global[elem_pid][the_pair].end(),
334  global_neighbor_to_i) == tpr_global[elem_pid][the_pair].end())
335  tpr_global[elem_pid][the_pair].push_back(global_neighbor_to_i);
336  }
337  }
338  }
339 
340  // flattening makes later manipulations a lot more concise. Store the result in
341  // _triples_to_receive
342  _triples_to_receive.clear();
343  for (const auto & kv : tpr_global)
344  {
345  const processor_id_type pid = kv.first;
346  _triples_to_receive[pid] = std::vector<dof_id_type>();
347  for (const auto & pr_vec : kv.second)
348  {
349  const dof_id_type i = pr_vec.first.first;
350  const dof_id_type j = pr_vec.first.second;
351  for (const auto & global_nd : pr_vec.second)
352  {
353  _triples_to_receive[pid].push_back(i);
354  _triples_to_receive[pid].push_back(j);
355  _triples_to_receive[pid].push_back(global_nd);
356  }
357  }
358  }
359 
360  _triples_to_send.clear();
361  auto triples_action_functor = [this](processor_id_type pid,
362  const std::vector<dof_id_type> & tts) {
363  _triples_to_send[pid] = tts;
364  };
365  Parallel::push_parallel_vector_data(this->comm(), _triples_to_receive, triples_action_functor);
366 
367  // _triples_to_send and _triples_to_receive have been built using global node IDs
368  // since all processors know about that. However, using global IDs means
369  // that every time we send/receive, we keep having to do things like
370  // _dkij_dvar[_connections.sequentialID(_triples_to_send[pid][i])][_connections.indexOfGlobalConnection(_triples_to_send[pid][i],
371  // _triples_to_send[pid][i + 1])] which is quite inefficient. So:
372  for (auto & kv : _triples_to_send)
373  {
374  const processor_id_type pid = kv.first;
375  const std::size_t num = kv.second.size();
376  for (std::size_t i = 0; i < num; i += 3)
377  {
379  _triples_to_send[pid][i], _triples_to_send[pid][i + 1]);
381  }
382  }
383  for (auto & kv : _triples_to_receive)
384  {
385  const processor_id_type pid = kv.first;
386  const std::size_t num = kv.second.size();
387  for (std::size_t i = 0; i < num; i += 3)
388  {
390  _triples_to_receive[pid][i], _triples_to_receive[pid][i + 1]);
392  }
393  }
394 }
395 
396 void
398 {
399  // Exchange u_nodal and k_ij
401 
402  // Exchange _du_dvar
403  std::map<processor_id_type, std::vector<std::vector<Real>>> du_dvar_to_send;
404  for (const auto & kv : _nodes_to_send)
405  {
406  const processor_id_type pid = kv.first;
407  du_dvar_to_send[pid] = std::vector<std::vector<Real>>();
408  for (const auto & nd : kv.second)
409  du_dvar_to_send[pid].push_back(_du_dvar[nd]);
410  }
411 
412  auto du_action_functor = [this](processor_id_type pid,
413  const std::vector<std::vector<Real>> & du_dvar_received) {
414  const std::size_t msg_size = du_dvar_received.size();
415  mooseAssert(
416  msg_size == _nodes_to_receive[pid].size(),
417  "Message size, "
418  << msg_size
419  << ", in du_dvar communication is incompatible with nodes_to_receive, which has size "
420  << _nodes_to_receive[pid].size());
421  for (unsigned i = 0; i < msg_size; ++i)
422  _du_dvar[_nodes_to_receive[pid][i]] = du_dvar_received[i];
423  };
424  Parallel::push_parallel_vector_data(this->comm(), du_dvar_to_send, du_action_functor);
425 
426  // Exchange _dkij_dvar
427  std::map<processor_id_type, std::vector<std::vector<Real>>> dkij_dvar_to_send;
428  for (const auto & kv : _triples_to_send)
429  {
430  const processor_id_type pid = kv.first;
431  dkij_dvar_to_send[pid] = std::vector<std::vector<Real>>();
432  const std::size_t num = kv.second.size();
433  for (std::size_t i = 0; i < num; i += 3)
434  {
435  const dof_id_type sequential_id = kv.second[i];
436  const unsigned index_to_seq = kv.second[i + 1];
437  const dof_id_type global_id = kv.second[i + 2];
438  dkij_dvar_to_send[pid].push_back(_dkij_dvar[sequential_id][index_to_seq][global_id]);
439  }
440  }
441 
442  auto dk_action_functor = [this](processor_id_type pid,
443  const std::vector<std::vector<Real>> & dkij_dvar_received) {
444  const std::size_t num = _triples_to_receive[pid].size();
445  mooseAssert(dkij_dvar_received.size() == num / 3,
446  "Message size, " << dkij_dvar_received.size()
447  << ", in dkij_dvar communication is incompatible with "
448  "triples_to_receive, which has size "
449  << _triples_to_receive[pid].size());
450  for (std::size_t i = 0; i < num; i += 3)
451  {
452  const dof_id_type sequential_id = _triples_to_receive[pid][i];
453  const unsigned index_to_seq = _triples_to_receive[pid][i + 1];
454  const dof_id_type global_id = _triples_to_receive[pid][i + 2];
455  for (unsigned pvar = 0; pvar < _num_vars; ++pvar)
456  _dkij_dvar[sequential_id][index_to_seq][global_id][pvar] += dkij_dvar_received[i / 3][pvar];
457  }
458  };
459  Parallel::push_parallel_vector_data(this->comm(), dkij_dvar_to_send, dk_action_functor);
460 }
PorousFlowAdvectiveFluxCalculatorBase::_triples_to_send
std::map< processor_id_type, std::vector< dof_id_type > > _triples_to_send
_triples_to_send[proc_id] indicates the dk(i, j)/du_nodal information that we will send to proc_id.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:158
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
PorousFlowAdvectiveFluxCalculatorBase::initialize
virtual void initialize() override
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:176
AdvectiveFluxCalculatorBase::threadJoin
virtual void threadJoin(const UserObject &uo) override
Definition: AdvectiveFluxCalculatorBase.C:251
PorousFlowAdvectiveFluxCalculatorBase::buildCommLists
virtual void buildCommLists() override
When using multiple processors, other processors will compute:
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:302
PorousFlowAdvectiveFluxCalculatorBase::_fluid_density_qp
const MaterialProperty< std::vector< Real > > & _fluid_density_qp
Fluid density for each phase (at the qp)
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:99
AdvectiveFluxCalculatorBase::execute
virtual void execute() override
Definition: AdvectiveFluxCalculatorBase.C:218
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
validParams< AdvectiveFluxCalculatorBase >
InputParameters validParams< AdvectiveFluxCalculatorBase >()
Definition: AdvectiveFluxCalculatorBase.C:19
PorousFlowAdvectiveFluxCalculatorBase::getdK_dvar
const std::map< dof_id_type, std::vector< Real > > & getdK_dvar(dof_id_type node_i, dof_id_type node_j) const
Returns, r, where r[global node k][a] = d(K[node_i][node_j])/d(porous_flow_variable[global node k][po...
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:239
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
libMesh::RealGradient
VectorValue< Real > RealGradient
Definition: GrainForceAndTorqueInterface.h:17
PorousFlowAdvectiveFluxCalculatorBase::computedU_dvar
virtual Real computedU_dvar(unsigned i, unsigned pvar) const =0
Compute d(u)/d(porous_flow_variable)
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
PorousFlowAdvectiveFluxCalculatorBase::_du_dvar_computed_by_thread
std::vector< bool > _du_dvar_computed_by_thread
Whether _du_dvar has been computed by the local thread.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:126
PorousFlowConnectedNodes::numNodes
std::size_t numNodes() const
number of nodes known by this class
Definition: PorousFlowConnectedNodes.C:31
PorousFlowAdvectiveFluxCalculatorBase::_phi
const VariablePhiValue & _phi
Kuzmin-Turek shape function.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:117
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
PorousFlowAdvectiveFluxCalculatorBase::executeOnElement
virtual void executeOnElement(dof_id_type global_i, dof_id_type global_j, unsigned local_i, unsigned local_j, unsigned qp) override
This is called by multiple times in execute() in a double loop over _current_elem's nodes (local_i an...
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:109
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
AdvectiveFluxCalculatorBase::_number_of_nodes
std::size_t _number_of_nodes
Number of nodes held by the _connections object.
Definition: AdvectiveFluxCalculatorBase.h:207
PorousFlowAdvectiveFluxCalculatorBase::_dflux_out_dvars
std::vector< std::map< dof_id_type, std::vector< Real > > > _dflux_out_dvars
_dflux_out_dvars[sequential_i][global_j][pvar] = d(flux_out[global version of sequential_i])/d(porous...
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:136
AdvectiveFluxCalculatorBase::buildCommLists
virtual void buildCommLists()
When using multiple processors, other processors will compute:
Definition: AdvectiveFluxCalculatorBase.C:850
PorousFlowAdvectiveFluxCalculatorBase::_grad_p
const MaterialProperty< std::vector< RealGradient > > & _grad_p
Gradient of the pore pressure in each phase.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:105
PorousFlowDictator
This holds maps between the nonlinear variables used in a PorousFlow simulation and the variable numb...
Definition: PorousFlowDictator.h:71
PorousFlowDictator::numPhases
unsigned int numPhases() const
The number of fluid phases.
Definition: PorousFlowDictator.C:105
PorousFlowAdvectiveFluxCalculatorBase::computeVelocity
virtual Real computeVelocity(unsigned i, unsigned j, unsigned qp) const override
Computes the transfer velocity between current node i and current node j at the current qp in the cur...
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:96
AdvectiveFluxCalculatorBase::timestepSetup
virtual void timestepSetup() override
Definition: AdvectiveFluxCalculatorBase.C:81
PorousFlowAdvectiveFluxCalculatorBase::_dfluid_density_qp_dvar
const MaterialProperty< std::vector< std::vector< Real > > > & _dfluid_density_qp_dvar
Derivative of the fluid density for each phase wrt PorousFlow variables (at the qp)
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:102
AdvectiveFluxCalculatorBase::finalize
virtual void finalize() override
Definition: AdvectiveFluxCalculatorBase.C:270
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
PorousFlowAdvectiveFluxCalculatorBase::threadJoin
virtual void threadJoin(const UserObject &uo) override
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:213
PorousFlowAdvectiveFluxCalculatorBase::getdFluxOut_dvars
const std::map< dof_id_type, std::vector< Real > > & getdFluxOut_dvars(unsigned node_id) const
Returns d(flux_out)/d(porous_flow_variables.
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:247
AdvectiveFluxCalculatorBase::_connections
PorousFlowConnectedNodes _connections
Holds the sequential and global nodal IDs, and info regarding mesh connections between them.
Definition: AdvectiveFluxCalculatorBase.h:204
PorousFlowAdvectiveFluxCalculatorBase::_dgrad_p_dgrad_var
const MaterialProperty< std::vector< std::vector< Real > > > & _dgrad_p_dgrad_var
Derivative of Grad porepressure in each phase wrt grad(PorousFlow variables)
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:108
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::initialize
virtual void initialize() override
Definition: AdvectiveFluxCalculatorBase.C:207
PorousFlowAdvectiveFluxCalculatorBase::_num_vars
const unsigned _num_vars
Number of PorousFlow variables.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:81
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
PorousFlowAdvectiveFluxCalculatorBase::_dgrad_p_dvar
const MaterialProperty< std::vector< std::vector< RealGradient > > > & _dgrad_p_dvar
Derivative of Grad porepressure in each phase wrt PorousFlow variables.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:111
PorousFlowAdvectiveFluxCalculatorBase::_dkij_dvar
std::vector< std::vector< std::map< dof_id_type, std::vector< Real > > > > _dkij_dvar
_dkij_dvar[sequential_i][j][global_k][a] = d(K[sequential_i][j])/d(porous_flow_variable[global_k][por...
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:133
PorousFlowAdvectiveFluxCalculatorBase::_gravity
const RealVectorValue _gravity
Gravity.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:84
AdvectiveFluxCalculatorBase
Base class to compute Advective fluxes.
Definition: AdvectiveFluxCalculatorBase.h:33
PorousFlowAdvectiveFluxCalculatorBase.h
PorousFlowAdvectiveFluxCalculatorBase::_dictator
const PorousFlowDictator & _dictator
PorousFlowDictator UserObject.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:78
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
PorousFlowAdvectiveFluxCalculatorBase::_grad_phi
const VariablePhiGradient & _grad_phi
grad(Kuzmin-Turek shape function)
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:120
PorousFlowAdvectiveFluxCalculatorBase::_phase
const unsigned int _phase
The phase.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:87
PorousFlowAdvectiveFluxCalculatorBase::_dpermeability_dgradvar
const MaterialProperty< std::vector< std::vector< RealTensorValue > > > & _dpermeability_dgradvar
d(permeabiity)/d(grad(PorousFlow variable))
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:96
PorousFlowAdvectiveFluxCalculatorBase
Base class to compute the advective flux of fluid in PorousFlow situations using the Kuzmin-Turek FEM...
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:32
PorousFlowConnectedNodes::globalIDs
const std::vector< dof_id_type > & globalIDs() const
Vector of all global node IDs (node numbers in the mesh)
Definition: PorousFlowConnectedNodes.C:163
PorousFlowAdvectiveFluxCalculatorBase::_dpermeability_dvar
const MaterialProperty< std::vector< RealTensorValue > > & _dpermeability_dvar
d(permeabiity)/d(PorousFlow variable)
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:93
PorousFlowAdvectiveFluxCalculatorBase::_du_dvar
std::vector< std::vector< Real > > _du_dvar
_du_dvar[sequential_i][a] = d(u[global version of sequential node i])/d(porous_flow_variable[a])
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:123
PorousFlowAdvectiveFluxCalculatorBase::exchangeGhostedInfo
virtual void exchangeGhostedInfo() override
Sends and receives multi-processor information regarding u_nodal and k_ij.
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:397
PorousFlowAdvectiveFluxCalculatorBase::_triples_to_receive
std::map< processor_id_type, std::vector< dof_id_type > > _triples_to_receive
_triples_to_receive[proc_id] indicates the dk(i, j)/du_nodal information that we will receive from pr...
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:147
PorousFlowAdvectiveFluxCalculatorBase::timestepSetup
virtual void timestepSetup() override
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:144
PorousFlowDictator::consistentFEType
bool consistentFEType() const
Whether the porous_flow_vars all have the same FEType or if no porous_flow_vars were provided.
Definition: PorousFlowDictator.C:167
PorousFlowAdvectiveFluxCalculatorBase::finalize
virtual void finalize() override
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:253
validParams< PorousFlowAdvectiveFluxCalculatorBase >
InputParameters validParams< PorousFlowAdvectiveFluxCalculatorBase >()
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:17
PorousFlowAdvectiveFluxCalculatorBase::_perm_derivs
const bool _perm_derivs
Flag to check whether permeabiity derivatives are non-zero.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:161
PorousFlowAdvectiveFluxCalculatorBase::execute
virtual void execute() override
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:193
PorousFlowAdvectiveFluxCalculatorBase::PorousFlowAdvectiveFluxCalculatorBase
PorousFlowAdvectiveFluxCalculatorBase(const InputParameters &parameters)
Definition: PorousFlowAdvectiveFluxCalculatorBase.C:45
PorousFlowAdvectiveFluxCalculatorBase::_permeability
const MaterialProperty< RealTensorValue > & _permeability
Permeability of porous material.
Definition: PorousFlowAdvectiveFluxCalculatorBase.h:90