https://mooseframework.inl.gov
SubChannel1PhaseProblem.h
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 #pragma once
11 #include "ExternalProblem.h"
12 #include "PostprocessorInterface.h"
13 #include "SubChannelApp.h"
14 #include "QuadSubChannelMesh.h"
15 #include "SolutionHandle.h"
17 #include <petscdm.h>
18 #include <petscdmda.h>
19 #include <petscksp.h>
20 #include <petscsys.h>
21 #include <petscvec.h>
22 #include <petscsnes.h>
23 
28 {
29 public:
31  virtual ~SubChannel1PhaseProblem();
32 
33  virtual void externalSolve() override;
34  virtual void syncSolutions(Direction direction) override;
35  virtual bool solverSystemConverged(const unsigned int) override;
36  virtual void initialSetup() override;
37 
38 protected:
41  {
42  Mat A;
43  Vec x;
44  };
45 
47  {
48  int i_ch;
51 
53  virtual Real computeFrictionFactor(FrictionStruct friction_args) = 0;
55  void computeWijFromSolve(int iblock);
57  void computeSumWij(int iblock);
59  void computeMdot(int iblock);
61  void computeWijPrime(int iblock);
63  virtual Real computeBeta(unsigned int i_gap, unsigned int iz) = 0;
65  void computeDP(int iblock);
67  void computeP(int iblock);
69  virtual void computeh(int iblock) = 0;
71  void computeT(int iblock);
73  void computeRho(int iblock);
75  void computeMu(int iblock);
77  void computeWijResidual(int iblock);
79  virtual Real computeAddedHeatPin(unsigned int i_ch, unsigned int iz) = 0;
81  Real computeAddedHeatDuct(unsigned int i_ch, unsigned int iz);
85  PetscErrorCode petscSnesSolver(int iblock,
86  const libMesh::DenseVector<Real> & solution,
89  friend PetscErrorCode formFunction(SNES snes, Vec x, Vec f, void * ctx);
90 
92  PetscErrorCode implicitPetscSolve(int iblock);
93 
95  virtual void initializeSolution() = 0;
96 
98  PetscScalar computeInterpolationCoefficients(PetscScalar Peclet = 0.0);
99  PetscScalar
100  computeInterpolatedValue(PetscScalar topValue, PetscScalar botValue, PetscScalar Peclet = 0.0);
101 
102  PetscErrorCode cleanUp();
105  unsigned int _n_blocks;
110  const Real _g_grav;
111  const Real & _kij;
112  unsigned int _n_cells;
113  unsigned int _n_gaps;
114  unsigned int _n_pins;
115  unsigned int _n_channels;
116  unsigned int _block_size;
121  std::vector<Real> _z_grid;
126  const bool _compute_density;
128  const bool _compute_viscosity;
130  const bool _compute_power;
132  const bool _pin_mesh_exist;
134  const bool _duct_mesh_exist;
142  const Real & _CT;
144  const Real & _P_tol;
146  const Real & _T_tol;
148  const int & _T_maxit;
150  const PetscReal & _rtol;
152  const PetscReal & _atol;
154  const PetscReal & _dtol;
156  const PetscInt & _maxit;
160  const bool _implicit_bool;
164  const bool _segregated_bool;
170  const bool _deformation;
171 
174  std::unique_ptr<SolutionHandle> _mdot_soln;
175  std::unique_ptr<SolutionHandle> _SumWij_soln;
176  std::unique_ptr<SolutionHandle> _P_soln;
177  std::unique_ptr<SolutionHandle> _DP_soln;
178  std::unique_ptr<SolutionHandle> _h_soln;
179  std::unique_ptr<SolutionHandle> _T_soln;
180  std::unique_ptr<SolutionHandle> _Tpin_soln;
181  std::unique_ptr<SolutionHandle> _Dpin_soln;
182  std::unique_ptr<SolutionHandle> _rho_soln;
183  std::unique_ptr<SolutionHandle> _mu_soln;
184  std::unique_ptr<SolutionHandle> _S_flow_soln;
185  std::unique_ptr<SolutionHandle> _w_perim_soln;
186  std::unique_ptr<SolutionHandle> _q_prime_soln;
187  std::unique_ptr<SolutionHandle> _q_prime_duct_soln; // Only used for ducted assemblies
188  std::unique_ptr<SolutionHandle> _Tduct_soln; // Only used for ducted assemblies
189  std::unique_ptr<SolutionHandle> _displacement_soln;
190 
192  inline PetscErrorCode createPetscVector(Vec & v, PetscInt n)
193  {
195  LibmeshPetscCall(VecCreate(PETSC_COMM_WORLD, &v));
196  LibmeshPetscCall(PetscObjectSetName((PetscObject)v, "Solution"));
197  LibmeshPetscCall(VecSetSizes(v, PETSC_DECIDE, n));
198  LibmeshPetscCall(VecSetFromOptions(v));
199  LibmeshPetscCall(VecZeroEntries(v));
200  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
201  }
202 
203  inline PetscErrorCode createPetscMatrix(Mat & M, PetscInt n, PetscInt m)
204  {
206  LibmeshPetscCall(MatCreate(PETSC_COMM_WORLD, &M));
207  LibmeshPetscCall(MatSetSizes(M, PETSC_DECIDE, PETSC_DECIDE, n, m));
208  LibmeshPetscCall(MatSetFromOptions(M));
209  LibmeshPetscCall(MatSetUp(M));
210  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
211  }
212 
213  template <class T>
214  PetscErrorCode populateVectorFromDense(Vec & x,
215  const T & solution,
216  const unsigned int first_axial_level,
217  const unsigned int last_axial_level,
218  const unsigned int cross_dimension);
219  template <class T>
220  PetscErrorCode populateDenseFromVector(const Vec & x,
221  T & solution,
222  const unsigned int first_axial_level,
223  const unsigned int last_axial_level,
224  const unsigned int cross_dimension);
225  template <class T>
226  PetscErrorCode populateVectorFromHandle(Vec & x,
227  const T & solution,
228  const unsigned int first_axial_level,
229  const unsigned int last_axial_level,
230  const unsigned int cross_dimension);
231  template <class T>
232  PetscErrorCode populateSolutionChan(const Vec & x,
233  T & solution,
234  const unsigned int first_axial_level,
235  const unsigned int last_axial_level,
236  const unsigned int cross_dimension);
237 
238  template <class T>
239  PetscErrorCode populateSolutionGap(const Vec & x,
240  T & solution,
241  const unsigned int first_axial_level,
242  const unsigned int last_axial_level,
243  const unsigned int cross_dimension);
244 
249  Vec _Wij_vec;
250  Vec _prod;
251  Vec _prodp;
257 
283 
301 
318  PetscInt _global_counter = 0;
319 
321  PetscScalar _added_K = 0.0;
322  PetscScalar _added_K_old = 1000.0;
323  PetscScalar _max_sumWij;
324  PetscScalar _max_sumWij_new;
325  PetscScalar _correction_factor = 1.0;
326 
327 public:
328  static InputParameters validParams();
329 };
330 
331 template <class T>
332 PetscErrorCode
334  T & loc_solution,
335  const unsigned int first_axial_level,
336  const unsigned int last_axial_level,
337  const unsigned int cross_dimension)
338 {
339  PetscScalar * xx;
340 
342  LibmeshPetscCall(VecGetArray(x, &xx));
343  for (unsigned int iz = first_axial_level; iz < last_axial_level + 1; iz++)
344  {
345  unsigned int iz_ind = iz - first_axial_level;
346  for (unsigned int i_l = 0; i_l < cross_dimension; i_l++)
347  {
348  loc_solution(i_l, iz) = xx[iz_ind * cross_dimension + i_l];
349  }
350  }
351  LibmeshPetscCall(VecRestoreArray(x, &xx));
352 
353  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
354 }
355 
356 template <class T>
357 PetscErrorCode
359  const T & loc_solution,
360  const unsigned int first_axial_level,
361  const unsigned int last_axial_level,
362  const unsigned int cross_dimension)
363 {
364  PetscScalar * xx;
365 
367  LibmeshPetscCall(VecGetArray(x, &xx));
368  for (unsigned int iz = first_axial_level; iz < last_axial_level + 1; iz++)
369  {
370  unsigned int iz_ind = iz - first_axial_level;
371  for (unsigned int i_l = 0; i_l < cross_dimension; i_l++)
372  {
373  auto * loc_node = _subchannel_mesh.getChannelNode(i_l, iz);
374  xx[iz_ind * cross_dimension + i_l] = loc_solution(loc_node);
375  }
376  }
377  LibmeshPetscCall(VecRestoreArray(x, &xx));
378 
379  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
380 }
381 
382 template <class T>
383 PetscErrorCode
385  const T & loc_solution,
386  const unsigned int first_axial_level,
387  const unsigned int last_axial_level,
388  const unsigned int cross_dimension)
389 {
390  PetscScalar * xx;
392  LibmeshPetscCall(VecGetArray(x, &xx));
393  for (unsigned int iz = first_axial_level; iz < last_axial_level; iz++)
394  {
395  unsigned int iz_ind = iz - first_axial_level;
396  for (unsigned int i_l = 0; i_l < cross_dimension; i_l++)
397  {
398  xx[iz_ind * cross_dimension + i_l] = loc_solution(i_l, iz);
399  }
400  }
401  LibmeshPetscCall(VecRestoreArray(x, &xx));
402  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
403 }
404 
405 template <class T>
406 PetscErrorCode
408  T & loc_solution,
409  const unsigned int first_axial_level,
410  const unsigned int last_axial_level,
411  const unsigned int cross_dimension)
412 {
413  PetscScalar * xx;
415  LibmeshPetscCall(VecGetArray(x, &xx));
416  Node * loc_node;
417  for (unsigned int iz = first_axial_level; iz < last_axial_level + 1; iz++)
418  {
419  unsigned int iz_ind = iz - first_axial_level;
420  for (unsigned int i_l = 0; i_l < cross_dimension; i_l++)
421  {
422  loc_node = _subchannel_mesh.getChannelNode(i_l, iz);
423  loc_solution.set(loc_node, xx[iz_ind * cross_dimension + i_l]);
424  }
425  }
426  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
427 }
428 
429 template <class T>
430 PetscErrorCode
432  T & loc_solution,
433  const unsigned int first_axial_level,
434  const unsigned int last_axial_level,
435  const unsigned int cross_dimension)
436 {
437  PetscScalar * xx;
439  LibmeshPetscCall(VecGetArray(x, &xx));
440  for (unsigned int iz = first_axial_level; iz < last_axial_level + 1; iz++)
441  {
442  unsigned int iz_ind = iz - first_axial_level;
443  for (unsigned int i_l = 0; i_l < cross_dimension; i_l++)
444  {
445  loc_solution(iz * cross_dimension + i_l) = xx[iz_ind * cross_dimension + i_l];
446  }
447  }
448  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
449 }
Mat _amc_sys_mdot_mat
Axial momentum system matrix.
const bool _pin_mesh_exist
Flag that informs if there is a pin mesh or not.
void computeRho(int iblock)
Computes Density per channel for block iblock.
Vec _amc_gravity_rhs
Axial momentum conservation - buoyancy force No implicit matrix.
std::unique_ptr< SolutionHandle > _Tduct_soln
void computeSumWij(int iblock)
Computes net diversion crossflow per channel for block iblock.
unsigned int _n_blocks
number of axial blocks
const bool _compute_power
Flag that informs if we need to solve the Enthalpy/Temperature equations or not.
std::unique_ptr< SolutionHandle > _T_soln
std::unique_ptr< SolutionHandle > _h_soln
const PetscReal & _dtol
The divergence tolerance for the ksp linear solver.
virtual Real computeBeta(unsigned int i_gap, unsigned int iz)=0
Computes turbulent mixing coefficient.
const bool _monolithic_thermal_bool
Thermal monolithic bool.
std::unique_ptr< SolutionHandle > _q_prime_duct_soln
Mat _cmc_friction_force_mat
Cross momentum conservation - friction force.
void computeP(int iblock)
Computes Pressure per channel for block iblock.
Real _TR
Flag that activates or deactivates the transient parts of the equations we solve by multiplication...
libMesh::DenseMatrix< Real > _Wij_old
const PostprocessorValue & _P_out
Outlet Pressure.
void computeT(int iblock)
Computes Temperature per channel for block iblock.
void computeMdot(int iblock)
Computes mass flow per channel for block iblock.
std::unique_ptr< SolutionHandle > _S_flow_soln
const bool _compute_density
Flag that activates or deactivates the calculation of density.
PetscErrorCode populateSolutionChan(const Vec &x, T &solution, const unsigned int first_axial_level, const unsigned int last_axial_level, const unsigned int cross_dimension)
Mat _amc_friction_force_mat
Axial momentum conservation - friction force.
PetscScalar _added_K
Added resistances for monolithic convergence.
Standard return structure for reusing in implicit/explicit formulations.
PetscErrorCode createPetscVector(Vec &v, PetscInt n)
Petsc Functions.
PetscFunctionBegin
const bool _duct_mesh_exist
Flag that informs if there is a pin mesh or not.
static InputParameters validParams()
virtual Real computeAddedHeatPin(unsigned int i_ch, unsigned int iz)=0
Computes added heat for channel i_ch and cell iz.
const Real & _T_tol
Convergence tolerance for the temperature loop in external solve.
void computeMu(int iblock)
Computes Viscosity per channel for block iblock.
PetscScalar computeInterpolatedValue(PetscScalar topValue, PetscScalar botValue, PetscScalar Peclet=0.0)
const bool _segregated_bool
Segregated solve.
Mat _cmc_sys_Wij_mat
Lateral momentum system matrix.
bool _converged
Variable that informs whether we exited external solve with a converged solution or not...
const bool _staggered_pressure_bool
Flag to define the usage of staggered or collocated pressure.
virtual void initializeSolution()=0
Function to initialize the solution & geometry fields.
Real computeAddedHeatDuct(unsigned int i_ch, unsigned int iz)
Function that computes the heat flux added by the duct.
Mat _mc_sumWij_mat
Matrices and vectors to be used in implicit assembly Mass conservation Mass conservation - sum of cro...
PetscErrorCode populateVectorFromDense(Vec &x, const T &solution, const unsigned int first_axial_level, const unsigned int last_axial_level, const unsigned int cross_dimension)
const Real & _CT
Turbulent modeling parameter used in axial momentum equation.
virtual Node * getChannelNode(unsigned int i_chan, unsigned iz) const =0
Get the subchannel mesh node for a given channel index and elevation index.
virtual void computeh(int iblock)=0
Computes Enthalpy per channel for block iblock.
std::unique_ptr< SolutionHandle > _rho_soln
Mat _amc_turbulent_cross_flows_mat
Mass conservation - density time derivative No implicit matrix.
PetscErrorCode populateVectorFromHandle(Vec &x, const T &solution, const unsigned int first_axial_level, const unsigned int last_axial_level, const unsigned int cross_dimension)
Mat _amc_advective_derivative_mat
Axial momentum conservation - advective (Eulerian) derivative.
Vec _hc_added_heat_rhs
Enthalpy conservation - source and sink.
std::vector< Real > _z_grid
axial location of nodes
const int & _T_maxit
Maximum iterations for the inner temperature loop.
Mat _amc_pressure_force_mat
Axial momentum conservation - pressure force.
void computeWijResidual(int iblock)
Computes Residual Matrix based on the lateral momentum conservation equation for block iblock...
const std::vector< double > x
Real f(Real x)
Test function for Brents method.
std::unique_ptr< SolutionHandle > _q_prime_soln
virtual void syncSolutions(Direction direction) override
PetscErrorCode populateSolutionGap(const Vec &x, T &solution, const unsigned int first_axial_level, const unsigned int last_axial_level, const unsigned int cross_dimension)
const PetscReal & _atol
The absolute convergence tolerance for the ksp linear solver.
PetscScalar computeInterpolationCoefficients(PetscScalar Peclet=0.0)
Functions that computes the interpolation scheme given the Peclet number.
Real PostprocessorValue
Real root(std::function< Real(Real)> const &f, Real x1, Real x2, Real tol=1.0e-12)
Finds the root of a function using Brent&#39;s method.
Definition: BrentsMethod.C:66
Mat _hc_cross_derivative_mat
Enthalpy conservation - cross flux derivative.
Mat _mc_axial_convection_mat
Mass conservation - axial convection.
Common class for single phase fluid properties.
auto Peclet(const T1 &volume_fraction, const T2 &cp, const T3 &rho, const T4 &vel, const T5 &D_h, const T6 &k)
Compute Peclet number.
Definition: Numerics.h:153
virtual Real computeFrictionFactor(FrictionStruct friction_args)=0
Returns friction factor.
const PetscInt & _maxit
The maximum number of iterations to use for the ksp linear solver.
PetscErrorCode createPetscMatrix(Mat &M, PetscInt n, PetscInt m)
friend PetscErrorCode formFunction(SNES snes, Vec x, Vec f, void *ctx)
This is the residual Vector function in a form compatible with the SNES PETC solvers.
const MooseEnum _interpolation_scheme
The interpolation method used in constructing the systems.
std::unique_ptr< SolutionHandle > _Dpin_soln
Base class for the 1-phase steady-state/transient subchannel solver.
libMesh::DenseMatrix< Real > _Wij_residual_matrix
PetscInt _global_counter
No implicit matrix.
std::unique_ptr< SolutionHandle > _mdot_soln
PetscErrorCode populateDenseFromVector(const Vec &x, T &solution, const unsigned int first_axial_level, const unsigned int last_axial_level, const unsigned int cross_dimension)
Real _dpz_error
average relative error in pressure drop of channels
const PetscReal & _rtol
The relative convergence tolerance, (relative decrease) for the ksp linear solver.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::unique_ptr< SolutionHandle > _displacement_soln
static const std::string v
Definition: NS.h:84
std::unique_ptr< SolutionHandle > _DP_soln
const bool _compute_viscosity
Flag that activates or deactivates the calculation of viscosity.
libMesh::DenseMatrix< Real > & _Wij
virtual void externalSolve() override
void computeWijPrime(int iblock)
Computes turbulent crossflow per gap for block iblock.
std::unique_ptr< SolutionHandle > _SumWij_soln
PetscErrorCode petscSnesSolver(int iblock, const libMesh::DenseVector< Real > &solution, libMesh::DenseVector< Real > &root)
Computes solution of nonlinear equation using snes and provided a residual in a formFunction.
libMesh::DenseVector< Real > residualFunction(int iblock, libMesh::DenseVector< Real > solution)
Computes Residual Vector based on the lateral momentum conservation equation for block iblock & updat...
const bool _verbose_subchannel
Boolean to printout information related to subchannel solve.
Mat _amc_cross_derivative_mat
Axial momentum conservation - cross flux derivative.
const Real & _P_tol
Convergence tolerance for the pressure loop in external solve.
virtual bool solverSystemConverged(const unsigned int) override
struct SubChannel1PhaseProblem::FrictionStruct _friction_args
Mat _hc_advective_derivative_mat
Enthalpy conservation - advective (Eulerian) derivative;.
Base class for subchannel meshes.
void computeWijFromSolve(int iblock)
Computes diversion crossflow per gap for block iblock.
const bool _implicit_bool
Flag to define the usage of a implicit or explicit solution.
Mat _amc_time_derivative_mat
Axial momentum conservation - time derivative.
Mat _hc_time_derivative_mat
Enthalpy Enthalpy conservation - time derivative.
PetscErrorCode implicitPetscSolve(int iblock)
Computes implicit solve using PetSc.
const bool _deformation
Flag that activates the effect of deformation (pin/duct) based on the auxvalues for displacement...
void computeDP(int iblock)
Computes Pressure Drop per channel for block iblock.
SubChannel1PhaseProblem(const InputParameters &params)
Mat _cmc_time_derivative_mat
Cross momentum Cross momentum conservation - time derivative.
libMesh::DenseMatrix< Real > _WijPrime
virtual void initialSetup() override
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
Mat _cmc_pressure_force_mat
Cross momentum conservation - pressure force.
std::unique_ptr< SolutionHandle > _P_soln
std::unique_ptr< SolutionHandle > _w_perim_soln
Mat _cmc_advective_derivative_mat
Cross momentum conservation - advective (Eulerian) derivative.
const SinglePhaseFluidProperties * _fp
Solutions handles and link to TH tables properties.
Mat _hc_sys_h_mat
System matrices.
std::unique_ptr< SolutionHandle > _mu_soln
std::unique_ptr< SolutionHandle > _Tpin_soln