www.mooseframework.org
AdvectiveFluxCalculatorBase.h
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 
10 #pragma once
11 
12 #include "ElementUserObject.h"
14 
16 
17 template <>
19 
33 class AdvectiveFluxCalculatorBase : public ElementUserObject
34 {
35 public:
36  AdvectiveFluxCalculatorBase(const InputParameters & parameters);
37 
38  virtual void timestepSetup() override;
39 
40  virtual void meshChanged() override;
41 
42  virtual void initialize() override;
43 
44  virtual void threadJoin(const UserObject & uo) override;
45 
46  virtual void finalize() override;
47 
48  virtual void execute() override;
49 
60  virtual void executeOnElement(
61  dof_id_type global_i, dof_id_type global_j, unsigned local_i, unsigned local_j, unsigned qp);
62 
68  Real getFluxOut(dof_id_type node_i) const;
69 
76  const std::map<dof_id_type, Real> & getdFluxOutdu(dof_id_type node_i) const;
77 
84  const std::vector<std::vector<Real>> & getdFluxOutdKjk(dof_id_type node_i) const;
85 
94  unsigned getValence(dof_id_type node_i) const;
95 
96 protected:
105  virtual void buildCommLists();
106 
111  virtual void exchangeGhostedInfo();
112 
121  virtual Real computeVelocity(unsigned i, unsigned j, unsigned qp) const = 0;
122 
127  virtual Real computeU(unsigned i) const = 0;
128 
131 
139  void limitFlux(Real a, Real b, Real & limited, Real & dlimited_db) const;
140 
149  Real rPlus(dof_id_type sequential_i,
150  std::vector<Real> & dlimited_du,
151  std::vector<Real> & dlimited_dk) const;
152 
161  Real rMinus(dof_id_type sequential_i,
162  std::vector<Real> & dlimited_du,
163  std::vector<Real> & dlimited_dk) const;
164 
170 
176  std::vector<std::vector<Real>> _kij;
177 
179  std::vector<Real> _flux_out;
180 
183  std::vector<std::map<dof_id_type, Real>> _dflux_out_du;
184 
190  std::vector<std::vector<std::vector<Real>>> _dflux_out_dKjk;
191 
195  std::vector<unsigned> _valence;
196 
198  std::vector<Real> _u_nodal;
199 
201  std::vector<bool> _u_nodal_computed_by_thread;
202 
205 
207  std::size_t _number_of_nodes;
208 
210  processor_id_type _my_pid;
211 
218  std::map<processor_id_type, std::vector<dof_id_type>> _nodes_to_receive;
219 
226  std::map<processor_id_type, std::vector<dof_id_type>> _nodes_to_send;
227 
235  std::map<processor_id_type, std::vector<std::pair<dof_id_type, dof_id_type>>> _pairs_to_receive;
236 
244  std::map<processor_id_type, std::vector<std::pair<dof_id_type, dof_id_type>>> _pairs_to_send;
245 
253 
255  enum class PQPlusMinusEnum
256  {
257  PPlus,
258  PMinus,
259  QPlus,
260  QMinus
261  };
262 
274  Real PQPlusMinus(dof_id_type sequential_i,
275  const PQPlusMinusEnum pq_plus_minus,
276  std::vector<Real> & derivs,
277  std::vector<Real> & dpq_dk) const;
278 
285  void zeroedConnection(std::map<dof_id_type, Real> & the_map, dof_id_type node_i) const;
286 
288  std::vector<std::vector<Real>> _dij;
290  std::vector<std::vector<Real>> _dDij_dKij;
292  std::vector<std::vector<Real>> _dDij_dKji;
294  std::vector<std::vector<Real>> _dDii_dKij;
296  std::vector<std::vector<Real>> _dDii_dKji;
297 
298  std::vector<std::vector<Real>> _lij;
299  std::vector<Real> _rP;
300  std::vector<Real> _rM;
301 
303  std::vector<std::vector<Real>> _drP;
305  std::vector<std::vector<Real>> _drM;
307  std::vector<std::vector<Real>> _drP_dk;
309  std::vector<std::vector<Real>> _drM_dk;
310 
312  std::vector<std::vector<Real>> _fa;
316  std::vector<std::vector<std::map<dof_id_type, Real>>> _dfa;
321  std::vector<std::vector<std::vector<Real>>> _dFij_dKik;
326  std::vector<std::vector<std::vector<Real>>> _dFij_dKjk;
327 };
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::execute
virtual void execute() override
Definition: AdvectiveFluxCalculatorBase.C:218
PorousFlowConnectedNodes
Class designed to hold node ID information and information about nodal connectivity.
Definition: PorousFlowConnectedNodes.h:44
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
validParams< AdvectiveFluxCalculatorBase >
InputParameters validParams< AdvectiveFluxCalculatorBase >()
Definition: AdvectiveFluxCalculatorBase.C:19
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum::MinMod
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
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
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
AdvectiveFluxCalculatorBase::FluxLimiterTypeEnum::MC
AdvectiveFluxCalculatorBase::buildCommLists
virtual void buildCommLists()
When using multiple processors, other processors will compute:
Definition: AdvectiveFluxCalculatorBase.C:850
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
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
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.h
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
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
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
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
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
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