12 #include "libmesh/string_to_enum.h"
13 #include "libmesh/parallel_sync.h"
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>(
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>(
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");
46 const InputParameters & parameters)
49 _num_vars(_dictator.numVariables()),
50 _gravity(getParam<RealVectorValue>(
"gravity")),
51 _phase(getParam<unsigned int>(
"phase")),
52 _permeability(getMaterialProperty<RealTensorValue>(
"PorousFlow_permeability_qp")),
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)),
72 _du_dvar_computed_by_thread(),
75 _triples_to_receive(),
77 _perm_derivs(_dictator.usePermDerivs())
81 "Phase number entered is greater than the number of phases specified in the "
82 "Dictator. Remember that indexing starts at 0");
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");
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.");
110 dof_id_type global_i, dof_id_type global_j,
unsigned local_i,
unsigned local_j,
unsigned qp)
117 for (
unsigned local_k = 0; local_k < _current_elem->n_nodes(); ++local_k)
119 const dof_id_type global_k = _current_elem->node_id(local_k);
120 for (
unsigned pvar = 0; pvar <
_num_vars; ++pvar)
122 RealVectorValue deriv =
132 for (
unsigned i = 0; i < LIBMESH_DIM; ++i)
137 _dkij_dvar[sequential_i][j][global_k][pvar] +=
138 _JxW[qp] * _coord[qp] * (-
_grad_phi[local_i][qp] * deriv *
_phi[local_j][qp]);
146 const bool resizing_was_needed =
155 if (resizing_was_needed)
160 _dflux_out_dvars.assign(num_nodes, std::map<dof_id_type, std::vector<Real>>());
162 for (dof_id_type sequential_i = 0; sequential_i < num_nodes; ++sequential_i)
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)
181 for (dof_id_type sequential_i = 0; sequential_i < num_nodes; ++sequential_i)
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)
200 for (
unsigned i = 0; i < _current_elem->n_nodes(); ++i)
202 const dof_id_type global_i = _current_elem->node_id(i);
206 for (
unsigned pvar = 0; pvar <
_num_vars; ++pvar)
217 static_cast<const PorousFlowAdvectiveFluxCalculatorBase &>(uo);
220 for (dof_id_type sequential_i = 0; sequential_i < num_nodes; ++sequential_i)
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];
232 for (dof_id_type sequential_i = 0; sequential_i <
_number_of_nodes; ++sequential_i)
238 const std::map<dof_id_type, std::vector<Real>> &
246 const std::map<dof_id_type, std::vector<Real>> &
264 const std::map<dof_id_type, Real> & dflux_out_du =
266 for (
const auto & node_du : dflux_out_du)
268 const dof_id_type j = node_du.first;
269 const Real dflux_out_du_j = node_du.second;
271 for (
unsigned pvar = 0; pvar <
_num_vars; ++pvar)
278 const std::vector<std::vector<Real>> & dflux_out_dKjk =
281 for (std::size_t index_j = 0; index_j < con_i.size(); ++index_j)
283 const dof_id_type node_j = con_i[index_j];
285 for (std::size_t index_k = 0; index_k < con_j.size(); ++index_k)
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)
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];
311 std::map<processor_id_type,
312 std::map<std::pair<dof_id_type, dof_id_type>, std::vector<dof_id_type>>>
314 for (
const auto & elem : _fe_problem.getEvaluableElementRange())
315 if (this->hasBlocks(elem->subdomain_id()))
317 const processor_id_type elem_pid = elem->processor_id();
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)
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>();
330 for (
const auto & global_neighbor_to_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);
343 for (
const auto & kv : tpr_global)
345 const processor_id_type pid = kv.first;
347 for (
const auto & pr_vec : kv.second)
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)
361 auto triples_action_functor = [
this](processor_id_type pid,
362 const std::vector<dof_id_type> & tts) {
365 Parallel::push_parallel_vector_data(this->comm(),
_triples_to_receive, triples_action_functor);
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)
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)
403 std::map<processor_id_type, std::vector<std::vector<Real>>> du_dvar_to_send;
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]);
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();
419 <<
", in du_dvar communication is incompatible with nodes_to_receive, which has size "
421 for (
unsigned i = 0; i < msg_size; ++i)
424 Parallel::push_parallel_vector_data(this->comm(), du_dvar_to_send, du_action_functor);
427 std::map<processor_id_type, std::vector<std::vector<Real>>> dkij_dvar_to_send;
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)
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]);
442 auto dk_action_functor = [
this](processor_id_type pid,
443 const std::vector<std::vector<Real>> & dkij_dvar_received) {
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 "
450 for (std::size_t i = 0; i < num; i += 3)
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];
459 Parallel::push_parallel_vector_data(this->comm(), dkij_dvar_to_send, dk_action_functor);