https://mooseframework.inl.gov
RhieChowMassFlux.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 
10 // MOOSE includes
11 #include "RhieChowMassFlux.h"
12 #include "SubProblem.h"
13 #include "MooseMesh.h"
14 #include "NS.h"
15 #include "VectorCompositeFunctor.h"
16 #include "PIMPLE.h"
17 #include "SIMPLE.h"
18 #include "PetscVectorReader.h"
19 #include "LinearSystem.h"
21 
22 // libMesh includes
23 #include "libmesh/mesh_base.h"
24 #include "libmesh/elem_range.h"
25 #include "libmesh/petsc_matrix.h"
26 
27 using namespace libMesh;
28 
29 registerMooseObject("NavierStokesApp", RhieChowMassFlux);
30 
33 {
36 
37  params.addClassDescription("Computes H/A and 1/A together with face mass fluxes for segregated "
38  "momentum-pressure equations using linear systems.");
39 
40  params.addRequiredParam<VariableName>(NS::pressure, "The pressure variable.");
41  params.addRequiredParam<VariableName>("u", "The x-component of velocity");
42  params.addParam<VariableName>("v", "The y-component of velocity");
43  params.addParam<VariableName>("w", "The z-component of velocity");
44  params.addRequiredParam<std::string>("p_diffusion_kernel",
45  "The diffusion kernel acting on the pressure.");
46 
47  params.addRequiredParam<MooseFunctorName>(NS::density, "Density functor");
48 
49  // We disable the execution of this, should only provide functions
50  // for the SIMPLE executioner
51  ExecFlagEnum & exec_enum = params.set<ExecFlagEnum>("execute_on", true);
52  exec_enum.addAvailableFlags(EXEC_NONE);
53  exec_enum = {EXEC_NONE};
54  params.suppressParameter<ExecFlagEnum>("execute_on");
55 
56  // Pressure projection
57  params.addParam<MooseEnum>("pressure_projection_method",
58  MooseEnum("standard consistent", "standard"),
59  "The method to use in the pressure projection for Ainv - "
60  "standard (SIMPLE) or consistent (SIMPLEC)");
61 
62  return params;
63 }
64 
66  : RhieChowFaceFluxProvider(params),
68  _moose_mesh(UserObject::_subproblem.mesh()),
69  _mesh(_moose_mesh.getMesh()),
70  _dim(blocksMaxDimension()),
71  _p(dynamic_cast<MooseLinearVariableFVReal *>(
72  &UserObject::_subproblem.getVariable(0, getParam<VariableName>(NS::pressure)))),
73  _vel(_dim, nullptr),
74  _HbyA_flux(_moose_mesh, blockIDs(), "HbyA_flux"),
75  _Ainv(_moose_mesh, blockIDs(), "Ainv"),
76  _face_mass_flux(
77  declareRestartableData<FaceCenteredMapFunctor<Real, std::unordered_map<dof_id_type, Real>>>(
78  "face_flux", _moose_mesh, blockIDs(), "face_values")),
79  _rho(getFunctor<Real>(NS::density)),
80  _pressure_projection_method(getParam<MooseEnum>("pressure_projection_method"))
81 {
82  if (!_p)
83  paramError(NS::pressure, "the pressure must be a MooseLinearVariableFVReal.");
84  checkBlocks(*_p);
85 
86  std::vector<std::string> vel_names = {"u", "v", "w"};
87  for (const auto i : index_range(_vel))
88  {
89  _vel[i] = dynamic_cast<MooseLinearVariableFVReal *>(
90  &UserObject::_subproblem.getVariable(0, getParam<VariableName>(vel_names[i])));
91 
92  if (!_vel[i])
93  paramError(vel_names[i], "the velocity must be a MOOSELinearVariableFVReal.");
94  checkBlocks(*_vel[i]);
95  }
96 
97  // Register the elemental/face functors which will be queried in the pressure equation
98  for (const auto tid : make_range(libMesh::n_threads()))
99  {
102  }
103 
104  if (!dynamic_cast<SIMPLE *>(getMooseApp().getExecutioner()) &&
105  !dynamic_cast<PIMPLE *>(getMooseApp().getExecutioner()))
106  mooseError(this->name(),
107  " should only be used with a linear segregated thermal-hydraulics solver!");
108 }
109 
110 void
112  const std::vector<LinearSystem *> & momentum_systems,
113  const LinearSystem & pressure_system,
114  const std::vector<unsigned int> & momentum_system_numbers)
115 {
116  _momentum_systems = momentum_systems;
117  _momentum_system_numbers = momentum_system_numbers;
118  _pressure_system = &pressure_system;
120 
122  for (auto & system : _momentum_systems)
123  {
124  _global_momentum_system_numbers.push_back(system->number());
125  _momentum_implicit_systems.push_back(dynamic_cast<LinearImplicitSystem *>(&system->system()));
126  }
127 
129 }
130 
131 void
133 {
134  _HbyA_flux.clear();
135  _Ainv.clear();
136  _face_mass_flux.clear();
138 }
139 
140 void
142 {
143  // We fetch the pressure diffusion kernel to ensure that the face flux correction
144  // is consistent with the pressure discretization in the Poisson equation.
145  std::vector<LinearFVFluxKernel *> flux_kernel;
146  auto base_query = _fe_problem.theWarehouse()
147  .query()
148  .template condition<AttribThread>(_tid)
149  .template condition<AttribSysNum>(_p->sys().number())
150  .template condition<AttribSystem>("LinearFVFluxKernel")
151  .template condition<AttribName>(getParam<std::string>("p_diffusion_kernel"))
152  .queryInto(flux_kernel);
153  if (flux_kernel.size() != 1)
154  paramError(
155  "p_diffusion_kernel",
156  "The kernel with the given name could not be found or multiple instances were identified.");
157  _p_diffusion_kernel = dynamic_cast<LinearFVAnisotropicDiffusion *>(flux_kernel[0]);
158  if (!_p_diffusion_kernel)
159  paramError("p_diffusion_kernel",
160  "The provided diffusion kernel should of type LinearFVAnisotropicDiffusion!");
161 }
162 
163 void
165 {
166  // We cache the cell volumes into a petsc vector for corrections here so we can use
167  // the optimized petsc operations for the normalization
169  for (const auto & elem_info : _fe_problem.mesh().elemInfoVector())
170  // We have to check this because the variable might not be defined on the given
171  // block
172  if (hasBlocks(elem_info->subdomain_id()))
173  {
174  const auto elem_dof = elem_info->dofIndices()[_global_pressure_system_number][0];
175  _cell_volumes->set(elem_dof, elem_info->volume() * elem_info->coordFactor());
176  }
177 
178  _cell_volumes->close();
179 
180  _flow_face_info.clear();
181  for (auto & fi : _fe_problem.mesh().faceInfo())
182  if (hasBlocks(fi->elemPtr()->subdomain_id()) ||
183  (fi->neighborPtr() && hasBlocks(fi->neighborPtr()->subdomain_id())))
184  _flow_face_info.push_back(fi);
185 }
186 
187 void
189 {
190  for (const auto & pair : _HbyA_flux)
191  _HbyA_flux[pair.first] = 0;
192 
193  for (const auto & pair : _Ainv)
194  _Ainv[pair.first] = 0;
195 }
196 
197 void
199 {
200  using namespace Moose::FV;
201 
202  const auto time_arg = Moose::currentState();
203 
204  // We loop through the faces and compute the resulting face fluxes from the
205  // initial conditions for velocity
206  for (auto & fi : _flow_face_info)
207  {
208  RealVectorValue density_times_velocity;
209 
210  // On internal face we do a regular interpolation with geometric weights
211  if (_vel[0]->isInternalFace(*fi))
212  {
213  const auto & elem_info = *fi->elemInfo();
214  const auto & neighbor_info = *fi->neighborInfo();
215 
216  Real elem_rho = _rho(makeElemArg(fi->elemPtr()), time_arg);
217  Real neighbor_rho = _rho(makeElemArg(fi->neighborPtr()), time_arg);
218 
219  for (const auto dim_i : index_range(_vel))
220  interpolate(InterpMethod::Average,
221  density_times_velocity(dim_i),
222  _vel[dim_i]->getElemValue(elem_info, time_arg) * elem_rho,
223  _vel[dim_i]->getElemValue(neighbor_info, time_arg) * neighbor_rho,
224  *fi,
225  true);
226  }
227  // On the boundary, we just take the boundary values
228  else
229  {
230  const Elem * const boundary_elem =
231  hasBlocks(fi->elemPtr()->subdomain_id()) ? fi->elemPtr() : fi->neighborPtr();
232 
233  const Moose::FaceArg boundary_face{
234  fi, Moose::FV::LimiterType::CentralDifference, true, false, boundary_elem, nullptr};
235 
236  const Real face_rho = _rho(boundary_face, time_arg);
237  for (const auto dim_i : index_range(_vel))
238  density_times_velocity(dim_i) =
239  face_rho * raw_value((*_vel[dim_i])(boundary_face, time_arg));
240  }
241 
242  _face_mass_flux[fi->id()] = density_times_velocity * fi->normal();
243  }
244 }
245 
246 Real
248 {
249  return _face_mass_flux.evaluate(&fi);
250 }
251 
252 Real
254 {
255  const Moose::FaceArg face_arg{&fi,
257  /*elem_is_upwind=*/true,
258  /*correct_skewness=*/false,
259  &fi.elem(),
260  /*state_limiter*/ nullptr};
261  const Real face_rho = _rho(face_arg, Moose::currentState());
262  return libmesh_map_find(_face_mass_flux, fi.id()) / face_rho;
263 }
264 
265 Real
267  const FaceInfo & fi,
268  const Moose::StateArg & time,
269  const THREAD_ID /*tid*/,
270  bool libmesh_dbg_var(subtract_mesh_velocity)) const
271 {
272  mooseAssert(!subtract_mesh_velocity, "RhieChowMassFlux does not support moving meshes yet!");
273 
275  mooseError("Interpolation methods other than Rhie-Chow are not supported!");
276  if (time.state != Moose::currentState().state)
277  mooseError("Older interpolation times are not supported!");
278 
279  return getVolumetricFaceFlux(fi);
280 }
281 
282 void
284 {
285  using namespace Moose::FV;
286 
287  const auto time_arg = Moose::currentState();
288 
289  // Petsc vector reader to make the repeated reading from the vector faster
291 
292  // We loop through the faces and compute the face fluxes using the pressure gradient
293  // and the momentum matrix/right hand side
294  for (auto & fi : _flow_face_info)
295  {
296  // Making sure the kernel knows which face we are on
298 
299  // We are setting this to 1.0 because we don't want to multiply the kernel contributions
300  // with the surface area yet. The surface area will be factored in in the advection kernels.
302 
303  Real p_grad_flux = 0.0;
304  if (_p->isInternalFace(*fi))
305  {
306  const auto & elem_info = *fi->elemInfo();
307  const auto & neighbor_info = *fi->neighborInfo();
308 
309  // Fetching the dof indices for the pressure variable
310  const auto elem_dof = elem_info.dofIndices()[_global_pressure_system_number][0];
311  const auto neighbor_dof = neighbor_info.dofIndices()[_global_pressure_system_number][0];
312 
313  // Fetching the values of the pressure for the element and the neighbor
314  const auto p_elem_value = p_reader(elem_dof);
315  const auto p_neighbor_value = p_reader(neighbor_dof);
316 
317  // Compute the elem matrix contributions for the face
318  const auto elem_matrix_contribution = _p_diffusion_kernel->computeElemMatrixContribution();
319  const auto neighbor_matrix_contribution =
321  const auto elem_rhs_contribution =
323 
324  // Compute the face flux from the matrix and right hand side contributions
325  p_grad_flux = (p_neighbor_value * neighbor_matrix_contribution +
326  p_elem_value * elem_matrix_contribution) -
327  elem_rhs_contribution;
328  }
329  else if (auto * bc_pointer = _p->getBoundaryCondition(*fi->boundaryIDs().begin()))
330  {
331  mooseAssert(fi->boundaryIDs().size() == 1, "We should only have one boundary on every face.");
332 
333  bc_pointer->setupFaceData(
334  fi, fi->faceType(std::make_pair(_p->number(), _global_pressure_system_number)));
335 
336  const ElemInfo & elem_info =
337  hasBlocks(fi->elemPtr()->subdomain_id()) ? *fi->elemInfo() : *fi->neighborInfo();
338  const auto p_elem_value = _p->getElemValue(elem_info, time_arg);
339  const auto matrix_contribution =
341  const auto rhs_contribution =
343 
344  // On the boundary, only the element side has a contribution
345  p_grad_flux = (p_elem_value * matrix_contribution - rhs_contribution);
346  }
347  // Compute the new face flux
348  _face_mass_flux[fi->id()] = -_HbyA_flux[fi->id()] + p_grad_flux;
349  }
350 }
351 
352 void
354 {
355  auto & pressure_gradient = _pressure_system->gradientContainer();
356 
357  // We set the dof value in the solution vector the same logic applies:
358  // u_C = -(H/A)_C - (1/A)_C*grad(p)_C where C is the cell index
359  for (const auto system_i : index_range(_momentum_implicit_systems))
360  {
361  auto working_vector = _Ainv_raw[system_i]->clone();
362  working_vector->pointwise_mult(*working_vector, *pressure_gradient[system_i]);
363  working_vector->add(*_HbyA_raw[system_i]);
364  working_vector->scale(-1.0);
365  (*_momentum_implicit_systems[system_i]->solution) = *working_vector;
366  _momentum_implicit_systems[system_i]->update();
367  _momentum_systems[system_i]->setSolution(
368  *_momentum_implicit_systems[system_i]->current_local_solution);
369  }
370 }
371 
372 void
374 {
375  // We loop through the faces and populate the coupling fields (face H/A and 1/H)
376  // with 0s for now. Pressure corrector solves will always come after the
377  // momentum source so we expect these fields to change before the actual solve.
378  for (auto & fi : _fe_problem.mesh().faceInfo())
379  {
380  _Ainv[fi->id()];
381  _HbyA_flux[fi->id()];
382  }
383 }
384 
385 void
387  const std::vector<std::unique_ptr<NumericVector<Number>>> & raw_hbya,
388  const std::vector<std::unique_ptr<NumericVector<Number>>> & raw_Ainv)
389 {
390  // We have the raw H/A and 1/A vectors in a petsc format. This function
391  // will create face functors from them
392  using namespace Moose::FV;
393  const auto time_arg = Moose::currentState();
394 
395  // Create the petsc vector readers for faster repeated access
396  std::vector<PetscVectorReader> hbya_reader;
397  for (const auto dim_i : index_range(raw_hbya))
398  hbya_reader.emplace_back(*raw_hbya[dim_i]);
399 
400  std::vector<PetscVectorReader> ainv_reader;
401  for (const auto dim_i : index_range(raw_Ainv))
402  ainv_reader.emplace_back(*raw_Ainv[dim_i]);
403 
404  // We loop through the faces and populate the coupling fields (face H/A and 1/H)
405  for (auto & fi : _flow_face_info)
406  {
407  Real face_rho = 0;
408  RealVectorValue face_hbya;
409 
410  // We do the lookup in advance
411  auto & Ainv = _Ainv[fi->id()];
412 
413  // If it is internal, we just interpolate (using geometric weights) to the face
414  if (_vel[0]->isInternalFace(*fi))
415  {
416  // Get the dof indices for the element and the neighbor
417  const auto & elem_info = *fi->elemInfo();
418  const auto & neighbor_info = *fi->neighborInfo();
419  const auto elem_dof = elem_info.dofIndices()[_global_momentum_system_numbers[0]][0];
420  const auto neighbor_dof = neighbor_info.dofIndices()[_global_momentum_system_numbers[0]][0];
421 
422  // Get the density values for the element and neighbor. We need this multiplication to make
423  // the coupling fields mass fluxes.
424  const Real elem_rho = _rho(makeElemArg(fi->elemPtr()), time_arg);
425  const Real neighbor_rho = _rho(makeElemArg(fi->neighborPtr()), time_arg);
426 
427  // Now we do the interpolation to the face
428  interpolate(Moose::FV::InterpMethod::Average, face_rho, elem_rho, neighbor_rho, *fi, true);
429  for (const auto dim_i : index_range(raw_hbya))
430  {
432  face_hbya(dim_i),
433  hbya_reader[dim_i](elem_dof),
434  hbya_reader[dim_i](neighbor_dof),
435  *fi,
436  true);
437  interpolate(InterpMethod::Average,
438  Ainv(dim_i),
439  elem_rho * ainv_reader[dim_i](elem_dof),
440  neighbor_rho * ainv_reader[dim_i](neighbor_dof),
441  *fi,
442  true);
443  }
444  }
445  else
446  {
447  const ElemInfo & elem_info =
448  hasBlocks(fi->elemPtr()->subdomain_id()) ? *fi->elemInfo() : *fi->neighborInfo();
449  const auto elem_dof = elem_info.dofIndices()[_global_momentum_system_numbers[0]][0];
450 
451  // If it is a Dirichlet BC, we use the dirichlet value the make sure the face flux
452  // is consistent
453  if (_vel[0]->isDirichletBoundaryFace(*fi))
454  {
455  const Moose::FaceArg boundary_face{
456  fi, Moose::FV::LimiterType::CentralDifference, true, false, elem_info.elem(), nullptr};
457  face_rho = _rho(boundary_face, Moose::currentState());
458 
459  for (const auto dim_i : make_range(_dim))
460  face_hbya(dim_i) =
461  -MetaPhysicL::raw_value((*_vel[dim_i])(boundary_face, Moose::currentState()));
462  }
463  // Otherwise we just do a one-term expansion (so we just use the element value)
464  else
465  {
466  const auto elem_dof = elem_info.dofIndices()[_global_momentum_system_numbers[0]][0];
467 
468  face_rho = _rho(makeElemArg(elem_info.elem()), time_arg);
469  for (const auto dim_i : make_range(_dim))
470  face_hbya(dim_i) = hbya_reader[dim_i](elem_dof);
471  }
472 
473  // We just do a one-term expansion for 1/A no matter what
474  const Real elem_rho = _rho(makeElemArg(elem_info.elem()), time_arg);
475  for (const auto dim_i : index_range(raw_Ainv))
476  Ainv(dim_i) = elem_rho * ainv_reader[dim_i](elem_dof);
477  }
478  // Lastly, we populate the face flux resulted by H/A
479  _HbyA_flux[fi->id()] = face_hbya * fi->normal() * face_rho;
480  }
481 }
482 
483 void
484 RhieChowMassFlux::computeHbyA(const bool with_updated_pressure, bool verbose)
485 {
486  if (verbose)
487  {
488  _console << "************************************" << std::endl;
489  _console << "Computing HbyA" << std::endl;
490  _console << "************************************" << std::endl;
491  }
493  "The momentum system shall be linked before calling this function!");
494 
495  auto & pressure_gradient = selectPressureGradient(with_updated_pressure);
496 
497  _HbyA_raw.clear();
498  _Ainv_raw.clear();
499  for (auto system_i : index_range(_momentum_systems))
500  {
501  LinearImplicitSystem * momentum_system = _momentum_implicit_systems[system_i];
502 
503  NumericVector<Number> & rhs = *(momentum_system->rhs);
504  NumericVector<Number> & current_local_solution = *(momentum_system->current_local_solution);
505  NumericVector<Number> & solution = *(momentum_system->solution);
506  PetscMatrix<Number> * mmat = dynamic_cast<PetscMatrix<Number> *>(momentum_system->matrix);
507  mooseAssert(mmat,
508  "The matrices used in the segregated INSFVRhieChow objects need to be convertable "
509  "to PetscMatrix!");
510 
511  if (verbose)
512  {
513  _console << "Matrix in rc object" << std::endl;
514  mmat->print();
515  }
516 
517  // First, we extract the diagonal and we will hold on to it for a little while
518  _Ainv_raw.push_back(current_local_solution.zero_clone());
519  NumericVector<Number> & Ainv = *(_Ainv_raw.back());
520 
521  mmat->get_diagonal(Ainv);
522 
523  if (verbose)
524  {
525  _console << "Velocity solution in H(u)" << std::endl;
526  solution.print();
527  }
528 
529  // Time to create H(u) = M_{offdiag} * u - b_{nonpressure}
530  _HbyA_raw.push_back(current_local_solution.zero_clone());
531  NumericVector<Number> & HbyA = *(_HbyA_raw.back());
532 
533  // We start with the matrix product part, we will do
534  // M*u - A*u for 2 reasons:
535  // 1, We assume A*u petsc operation is faster than setting the matrix diagonal to 0
536  // 2, In PISO loops we need to reuse the matrix so we can't just set the diagonals to 0
537 
538  // We create a working vector to ease some of the operations, we initialize its values
539  // with the current solution values to have something for the A*u term
540  auto working_vector = momentum_system->current_local_solution->zero_clone();
541  PetscVector<Number> * working_vector_petsc =
542  dynamic_cast<PetscVector<Number> *>(working_vector.get());
543  mooseAssert(working_vector_petsc,
544  "The vectors used in the RhieChowMassFlux objects need to be convertable "
545  "to PetscVectors!");
546 
547  mmat->vector_mult(HbyA, solution);
548  working_vector_petsc->pointwise_mult(Ainv, solution);
549  HbyA.add(-1.0, *working_vector_petsc);
550 
551  if (verbose)
552  {
553  _console << " H(u)" << std::endl;
554  HbyA.print();
555  }
556 
557  // We continue by adding the momentum right hand side contributions
558  HbyA.add(-1.0, rhs);
559 
560  // Unfortunately, the pressure forces are included in the momentum RHS
561  // so we have to correct them back
562  working_vector_petsc->pointwise_mult(*pressure_gradient[system_i], *_cell_volumes);
563  HbyA.add(-1.0, *working_vector_petsc);
564 
565  if (verbose)
566  {
567  _console << "total RHS" << std::endl;
568  rhs.print();
569  _console << "pressure RHS" << std::endl;
570  pressure_gradient[system_i]->print();
571  _console << " H(u)-rhs-relaxation_source" << std::endl;
572  HbyA.print();
573  }
574 
575  // It is time to create element-wise 1/A-s based on the the diagonal of the momentum matrix
576  *working_vector_petsc = 1.0;
577  Ainv.pointwise_divide(*working_vector_petsc, Ainv);
578 
579  // Create 1/A*(H(u)-RHS)
580  HbyA.pointwise_mult(HbyA, Ainv);
581 
582  if (verbose)
583  {
584  _console << " (H(u)-rhs)/A" << std::endl;
585  HbyA.print();
586  }
587 
588  if (_pressure_projection_method == "consistent")
589  {
590 
591  // Consistent Corrections to SIMPLE
592  // 1. Ainv_old = 1/a_p <- Ainv = 1/(a_p + \sum_n a_n)
593  // 2. H(u) <- H(u*) + H(u') = H(u*) - (Ainv - Ainv_old) * grad(p) * Vc
594 
595  if (verbose)
596  _console << "Performing SIMPLEC projection." << std::endl;
597 
598  // Lambda function to calculate the sum of diagonal and neighbor coefficients
599  auto get_row_sum = [mmat](NumericVector<Number> & sum_vector)
600  {
601  // Ensure the sum_vector is zeroed out
602  sum_vector.zero();
603 
604  // Local row size
605  const auto local_size = mmat->local_m();
606 
607  for (const auto row_i : make_range(local_size))
608  {
609  // Get all non-zero components of the row of the matrix
610  const auto global_index = mmat->row_start() + row_i;
611  std::vector<numeric_index_type> indices;
612  std::vector<Real> values;
613  mmat->get_row(global_index, indices, values);
614 
615  // Sum row elements (no absolute values)
616  const Real row_sum = std::accumulate(values.cbegin(), values.cend(), 0.0);
617 
618  // Add the sum of diagonal and elements to the sum_vector
619  sum_vector.add(global_index, row_sum);
620  }
621  sum_vector.close();
622  };
623 
624  // Create a temporary vector to store the sum of diagonal and neighbor coefficients
625  auto row_sum = current_local_solution.zero_clone();
626  get_row_sum(*row_sum);
627 
628  // Create vector with new inverse projection matrix
629  auto Ainv_full = current_local_solution.zero_clone();
630  *working_vector_petsc = 1.0;
631  Ainv_full->pointwise_divide(*working_vector_petsc, *row_sum);
632  const auto Ainv_full_old = Ainv_full->clone();
633 
634  // Correct HbyA
635  Ainv_full->add(-1.0, Ainv);
636  working_vector_petsc->pointwise_mult(*Ainv_full, *pressure_gradient[system_i]);
637  working_vector_petsc->pointwise_mult(*working_vector_petsc, *_cell_volumes);
638  HbyA.add(-1.0, *working_vector_petsc);
639 
640  // Correct Ainv
641  Ainv = *Ainv_full_old;
642  }
643 
644  Ainv.pointwise_mult(Ainv, *_cell_volumes);
645  }
646 
647  // We fill the 1/A and H/A functors
649 
650  if (verbose)
651  {
652  _console << "************************************" << std::endl;
653  _console << "DONE Computing HbyA " << std::endl;
654  _console << "************************************" << std::endl;
655  }
656 }
657 
658 std::vector<std::unique_ptr<NumericVector<Number>>> &
659 RhieChowMassFlux::selectPressureGradient(const bool updated_pressure)
660 {
661  if (updated_pressure)
662  {
663  _grad_p_current.clear();
664  for (const auto & component : _pressure_system->gradientContainer())
665  _grad_p_current.push_back(component->clone());
666  }
667 
668  return _grad_p_current;
669 }
std::vector< LinearSystem * > _momentum_systems
Pointers to the linear system(s) in moose corresponding to the momentum equation(s) ...
const MooseEnum _pressure_projection_method
Enumerator for the method used for pressure projection.
virtual void initialize() override
void setupMeshInformation()
Compute the cell volumes on the mesh.
User object responsible for determining the face fluxes using the Rhie-Chow interpolation in a segreg...
unsigned int n_threads()
A functor whose evaluation relies on querying a map where the keys are face info ids and the values c...
const std::set< BoundaryID > & boundaryIDs() const
virtual numeric_index_type local_m() const final
T & getMesh(MooseMesh &mesh)
function to cast mesh
Definition: SCM.h:35
void checkBlocks(const VarType &var) const
Check the block consistency between the passed in var and us.
const std::vector< const ElemInfo *> & elemInfoVector() const
std::vector< libMesh::LinearImplicitSystem * > _momentum_implicit_systems
Pointers to the momentum equation implicit system(s) from libmesh.
unsigned int number() const
const std::vector< std::unique_ptr< NumericVector< Number > > > & gradientContainer() const
Real getMassFlux(const FaceInfo &fi) const
Get the face velocity times density (used in advection terms)
void vector_mult(NumericVector< T > &dest, const NumericVector< T > &arg) const
virtual Real computeElemMatrixContribution() override
static const std::string component
Definition: NS.h:153
const Elem & elem() const
virtual std::unique_ptr< NumericVector< T > > zero_clone() const=0
const ElemInfo * neighborInfo() const
const ExecFlagType EXEC_NONE
const Elem * elem() const
registerMooseObject("NavierStokesApp", RhieChowMassFlux)
void computeFaceMassFlux()
Update the values of the face velocities in the containers.
void addFunctor(const std::string &name, const Moose::FunctorBase< T > &functor, const THREAD_ID tid)
MeshBase & mesh
static const std::string density
Definition: NS.h:33
auto raw_value(const Eigen::Map< T > &in)
virtual void setupFaceData(const FaceInfo *face_info)
NumericVector< Number > * rhs
void addAvailableFlags(const ExecFlagType &flag, Args... flags)
const ElemInfo * elemInfo() const
std::vector< std::unique_ptr< NumericVector< Number > > > _Ainv_raw
We hold on to the cell-based 1/A vectors so that we can easily reconstruct the cell velocities as wel...
const LinearSystem * _pressure_system
Pointer to the pressure system.
LinearFVBoundaryCondition * getBoundaryCondition(const BoundaryID bd_id) const
FaceCenteredMapFunctor< RealVectorValue, std::unordered_map< dof_id_type, RealVectorValue > > _Ainv
A map functor from faces to $(1/A)_f$.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
std::vector< unsigned int > _momentum_system_numbers
Numbers of the momentum system(s)
virtual void initialSetup() override
virtual void meshChanged() override
virtual void pointwise_divide(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
virtual const std::string & name() const
virtual Real computeBoundaryMatrixContribution(const LinearFVBoundaryCondition &bc) override
MooseApp & getMooseApp() const
std::vector< unsigned int > _global_momentum_system_numbers
Global numbers of the momentum system(s)
static InputParameters validParams()
SubProblem & _subproblem
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
bool isInternalFace(const FaceInfo &) const
void computeHbyA(const bool with_updated_pressure, const bool verbose)
Computes the inverse of the diagonal (1/A) of the system matrix plus the H/A components for the press...
const std::vector< const FaceInfo *> & faceInfo() const
ValueType evaluate(const FaceInfo *const fi) const
Evaluate the face functor using a FaceInfo argument.
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
FaceCenteredMapFunctor< Real, std::unordered_map< dof_id_type, Real > > _HbyA_flux
A map functor from faces to $HbyA_{ij} = (A_{offdiag}*{(predicted~velocity)} - {Source})_{ij}/A_{ij}$...
const Elem * neighborPtr() const
virtual const NumericVector< Number > *const & currentSolution() const override final
TheWarehouse & theWarehouse() const
std::vector< const MooseLinearVariableFVReal * > _vel
The thread 0 copy of the x-velocity variable.
const MooseLinearVariableFVReal *const _p
The thread 0 copy of the pressure variable.
unsigned int _global_pressure_system_number
Global number of the pressure system.
RhieChowMassFlux(const InputParameters &params)
void linkMomentumPressureSystems(const std::vector< LinearSystem *> &momentum_systems, const LinearSystem &pressure_system, const std::vector< unsigned int > &momentum_system_numbers)
Update the momentum system-related information.
const Moose::Functor< Real > & _rho
Functor describing the density of the fluid.
static InputParameters validParams()
virtual Real computeElemRightHandSideContribution() override
std::unique_ptr< NumericVector< Number > > solution
void computeCellVelocity()
Update the cell values of the velocity variables.
virtual Real computeNeighborMatrixContribution() override
Real getVolumetricFaceFlux(const FaceInfo &fi) const
Get the volumetric face flux (used in advection terms)
void setCurrentFaceArea(const Real area)
virtual void print(std::ostream &os=libMesh::out) const
virtual Real computeBoundaryRHSContribution(const LinearFVBoundaryCondition &bc) override
const Point & normal() const
void paramError(const std::string &param, Args... args) const
unsigned int number() const
virtual const MooseVariableFieldBase & getVariable(const THREAD_ID tid, const std::string &var_name, Moose::VarKindType expected_var_type=Moose::VarKindType::VAR_ANY, Moose::VarFieldType expected_var_field_type=Moose::VarFieldType::VAR_FIELD_ANY) const=0
virtual void get_row(numeric_index_type i, std::vector< numeric_index_type > &indices, std::vector< T > &values) const override
std::vector< std::unique_ptr< NumericVector< Number > > > & selectPressureGradient(const bool updated_pressure)
Select the right pressure gradient field and return a reference to the container. ...
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2)=0
const Elem * elemPtr() const
const std::vector< std::vector< dof_id_type > > & dofIndices() const
virtual void get_diagonal(NumericVector< T > &dest) const override
FaceCenteredMapFunctor< Real, std::unordered_map< dof_id_type, Real > > & _face_mass_flux
A map functor from faces to mass fluxes which are used in the advection terms.
std::vector< std::unique_ptr< NumericVector< Number > > > _grad_p_current
for a PISO iteration we need to hold on to the original pressure gradient field.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
subdomain_id_type subdomain_id() const
FEProblemBase & _fe_problem
SparseMatrix< Number > * matrix
Query query()
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
static const std::string pressure
Definition: NS.h:56
const THREAD_ID _tid
IntRange< T > make_range(T beg, T end)
virtual MooseMesh & mesh() override
dof_id_type id() const
void mooseError(Args &&... args) const
virtual numeric_index_type row_start() const override
std::unique_ptr< NumericVector< Number > > current_local_solution
static InputParameters validParams()
const ConsoleStream _console
std::vector< const FaceInfo * > _flow_face_info
The subset of the FaceInfo objects that actually cover the subdomains which the flow field is defined...
bool hasBlocks(const SubdomainName &name) const
void populateCouplingFunctors(const std::vector< std::unique_ptr< NumericVector< Number >>> &raw_hbya, const std::vector< std::unique_ptr< NumericVector< Number >>> &raw_Ainv)
Populate the face values of the H/A and 1/A fields.
virtual void add(const numeric_index_type i, const T value)=0
void interpolate(InterpMethod m, T &result, const T2 &value1, const T3 &value2, const FaceInfo &fi, const bool one_is_elem)
StateArg currentState()
LinearFVAnisotropicDiffusion * _p_diffusion_kernel
Pointer to the pressure diffusion term in the pressure Poisson equation.
auto index_range(const T &sizable)
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
std::vector< std::unique_ptr< NumericVector< Number > > > _HbyA_raw
We hold on to the cell-based HbyA vectors so that we can easily reconstruct the cell velocities as we...
virtual void pointwise_divide(const NumericVector< T > &vec1, const NumericVector< T > &vec2)=0
void initCouplingField()
Initialize the coupling fields (HbyA and Ainv)
unsigned int state
virtual System & system() override
unsigned int THREAD_ID
uint8_t dof_id_type
void initFaceMassFlux()
Initialize the container for face velocities.
VarFaceNeighbors faceType(const std::pair< unsigned int, unsigned int > &var_sys) const
std::unique_ptr< NumericVector< Number > > _cell_volumes
We will hold a vector of cell volumes to make sure we can do volume corrections rapidly.
const unsigned int _dim
The dimension of the mesh, e.g. 3 for hexes and tets, 2 for quads and tris.