https://mooseframework.inl.gov
WCNSLinearFVScalarTransportPhysics.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "WCNSFVFlowPhysicsBase.h"
12 #include "NS.h"
13 
16 
19 {
21  params.addClassDescription("Define the Navier Stokes weakly-compressible scalar field transport "
22  "equation(s) using the linear finite volume discretization");
23  params.addParam<bool>("use_nonorthogonal_correction",
24  true,
25  "If the nonorthogonal correction should be used when computing the normal "
26  "gradient, notably in the diffusion term.");
27 
28  // Not supported
29  params.suppressParameter<MooseEnum>("preconditioning");
30 
31  return params;
32 }
33 
35  const InputParameters & parameters)
37 {
39  _flow_equations_physics->paramError("porous_medium_treatment",
40  "Porous media scalar advection is currently unimplemented");
41 }
42 
43 void
45 {
46  // For compatibility with Modules/NavierStokesFV syntax
48  return;
49 
50  auto params = getFactory().getValidParams("MooseLinearVariableFVReal");
51  assignBlocks(params, _blocks);
52 
53  for (const auto name_i : index_range(_passive_scalar_names))
54  {
55  // Dont add if the user already defined the variable
56  if (!shouldCreateVariable(_passive_scalar_names[name_i], _blocks, /*error if aux*/ true))
57  {
58  reportPotentiallyMissedParameters({"system_names", "passive_scalar_scaling"},
59  "MooseLinearVariableFVReal");
60  continue;
61  }
62 
63  params.set<SolverSystemName>("solver_sys") = getSolverSystem(name_i);
64  if (isParamValid("passive_scalar_scaling"))
65  params.set<std::vector<Real>>("scaling") = {
66  getParam<std::vector<Real>>("passive_scalar_scaling")[name_i]};
67 
68  getProblem().addVariable("MooseLinearVariableFVReal", _passive_scalar_names[name_i], params);
69  }
70 }
71 
72 void
74 {
75  std::string kernel_type = "LinearFVTimeDerivative";
76  InputParameters params = getFactory().getValidParams(kernel_type);
77  assignBlocks(params, _blocks);
78 
79  for (const auto & vname : _passive_scalar_names)
80  {
81  params.set<LinearVariableName>("variable") = vname;
82  if (shouldCreateTimeDerivative(vname, _blocks, /*error if already defined */ false))
83  getProblem().addLinearFVKernel(kernel_type, prefix() + "ins_" + vname + "_time", params);
84  }
85 }
86 
87 void
89 {
90  const std::string kernel_type = "LinearFVScalarAdvection";
91  InputParameters params = getFactory().getValidParams(kernel_type);
92 
93  assignBlocks(params, _blocks);
94  params.set<UserObjectName>("rhie_chow_user_object") = _flow_equations_physics->rhieChowUOName();
95  params.set<MooseEnum>("advected_interp_method") =
96  getParam<MooseEnum>("passive_scalar_advection_interpolation");
97  setSlipVelocityParams(params);
98 
99  for (const auto & vname : _passive_scalar_names)
100  {
101  params.set<LinearVariableName>("variable") = vname;
102  getProblem().addLinearFVKernel(kernel_type, prefix() + "ins_" + vname + "_advection", params);
103  }
104 }
105 
106 void
108 {
109  // Direct specification of diffusion term
110  const auto passive_scalar_diffusivities =
111  getParam<std::vector<MooseFunctorName>>("passive_scalar_diffusivity");
112 
113  if (passive_scalar_diffusivities.size())
114  {
115  const std::string kernel_type = "LinearFVDiffusion";
116  InputParameters params = getFactory().getValidParams(kernel_type);
117  assignBlocks(params, _blocks);
118  params.set<bool>("use_nonorthogonal_correction") =
119  getParam<bool>("use_nonorthogonal_correction");
120  for (const auto name_i : index_range(_passive_scalar_names))
121  {
122  params.set<LinearVariableName>("variable") = _passive_scalar_names[name_i];
123  params.set<MooseFunctorName>("diffusion_coeff") = passive_scalar_diffusivities[name_i];
125  kernel_type, prefix() + "ins_" + _passive_scalar_names[name_i] + "_diffusion", params);
126  }
127  }
128 }
129 
130 void
132 {
133  const std::string kernel_type = "LinearFVSource";
134  InputParameters params = getFactory().getValidParams(kernel_type);
135  assignBlocks(params, _blocks);
136 
137  for (const auto scalar_i : index_range(_passive_scalar_names))
138  {
139  params.set<LinearVariableName>("variable") = _passive_scalar_names[scalar_i];
140 
141  if (_passive_scalar_sources.size())
142  {
143  // Added for backward compatibility with former Modules/NavierStokesFV syntax
144  params.set<MooseFunctorName>("source_density") = _passive_scalar_sources[scalar_i];
146  kernel_type, prefix() + "ins_" + _passive_scalar_names[scalar_i] + "_source", params);
147  }
148 
150  for (const auto i : index_range(_passive_scalar_coupled_sources[scalar_i]))
151  {
152  params.set<MooseFunctorName>("source_density") =
153  _passive_scalar_coupled_sources[scalar_i][i];
154  if (_passive_scalar_sources_coef.size())
155  params.set<Real>("scaling_factor") = _passive_scalar_sources_coef[scalar_i][i];
156 
157  getProblem().addLinearFVKernel(kernel_type,
158  prefix() + "ins_" + _passive_scalar_names[scalar_i] +
159  "_coupled_source_" + std::to_string(i),
160  params);
161  }
162  }
163 }
164 
165 void
167 {
168  const auto & inlet_boundaries = _flow_equations_physics->getInletBoundaries();
169  if (inlet_boundaries.empty())
170  return;
171 
172  // Boundary checks
173  // TODO: once we have vectors of MooseEnum, we could use the same templated check for types and
174  // functors
175  if (inlet_boundaries.size() * _passive_scalar_names.size() != _passive_scalar_inlet_types.size())
176  paramError(
177  "passive_scalar_inlet_types",
178  "The number of scalar inlet types (" + std::to_string(_passive_scalar_inlet_types.size()) +
179  ") is not equal to the number of inlet boundaries (" +
180  std::to_string(inlet_boundaries.size()) + ") times the number of passive scalars (" +
181  std::to_string(_passive_scalar_names.size()) + ")");
183  paramError("passive_scalar_inlet_functors",
184  "The number of groups of inlet functors (" +
185  std::to_string(_passive_scalar_inlet_functors.size()) +
186  ") is not equal to the number of passive scalars (" +
187  std::to_string(_passive_scalar_names.size()) + ")");
188 
189  for (const auto name_i : index_range(_passive_scalar_names))
190  {
191  if (inlet_boundaries.size() != _passive_scalar_inlet_functors[name_i].size())
192  paramError("passive_scalar_inlet_functors",
193  "The number of inlet boundary functors for scalar '" +
194  _passive_scalar_names[name_i] +
195  "' does not match the number of inlet boundaries (" +
196  std::to_string(_passive_scalar_inlet_functors[name_i].size()) + ")");
197 
198  unsigned int num_inlets = inlet_boundaries.size();
199  for (unsigned int bc_ind = 0; bc_ind < num_inlets; ++bc_ind)
200  {
201  if (_passive_scalar_inlet_types[name_i * num_inlets + bc_ind] == "fixed-value")
202  {
203  const std::string bc_type = "LinearFVAdvectionDiffusionFunctorDirichletBC";
204  InputParameters params = getFactory().getValidParams(bc_type);
205  params.set<LinearVariableName>("variable") = _passive_scalar_names[name_i];
206  params.set<MooseFunctorName>("functor") = _passive_scalar_inlet_functors[name_i][bc_ind];
207  params.set<std::vector<BoundaryName>>("boundary") = {inlet_boundaries[bc_ind]};
208 
210  bc_type, _passive_scalar_names[name_i] + "_" + inlet_boundaries[bc_ind], params);
211  }
212  else if (_passive_scalar_inlet_types[name_i * num_inlets + bc_ind] == "flux-mass" ||
213  _passive_scalar_inlet_types[name_i * num_inlets + bc_ind] == "flux-velocity")
214  {
215  mooseError("Flux boundary conditions not supported at this time using the linear finite "
216  "volume discretization");
217  }
218  }
219  }
220 }
221 
222 void
224 {
225  const auto & outlet_boundaries = _flow_equations_physics->getOutletBoundaries();
226  if (outlet_boundaries.empty())
227  return;
228 
229  for (const auto & outlet_bdy : outlet_boundaries)
230  {
231  const std::string bc_type = "LinearFVAdvectionDiffusionOutflowBC";
232  InputParameters params = getFactory().getValidParams(bc_type);
233  params.set<std::vector<BoundaryName>>("boundary") = {outlet_bdy};
234  params.set<bool>("use_two_term_expansion") =
235  getParam<bool>("passive_scalar_two_term_bc_expansion");
236 
237  for (const auto name_i : index_range(_passive_scalar_names))
238  {
239  params.set<LinearVariableName>("variable") = _passive_scalar_names[name_i];
240  getProblem().addLinearFVBC(bc_type, _passive_scalar_names[name_i] + "_" + outlet_bdy, params);
241  }
242  }
243 }
std::string prefix() const
void assignBlocks(InputParameters &params, const std::vector< SubdomainName > &blocks) const
bool shouldCreateVariable(const VariableName &var_name, const std::vector< SubdomainName > &blocks, const bool error_if_aux)
Factory & getFactory()
const std::vector< BoundaryName > & getOutletBoundaries() const
Get the outlet boundaries.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const std::vector< BoundaryName > & getInletBoundaries() const
Get the inlet boundaries.
T & set(const std::string &name, bool quiet_mode=false)
registerNavierStokesPhysicsBaseTasks("NavierStokesApp", WCNSLinearFVScalarTransportPhysics)
InputParameters getValidParams(const std::string &name) const
unsigned int size() const
Creates all the objects needed to solve the Navier Stokes scalar transport equations using the linear...
std::vector< NonlinearVariableName > _passive_scalar_names
Names of the passive scalar variables.
bool shouldCreateTimeDerivative(const VariableName &var_name, const std::vector< SubdomainName > &blocks, const bool error_if_already_defined) const
std::vector< SubdomainName > _blocks
void suppressParameter(const std::string &name)
bool isParamValid(const std::string &name) const
std::vector< std::vector< MooseFunctorName > > _passive_scalar_inlet_functors
Functors describing the inlet boundary values. See passive_scalar_inlet_types for what the functors a...
virtual FEProblemBase & getProblem()
virtual UserObjectName rhieChowUOName() const =0
Return the name of the Rhie Chow user object.
virtual void addScalarInletBC() override
Functions adding boundary conditions for the incompressible simulation.
const SolverSystemName & getSolverSystem(unsigned int variable_index) const
std::vector< std::vector< Real > > _passive_scalar_sources_coef
Coefficients multiplying for the passive scalar sources. Inner indexing is scalar variable index...
virtual void addScalarTimeKernels() override
Functions adding kernels for the incompressible / weakly-compressible scalar transport equation If th...
const bool _porous_medium_treatment
Switch to show if porous medium treatment is requested or not.
std::vector< std::vector< MooseFunctorName > > _passive_scalar_coupled_sources
Functors for the passive scalar (coupled) sources. Inner indexing is scalar variable index...
void paramError(const std::string &param, Args... args) const
virtual void addVariable(const std::string &var_type, const std::string &var_name, InputParameters &params)
virtual void setSlipVelocityParams(InputParameters &) const
virtual void addLinearFVKernel(const std::string &kernel_name, const std::string &name, InputParameters &parameters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void addLinearFVBC(const std::string &fv_bc_name, const std::string &name, InputParameters &parameters)
Creates all the objects needed to solve the Navier Stokes scalar transport equations.
void mooseError(Args &&... args) const
std::vector< MooseFunctorName > _passive_scalar_sources
Functors for the passive scalar sources. Indexing is scalar variable index.
void addClassDescription(const std::string &doc_string)
void reportPotentiallyMissedParameters(const std::vector< std::string > &param_names, const std::string &object_type) const
WCNSLinearFVScalarTransportPhysics(const InputParameters &parameters)
MultiMooseEnum _passive_scalar_inlet_types
Passive scalar inlet boundary types.
const bool _has_scalar_equation
A boolean to help compatibility with the old Modules/NavierStokesFV syntax or to deliberately skip ad...
virtual void addFVKernel(const std::string &kernel_name, const std::string &name, InputParameters &parameters)
const WCNSFVFlowPhysicsBase * _flow_equations_physics
Flow physics.
auto index_range(const T &sizable)
registerWCNSFVScalarTransportBaseTasks("NavierStokesApp", WCNSLinearFVScalarTransportPhysics)
virtual void addScalarSourceKernels() override
Equivalent of NSFVAction addScalarCoupledSourceKernels.