https://mooseframework.inl.gov
SubChannel1PhaseProblem.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 "SystemBase.h"
12 #include "libmesh/petsc_vector.h"
13 #include "libmesh/dense_matrix.h"
14 #include "libmesh/dense_vector.h"
15 #include <iostream>
16 #include <cmath>
17 #include "AuxiliarySystem.h"
18 #include "SCM.h"
19 
20 struct Ctx
21 {
22  int iblock;
24 };
25 
26 PetscErrorCode
27 formFunction(SNES, Vec x, Vec f, void * ctx)
28 {
29  const PetscScalar * xx;
30  PetscScalar * ff;
31  PetscInt size;
32 
34  Ctx * cc = static_cast<Ctx *>(ctx);
35  LibmeshPetscCallQ(VecGetSize(x, &size));
36 
37  libMesh::DenseVector<Real> solution_seed(size, 0.0);
38  LibmeshPetscCallQ(VecGetArrayRead(x, &xx));
39  for (PetscInt i = 0; i < size; i++)
40  solution_seed(i) = xx[i];
41 
42  LibmeshPetscCallQ(VecRestoreArrayRead(x, &xx));
43 
44  libMesh::DenseVector<Real> Wij_residual_vector =
45  cc->schp->residualFunction(cc->iblock, solution_seed);
46 
47  LibmeshPetscCallQ(VecGetArray(f, &ff));
48  for (int i = 0; i < size; i++)
49  ff[i] = Wij_residual_vector(i);
50 
51  LibmeshPetscCallQ(VecRestoreArray(f, &ff));
52  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
53 }
54 
57 {
58  MooseEnum schemes("upwind downwind central_difference exponential", "central_difference");
59  MooseEnum gravity_direction("counter_flow co_flow none", "counter_flow");
62  params.addClassDescription("Base class of the subchannel solvers");
63  params.addRequiredParam<unsigned int>("n_blocks", "The number of blocks in the axial direction");
64  params.addRequiredParam<Real>("CT", "Turbulent modeling parameter");
65  params.addParam<Real>("P_tol", 1e-6, "Pressure tolerance");
66  params.addParam<Real>("T_tol", 1e-6, "Temperature tolerance");
67  params.addParam<int>("T_maxit", 100, "Maximum number of iterations for inner temperature loop");
68  params.addParam<PetscReal>("rtol", 1e-6, "Relative tolerance for ksp solver");
69  params.addParam<PetscReal>("atol", 1e-6, "Absolute tolerance for ksp solver");
70  params.addParam<PetscReal>("dtol", 1e5, "Divergence tolerance or ksp solver");
71  params.addParam<PetscInt>("maxit", 1e4, "Maximum number of iterations for ksp solver");
72  params.addParam<MooseEnum>("interpolation_scheme",
73  schemes,
74  "Interpolation scheme used for the method. Default is exponential");
75  params.addParam<MooseEnum>(
76  "gravity", gravity_direction, "Direction of gravity. Default is counter_flow");
77  params.addParam<bool>(
78  "implicit", false, "Boolean to define the use of explicit or implicit solution.");
79  params.addParam<bool>(
80  "staggered_pressure", false, "Boolean to define the use of explicit or implicit solution.");
81  params.addParam<bool>(
82  "segregated", true, "Boolean to define whether to use a segregated solution.");
83  params.addParam<bool>(
84  "monolithic_thermal", false, "Boolean to define whether to use thermal monolithic solve.");
85  params.addParam<bool>(
86  "verbose_subchannel", false, "Boolean to print out information related to subchannel solve.");
87  params.addParam<bool>(
88  "deformation",
89  false,
90  "Boolean that activates the deformation effect based on values for: displacement, Dpin");
91  params.addRequiredParam<bool>("compute_density", "Flag that enables the calculation of density");
92  params.addRequiredParam<bool>("compute_viscosity",
93  "Flag that enables the calculation of viscosity");
94  params.addRequiredParam<bool>(
95  "compute_power",
96  "Flag that informs whether we solve the Enthalpy/Temperature equations or not");
97  params.addRequiredParam<PostprocessorName>(
98  "P_out", "The postprocessor (or scalar) that provides the value of outlet pressure [Pa]");
99  params.addRequiredParam<UserObjectName>("fp", "Fluid properties user object name");
100  return params;
101 }
102 
104  : ExternalProblem(params),
106  _subchannel_mesh(SCM::getMesh<SubChannelMesh>(_mesh)),
107  _n_blocks(getParam<unsigned int>("n_blocks")),
108  _Wij(declareRestartableData<libMesh::DenseMatrix<Real>>("Wij")),
109  _g_grav(9.81),
110  _kij(_subchannel_mesh.getKij()),
111  _one(1.0),
112  _compute_density(getParam<bool>("compute_density")),
113  _compute_viscosity(getParam<bool>("compute_viscosity")),
114  _compute_power(getParam<bool>("compute_power")),
115  _pin_mesh_exist(_subchannel_mesh.pinMeshExist()),
116  _duct_mesh_exist(_subchannel_mesh.ductMeshExist()),
117  _P_out(getPostprocessorValue("P_out")),
118  _CT(getParam<Real>("CT")),
119  _P_tol(getParam<Real>("P_tol")),
120  _T_tol(getParam<Real>("T_tol")),
121  _T_maxit(getParam<int>("T_maxit")),
122  _rtol(getParam<PetscReal>("rtol")),
123  _atol(getParam<PetscReal>("atol")),
124  _dtol(getParam<PetscReal>("dtol")),
125  _maxit(getParam<PetscInt>("maxit")),
126  _interpolation_scheme(getParam<MooseEnum>("interpolation_scheme")),
127  _gravity_direction(getParam<MooseEnum>("gravity")),
128  _dir_grav(computeGravityDir(_gravity_direction)),
129  _implicit_bool(getParam<bool>("implicit")),
130  _staggered_pressure_bool(getParam<bool>("staggered_pressure")),
131  _segregated_bool(getParam<bool>("segregated")),
132  _monolithic_thermal_bool(getParam<bool>("monolithic_thermal")),
133  _verbose_subchannel(getParam<bool>("verbose_subchannel")),
134  _deformation(getParam<bool>("deformation")),
135  _fp(nullptr),
136  _Tpin_soln(nullptr),
137  _duct_heat_flux_soln(nullptr),
138  _Tduct_soln(nullptr)
139 {
140  // Require only one process
141  if (n_processors() > 1)
142  mooseError("Cannot use more than one MPI process.");
149  // Turbulent crossflow (stuff that live on the gaps)
150  if (!_app.isRestarting() && !_app.isRecovering())
151  {
152  _Wij.resize(_n_gaps, _n_cells + 1);
153  _Wij.zero();
154  }
156  _Wij_old.zero();
158  _WijPrime.zero();
161  _converged = true;
162 
163  // Mass conservation components
164  LibmeshPetscCall(
166  LibmeshPetscCall(createPetscVector(_Wij_vec, _block_size * _n_gaps));
167  LibmeshPetscCall(createPetscVector(_prod, _block_size * _n_channels));
168  LibmeshPetscCall(createPetscVector(_prodp, _block_size * _n_channels));
169  LibmeshPetscCall(createPetscMatrix(
172 
173  // Axial momentum conservation components
174  LibmeshPetscCall(createPetscMatrix(
177  LibmeshPetscCall(createPetscMatrix(
180  LibmeshPetscCall(createPetscMatrix(
183  LibmeshPetscCall(createPetscMatrix(
186  LibmeshPetscCall(createPetscMatrix(
190  LibmeshPetscCall(createPetscMatrix(
193  LibmeshPetscCall(
196 
197  // Lateral momentum conservation components
198  LibmeshPetscCall(
201  LibmeshPetscCall(createPetscMatrix(
204  LibmeshPetscCall(
207  LibmeshPetscCall(
210  LibmeshPetscCall(
213 
214  // Energy conservation components
215  LibmeshPetscCall(createPetscMatrix(
218  LibmeshPetscCall(createPetscMatrix(
221  LibmeshPetscCall(createPetscMatrix(
225  LibmeshPetscCall(
228 
229  if ((_n_blocks == _n_cells) && _implicit_bool)
230  {
231  mooseError(name(),
232  ": When implicit number of blocks can't be equal to number of cells. This will "
233  "cause problems with the subchannel interpolation scheme.");
234  }
235 }
236 
237 void
239 {
241  _fp = &getUserObject<SinglePhaseFluidProperties>(getParam<UserObjectName>("fp"));
242  _mdot_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::MASS_FLOW_RATE));
243  _SumWij_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::SUM_CROSSFLOW));
244  _P_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PRESSURE));
245  _DP_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PRESSURE_DROP));
246  _h_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::ENTHALPY));
247  _T_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::TEMPERATURE));
248  if (_pin_mesh_exist)
249  {
250  _Tpin_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PIN_TEMPERATURE));
251  _Dpin_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PIN_DIAMETER));
252  }
253  _rho_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DENSITY));
254  _mu_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::VISCOSITY));
255  _S_flow_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::SURFACE_AREA));
256  _w_perim_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::WETTED_PERIMETER));
257  _q_prime_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::LINEAR_HEAT_RATE));
259  std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DISPLACEMENT));
260  if (_duct_mesh_exist)
261  {
263  std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DUCT_HEAT_FLUX));
264  _Tduct_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DUCT_TEMPERATURE));
265  }
266 }
267 
269 {
270  PetscErrorCode ierr = cleanUp();
271  if (ierr)
272  mooseError(name(), ": Error in memory cleanup");
273 }
274 
275 PetscErrorCode
277 {
279  // We need to clean up the petsc matrices/vectors
280  // Mass conservation components
281  LibmeshPetscCall(MatDestroy(&_mc_sumWij_mat));
282  LibmeshPetscCall(VecDestroy(&_Wij_vec));
283  LibmeshPetscCall(VecDestroy(&_prod));
284  LibmeshPetscCall(VecDestroy(&_prodp));
285  LibmeshPetscCall(MatDestroy(&_mc_axial_convection_mat));
286  LibmeshPetscCall(VecDestroy(&_mc_axial_convection_rhs));
287 
288  // Axial momentum conservation components
289  LibmeshPetscCall(MatDestroy(&_amc_turbulent_cross_flows_mat));
290  LibmeshPetscCall(VecDestroy(&_amc_turbulent_cross_flows_rhs));
291  LibmeshPetscCall(MatDestroy(&_amc_time_derivative_mat));
292  LibmeshPetscCall(VecDestroy(&_amc_time_derivative_rhs));
293  LibmeshPetscCall(MatDestroy(&_amc_advective_derivative_mat));
294  LibmeshPetscCall(VecDestroy(&_amc_advective_derivative_rhs));
295  LibmeshPetscCall(MatDestroy(&_amc_cross_derivative_mat));
296  LibmeshPetscCall(VecDestroy(&_amc_cross_derivative_rhs));
297  LibmeshPetscCall(MatDestroy(&_amc_friction_force_mat));
298  LibmeshPetscCall(VecDestroy(&_amc_friction_force_rhs));
299  LibmeshPetscCall(VecDestroy(&_amc_gravity_rhs));
300  LibmeshPetscCall(MatDestroy(&_amc_pressure_force_mat));
301  LibmeshPetscCall(VecDestroy(&_amc_pressure_force_rhs));
302  LibmeshPetscCall(MatDestroy(&_amc_sys_mdot_mat));
303  LibmeshPetscCall(VecDestroy(&_amc_sys_mdot_rhs));
304 
305  // Lateral momentum conservation components
306  LibmeshPetscCall(MatDestroy(&_cmc_time_derivative_mat));
307  LibmeshPetscCall(VecDestroy(&_cmc_time_derivative_rhs));
308  LibmeshPetscCall(MatDestroy(&_cmc_advective_derivative_mat));
309  LibmeshPetscCall(VecDestroy(&_cmc_advective_derivative_rhs));
310  LibmeshPetscCall(MatDestroy(&_cmc_friction_force_mat));
311  LibmeshPetscCall(VecDestroy(&_cmc_friction_force_rhs));
312  LibmeshPetscCall(MatDestroy(&_cmc_pressure_force_mat));
313  LibmeshPetscCall(VecDestroy(&_cmc_pressure_force_rhs));
314  LibmeshPetscCall(MatDestroy(&_cmc_sys_Wij_mat));
315  LibmeshPetscCall(VecDestroy(&_cmc_sys_Wij_rhs));
316 
317  // Energy conservation components
318  LibmeshPetscCall(MatDestroy(&_hc_time_derivative_mat));
319  LibmeshPetscCall(VecDestroy(&_hc_time_derivative_rhs));
320  LibmeshPetscCall(MatDestroy(&_hc_advective_derivative_mat));
321  LibmeshPetscCall(VecDestroy(&_hc_advective_derivative_rhs));
322  LibmeshPetscCall(MatDestroy(&_hc_cross_derivative_mat));
323  LibmeshPetscCall(VecDestroy(&_hc_cross_derivative_rhs));
324  LibmeshPetscCall(VecDestroy(&_hc_added_heat_rhs));
325  LibmeshPetscCall(MatDestroy(&_hc_sys_h_mat));
326  LibmeshPetscCall(VecDestroy(&_hc_sys_h_rhs));
327 
328  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
329 }
330 
331 bool
333 {
334  return _converged;
335 }
336 
337 PetscScalar
339 {
340  switch (_interpolation_scheme)
341  {
342  case 0: // upwind interpolation
343  return 1.0;
344  case 1: // downwind interpolation
345  return 0.0;
346  case 2: // central_difference interpolation
347  return 0.5;
348  case 3: // exponential interpolation (Peclet limited)
349  return ((Peclet - 1.0) * std::exp(Peclet) + 1) / (Peclet * (std::exp(Peclet) - 1.) + 1e-10);
350  default:
351  mooseError(name(),
352  ": Interpolation scheme should be a string: upwind, downwind, central_difference, "
353  "exponential");
354  }
355 }
356 
357 PetscScalar
359  PetscScalar botValue,
360  PetscScalar Peclet)
361 {
363  return alpha * botValue + (1.0 - alpha) * topValue;
364 }
365 
366 void
368 {
369  unsigned int last_node = (iblock + 1) * _block_size;
370  unsigned int first_node = iblock * _block_size + 1;
371  // Initial guess, port crossflow of block (iblock) into a vector that will act as my initial guess
372  libMesh::DenseVector<Real> solution_seed(_n_gaps * _block_size, 0.0);
373  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
374  {
375  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
376  {
377  int i = _n_gaps * (iz - first_node) + i_gap; // column wise transfer
378  solution_seed(i) = _Wij(i_gap, iz);
379  }
380  }
381 
382  // Solving the combined lateral momentum equation for Wij using a PETSc solver and update vector
383  // root
385  LibmeshPetscCall(petscSnesSolver(iblock, solution_seed, root));
386 
387  // Assign the solution to the cross-flow matrix
388  int i = 0;
389  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
390  {
391  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
392  {
393  _Wij(i_gap, iz) = root(i);
394  i++;
395  }
396  }
397 }
398 
399 void
401 {
402  unsigned int last_node = (iblock + 1) * _block_size;
403  unsigned int first_node = iblock * _block_size + 1;
404  // Add to solution vector if explicit
405  if (!_implicit_bool)
406  {
407  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
408  {
409  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
410  {
411  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
412  Real sumWij = 0.0;
413  // Calculate sum of crossflow into channel i from channels j around i
414  unsigned int counter = 0;
415  for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
416  {
417  sumWij += _subchannel_mesh.getCrossflowSign(i_ch, counter) * _Wij(i_gap, iz);
418  counter++;
419  }
420  // The net crossflow coming out of cell i [kg/sec]
421  _SumWij_soln->set(node_out, sumWij);
422  }
423  }
424  }
425  // Add to matrix if explicit
426  else
427  {
428  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
429  {
430  unsigned int iz_ind = iz - first_node;
431  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
432  {
433  // Calculate sum of crossflow into channel i from channels j around i
434  unsigned int counter = 0;
435  for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
436  {
437  PetscInt row = i_ch + _n_channels * iz_ind;
438  PetscInt col = i_gap + _n_gaps * iz_ind;
439  PetscScalar value = _subchannel_mesh.getCrossflowSign(i_ch, counter);
440  LibmeshPetscCall(MatSetValues(_mc_sumWij_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
441  counter++;
442  }
443  }
444  }
445  LibmeshPetscCall(MatAssemblyBegin(_mc_sumWij_mat, MAT_FINAL_ASSEMBLY));
446  LibmeshPetscCall(MatAssemblyEnd(_mc_sumWij_mat, MAT_FINAL_ASSEMBLY));
447  if (_segregated_bool)
448  {
449  Vec loc_prod;
450  Vec loc_Wij;
451  LibmeshPetscCall(VecDuplicate(_amc_sys_mdot_rhs, &loc_prod));
452  LibmeshPetscCall(VecDuplicate(_Wij_vec, &loc_Wij));
454  loc_Wij, _Wij, first_node, last_node, _n_gaps));
455  LibmeshPetscCall(MatMult(_mc_sumWij_mat, loc_Wij, loc_prod));
456  LibmeshPetscCall(populateSolutionChan<SolutionHandle>(
457  loc_prod, *_SumWij_soln, first_node, last_node, _n_channels));
458  LibmeshPetscCall(VecDestroy(&loc_prod));
459  LibmeshPetscCall(VecDestroy(&loc_Wij));
460  }
461  }
462 }
463 
464 void
466 {
467  unsigned int last_node = (iblock + 1) * _block_size;
468  unsigned int first_node = iblock * _block_size + 1;
469  if (!_implicit_bool)
470  {
471  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
472  {
473  auto dz = _z_grid[iz] - _z_grid[iz - 1];
474  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
475  {
476  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
477  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
478  auto volume = dz * (*_S_flow_soln)(node_in);
479  auto time_term = _TR * ((*_rho_soln)(node_out)-_rho_soln->old(node_out)) * volume / _dt;
480  // Wij positive out of i into j;
481  auto mdot_out = (*_mdot_soln)(node_in) - (*_SumWij_soln)(node_out)-time_term;
482  if (mdot_out < 0)
483  {
484  _console << "Wij = : " << _Wij << "\n";
485  mooseError(name(),
486  " : Calculation of negative mass flow mdot_out = : ",
487  mdot_out,
488  " Axial Level= : ",
489  iz,
490  " - Implicit solves are required for recirculating flow.");
491  }
492  _mdot_soln->set(node_out, mdot_out); // kg/sec
493  }
494  }
495  }
496  else
497  {
498  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
499  {
500  auto dz = _z_grid[iz] - _z_grid[iz - 1];
501  auto iz_ind = iz - first_node;
502  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
503  {
504  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
505  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
506  auto volume = dz * (*_S_flow_soln)(node_in);
507 
508  // Adding time derivative to the RHS
509  auto time_term = _TR * ((*_rho_soln)(node_out)-_rho_soln->old(node_out)) * volume / _dt;
510  PetscInt row_vec = i_ch + _n_channels * iz_ind;
511  PetscScalar value_vec = -1.0 * time_term;
512  LibmeshPetscCall(
513  VecSetValues(_mc_axial_convection_rhs, 1, &row_vec, &value_vec, INSERT_VALUES));
514 
515  // Imposing bottom boundary condition or adding of diagonal elements
516  if (iz == first_node)
517  {
518  PetscScalar value_vec = (*_mdot_soln)(node_in);
519  PetscInt row_vec = i_ch + _n_channels * iz_ind;
520  LibmeshPetscCall(
521  VecSetValues(_mc_axial_convection_rhs, 1, &row_vec, &value_vec, ADD_VALUES));
522  }
523  else
524  {
525  PetscInt row = i_ch + _n_channels * iz_ind;
526  PetscInt col = i_ch + _n_channels * (iz_ind - 1);
527  PetscScalar value = -1.0;
528  LibmeshPetscCall(
529  MatSetValues(_mc_axial_convection_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
530  }
531 
532  // Adding diagonal elements
533  PetscInt row = i_ch + _n_channels * iz_ind;
534  PetscInt col = i_ch + _n_channels * iz_ind;
535  PetscScalar value = 1.0;
536  LibmeshPetscCall(
537  MatSetValues(_mc_axial_convection_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
538 
539  // Adding cross flows RHS
540  if (_segregated_bool)
541  {
542  PetscScalar value_vec_2 = -1.0 * (*_SumWij_soln)(node_out);
543  PetscInt row_vec_2 = i_ch + _n_channels * iz_ind;
544  LibmeshPetscCall(
545  VecSetValues(_mc_axial_convection_rhs, 1, &row_vec_2, &value_vec_2, ADD_VALUES));
546  }
547  }
548  }
549  LibmeshPetscCall(MatAssemblyBegin(_mc_axial_convection_mat, MAT_FINAL_ASSEMBLY));
550  LibmeshPetscCall(MatAssemblyEnd(_mc_axial_convection_mat, MAT_FINAL_ASSEMBLY));
551 
552  if (_segregated_bool)
553  {
554  KSP ksploc;
555  PC pc;
556  Vec sol;
557  LibmeshPetscCall(VecDuplicate(_mc_axial_convection_rhs, &sol));
558  LibmeshPetscCall(KSPCreate(PETSC_COMM_WORLD, &ksploc));
559  LibmeshPetscCall(KSPSetOperators(ksploc, _mc_axial_convection_mat, _mc_axial_convection_mat));
560  LibmeshPetscCall(KSPGetPC(ksploc, &pc));
561  LibmeshPetscCall(PCSetType(pc, PCJACOBI));
562  LibmeshPetscCall(KSPSetTolerances(ksploc, _rtol, _atol, _dtol, _maxit));
563  LibmeshPetscCall(KSPSetFromOptions(ksploc));
564  LibmeshPetscCall(KSPSolve(ksploc, _mc_axial_convection_rhs, sol));
565  LibmeshPetscCall(populateSolutionChan<SolutionHandle>(
566  sol, *_mdot_soln, first_node, last_node, _n_channels));
567  LibmeshPetscCall(VecZeroEntries(_mc_axial_convection_rhs));
568  LibmeshPetscCall(KSPDestroy(&ksploc));
569  LibmeshPetscCall(VecDestroy(&sol));
570  }
571  }
572 }
573 
574 void
576 {
577  unsigned int last_node = (iblock + 1) * _block_size;
578  unsigned int first_node = iblock * _block_size + 1;
579  if (!_implicit_bool)
580  {
581  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
582  {
583  auto k_grid = _subchannel_mesh.getKGrid();
584  auto dz = _z_grid[iz] - _z_grid[iz - 1];
585  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
586  {
587  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
588  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
589  auto rho_in = (*_rho_soln)(node_in);
590  auto rho_out = (*_rho_soln)(node_out);
591  auto mu_in = (*_mu_soln)(node_in);
592  auto S = (*_S_flow_soln)(node_in);
593  auto w_perim = (*_w_perim_soln)(node_in);
594  // hydraulic diameter in the i direction
595  auto Dh_i = 4.0 * S / w_perim;
596  auto time_term = _TR * ((*_mdot_soln)(node_out)-_mdot_soln->old(node_out)) * dz / _dt -
597  dz * 2.0 * (*_mdot_soln)(node_out) * (rho_out - _rho_soln->old(node_out)) /
598  rho_in / _dt;
599  auto mass_term1 =
600  std::pow((*_mdot_soln)(node_out), 2.0) * (1.0 / S / rho_out - 1.0 / S / rho_in);
601  auto mass_term2 = -2.0 * (*_mdot_soln)(node_out) * (*_SumWij_soln)(node_out) / S / rho_in;
602  auto crossflow_term = 0.0;
603  auto turbulent_term = 0.0;
604  unsigned int counter = 0;
605  for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
606  {
607  auto chans = _subchannel_mesh.getGapChannels(i_gap);
608  unsigned int ii_ch = chans.first;
609  unsigned int jj_ch = chans.second;
610  auto * node_in_i = _subchannel_mesh.getChannelNode(ii_ch, iz - 1);
611  auto * node_in_j = _subchannel_mesh.getChannelNode(jj_ch, iz - 1);
612  auto * node_out_i = _subchannel_mesh.getChannelNode(ii_ch, iz);
613  auto * node_out_j = _subchannel_mesh.getChannelNode(jj_ch, iz);
614  auto rho_i = (*_rho_soln)(node_in_i);
615  auto rho_j = (*_rho_soln)(node_in_j);
616  auto Si = (*_S_flow_soln)(node_in_i);
617  auto Sj = (*_S_flow_soln)(node_in_j);
618  Real u_star = 0.0;
619  // figure out donor axial velocity
620  if (_Wij(i_gap, iz) > 0.0)
621  u_star = (*_mdot_soln)(node_out_i) / Si / rho_i;
622  else
623  u_star = (*_mdot_soln)(node_out_j) / Sj / rho_j;
624 
625  crossflow_term +=
626  _subchannel_mesh.getCrossflowSign(i_ch, counter) * _Wij(i_gap, iz) * u_star;
627 
628  turbulent_term += _WijPrime(i_gap, iz) * (2 * (*_mdot_soln)(node_out) / rho_in / S -
629  (*_mdot_soln)(node_out_j) / Sj / rho_j -
630  (*_mdot_soln)(node_out_i) / Si / rho_i);
631  counter++;
632  }
633  turbulent_term *= _CT;
634  auto Re = (((*_mdot_soln)(node_in) / S) * Dh_i / mu_in);
635  _friction_args.Re = Re;
636  _friction_args.i_ch = i_ch;
637  _friction_args.S = S;
638  _friction_args.w_perim = w_perim;
641  auto ki = 0.0;
642  if ((*_mdot_soln)(node_out) >= 0)
643  ki = k_grid[i_ch][iz - 1];
644  else
645  ki = k_grid[i_ch][iz];
646  auto friction_term = (fi * dz / Dh_i + ki) * 0.5 *
647  (*_mdot_soln)(node_out)*std::abs((*_mdot_soln)(node_out)) /
648  (S * (*_rho_soln)(node_out));
649  auto gravity_term = _dir_grav * _g_grav * (*_rho_soln)(node_out)*dz * S;
650  auto DP = std::pow(S, -1.0) * (time_term + mass_term1 + mass_term2 + crossflow_term +
651  turbulent_term + friction_term + gravity_term); // Pa
652  _DP_soln->set(node_out, DP);
653  }
654  }
655  }
656  else
657  {
658  LibmeshPetscCall(MatZeroEntries(_amc_time_derivative_mat));
659  LibmeshPetscCall(MatZeroEntries(_amc_advective_derivative_mat));
660  LibmeshPetscCall(MatZeroEntries(_amc_cross_derivative_mat));
661  LibmeshPetscCall(MatZeroEntries(_amc_friction_force_mat));
662  LibmeshPetscCall(VecZeroEntries(_amc_time_derivative_rhs));
663  LibmeshPetscCall(VecZeroEntries(_amc_advective_derivative_rhs));
664  LibmeshPetscCall(VecZeroEntries(_amc_cross_derivative_rhs));
665  LibmeshPetscCall(VecZeroEntries(_amc_friction_force_rhs));
666  LibmeshPetscCall(VecZeroEntries(_amc_gravity_rhs));
667  LibmeshPetscCall(MatZeroEntries(_amc_sys_mdot_mat));
668  LibmeshPetscCall(VecZeroEntries(_amc_sys_mdot_rhs));
669  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
670  {
671  auto k_grid = _subchannel_mesh.getKGrid();
672  auto dz = _z_grid[iz] - _z_grid[iz - 1];
673  auto iz_ind = iz - first_node;
674  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
675  {
676  // inlet and outlet nodes
677  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
678  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
679 
680  // interpolation weight coefficient
681  PetscScalar Pe = 0.5;
682  if (_interpolation_scheme == 3)
683  {
684  // Compute the Peclet number
685  auto S_in = (*_S_flow_soln)(node_in);
686  auto S_out = (*_S_flow_soln)(node_out);
687  auto S_interp = computeInterpolatedValue(S_out, S_in, 0.5);
688  auto w_perim_in = (*_w_perim_soln)(node_in);
689  auto w_perim_out = (*_w_perim_soln)(node_out);
690  auto w_perim_interp = this->computeInterpolatedValue(w_perim_out, w_perim_in, 0.5);
691  auto mdot_loc =
692  this->computeInterpolatedValue((*_mdot_soln)(node_out), (*_mdot_soln)(node_in), 0.5);
693  auto mu_in = (*_mu_soln)(node_in);
694  auto mu_out = (*_mu_soln)(node_out);
695  auto mu_interp = this->computeInterpolatedValue(mu_out, mu_in, 0.5);
696  auto Dh_i = 4.0 * S_interp / w_perim_interp;
697  // Compute friction
698  auto Re = ((mdot_loc / S_interp) * Dh_i / mu_interp);
699  _friction_args.Re = Re;
700  _friction_args.i_ch = i_ch;
701  _friction_args.S = S_interp;
702  _friction_args.w_perim = w_perim_interp;
705  auto ki = 0.0;
706  if ((*_mdot_soln)(node_out) >= 0)
707  ki = k_grid[i_ch][iz - 1];
708  else
709  ki = k_grid[i_ch][iz];
710  Pe = 1.0 / ((fi * dz / Dh_i + ki) * 0.5) * mdot_loc / std::abs(mdot_loc);
711  }
713 
714  // inlet, outlet, and interpolated density
715  auto rho_in = (*_rho_soln)(node_in);
716  auto rho_out = (*_rho_soln)(node_out);
717  auto rho_interp = computeInterpolatedValue(rho_out, rho_in, Pe);
718 
719  // inlet, outlet, and interpolated viscosity
720  auto mu_in = (*_mu_soln)(node_in);
721  auto mu_out = (*_mu_soln)(node_out);
722  auto mu_interp = computeInterpolatedValue(mu_out, mu_in, Pe);
723 
724  // inlet, outlet, and interpolated axial surface area
725  auto S_in = (*_S_flow_soln)(node_in);
726  auto S_out = (*_S_flow_soln)(node_out);
727  auto S_interp = computeInterpolatedValue(S_out, S_in, Pe);
728 
729  // inlet, outlet, and interpolated wetted perimeter
730  auto w_perim_in = (*_w_perim_soln)(node_in);
731  auto w_perim_out = (*_w_perim_soln)(node_out);
732  auto w_perim_interp = computeInterpolatedValue(w_perim_out, w_perim_in, Pe);
733 
734  // hydraulic diameter in the i direction
735  auto Dh_i = 4.0 * S_interp / w_perim_interp;
736 
738  if (iz == first_node)
739  {
740  PetscScalar value_vec_tt = -1.0 * _TR * alpha * (*_mdot_soln)(node_in)*dz / _dt;
741  PetscInt row_vec_tt = i_ch + _n_channels * iz_ind;
742  LibmeshPetscCall(
743  VecSetValues(_amc_time_derivative_rhs, 1, &row_vec_tt, &value_vec_tt, ADD_VALUES));
744  }
745  else
746  {
747  PetscInt row_tt = i_ch + _n_channels * iz_ind;
748  PetscInt col_tt = i_ch + _n_channels * (iz_ind - 1);
749  PetscScalar value_tt = _TR * alpha * dz / _dt;
750  LibmeshPetscCall(MatSetValues(
751  _amc_time_derivative_mat, 1, &row_tt, 1, &col_tt, &value_tt, INSERT_VALUES));
752  }
753 
754  // Adding diagonal elements
755  PetscInt row_tt = i_ch + _n_channels * iz_ind;
756  PetscInt col_tt = i_ch + _n_channels * iz_ind;
757  PetscScalar value_tt = _TR * (1.0 - alpha) * dz / _dt;
758  LibmeshPetscCall(MatSetValues(
759  _amc_time_derivative_mat, 1, &row_tt, 1, &col_tt, &value_tt, INSERT_VALUES));
760 
761  // Adding RHS elements
762  PetscScalar mdot_old_interp =
763  computeInterpolatedValue(_mdot_soln->old(node_out), _mdot_soln->old(node_in), Pe);
764  PetscScalar value_vec_tt = _TR * mdot_old_interp * dz / _dt;
765  PetscInt row_vec_tt = i_ch + _n_channels * iz_ind;
766  LibmeshPetscCall(
767  VecSetValues(_amc_time_derivative_rhs, 1, &row_vec_tt, &value_vec_tt, ADD_VALUES));
768 
770  if (iz == first_node)
771  {
772  PetscScalar value_vec_at = std::pow((*_mdot_soln)(node_in), 2.0) / (S_in * rho_in);
773  PetscInt row_vec_at = i_ch + _n_channels * iz_ind;
774  LibmeshPetscCall(VecSetValues(
775  _amc_advective_derivative_rhs, 1, &row_vec_at, &value_vec_at, ADD_VALUES));
776  }
777  else
778  {
779  PetscInt row_at = i_ch + _n_channels * iz_ind;
780  PetscInt col_at = i_ch + _n_channels * (iz_ind - 1);
781  PetscScalar value_at = -1.0 * std::abs((*_mdot_soln)(node_in)) / (S_in * rho_in);
782  LibmeshPetscCall(MatSetValues(
783  _amc_advective_derivative_mat, 1, &row_at, 1, &col_at, &value_at, INSERT_VALUES));
784  }
785 
786  // Adding diagonal elements
787  PetscInt row_at = i_ch + _n_channels * iz_ind;
788  PetscInt col_at = i_ch + _n_channels * iz_ind;
789  PetscScalar value_at = std::abs((*_mdot_soln)(node_out)) / (S_out * rho_out);
790  LibmeshPetscCall(MatSetValues(
791  _amc_advective_derivative_mat, 1, &row_at, 1, &col_at, &value_at, INSERT_VALUES));
792 
794  unsigned int counter = 0;
795  unsigned int cross_index = iz; // iz-1;
796  for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
797  {
798  auto chans = _subchannel_mesh.getGapChannels(i_gap);
799  unsigned int ii_ch = chans.first;
800  unsigned int jj_ch = chans.second;
801  auto * node_in_i = _subchannel_mesh.getChannelNode(ii_ch, iz - 1);
802  auto * node_in_j = _subchannel_mesh.getChannelNode(jj_ch, iz - 1);
803  auto * node_out_i = _subchannel_mesh.getChannelNode(ii_ch, iz);
804  auto * node_out_j = _subchannel_mesh.getChannelNode(jj_ch, iz);
805  auto rho_i =
806  computeInterpolatedValue((*_rho_soln)(node_out_i), (*_rho_soln)(node_in_i), Pe);
807  auto rho_j =
808  computeInterpolatedValue((*_rho_soln)(node_out_j), (*_rho_soln)(node_in_j), Pe);
809  auto S_i =
810  computeInterpolatedValue((*_S_flow_soln)(node_out_i), (*_S_flow_soln)(node_in_i), Pe);
811  auto S_j =
812  computeInterpolatedValue((*_S_flow_soln)(node_out_j), (*_S_flow_soln)(node_in_j), Pe);
813  auto u_star = 0.0;
814  // figure out donor axial velocity
815  if (_Wij(i_gap, cross_index) > 0.0)
816  {
817  if (iz == first_node)
818  {
819  u_star = (*_mdot_soln)(node_in_i) / S_i / rho_i;
820  PetscScalar value_vec_ct = -1.0 * alpha *
821  _subchannel_mesh.getCrossflowSign(i_ch, counter) *
822  _Wij(i_gap, cross_index) * u_star;
823  PetscInt row_vec_ct = i_ch + _n_channels * iz_ind;
824  LibmeshPetscCall(VecSetValues(
825  _amc_cross_derivative_rhs, 1, &row_vec_ct, &value_vec_ct, ADD_VALUES));
826  }
827  else
828  {
829  PetscScalar value_ct = alpha * _subchannel_mesh.getCrossflowSign(i_ch, counter) *
830  _Wij(i_gap, cross_index) / S_i / rho_i;
831  PetscInt row_ct = i_ch + _n_channels * iz_ind;
832  PetscInt col_ct = ii_ch + _n_channels * (iz_ind - 1);
833  LibmeshPetscCall(MatSetValues(
834  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
835  }
836  PetscScalar value_ct = (1.0 - alpha) *
837  _subchannel_mesh.getCrossflowSign(i_ch, counter) *
838  _Wij(i_gap, cross_index) / S_i / rho_i;
839  PetscInt row_ct = i_ch + _n_channels * iz_ind;
840  PetscInt col_ct = ii_ch + _n_channels * iz_ind;
841  LibmeshPetscCall(MatSetValues(
842  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
843  }
844  else if (_Wij(i_gap, cross_index) < 0.0) // _Wij=0 operations not necessary
845  {
846  if (iz == first_node)
847  {
848  u_star = (*_mdot_soln)(node_in_j) / S_j / rho_j;
849  PetscScalar value_vec_ct = -1.0 * alpha *
850  _subchannel_mesh.getCrossflowSign(i_ch, counter) *
851  _Wij(i_gap, cross_index) * u_star;
852  PetscInt row_vec_ct = i_ch + _n_channels * iz_ind;
853  LibmeshPetscCall(VecSetValues(
854  _amc_cross_derivative_rhs, 1, &row_vec_ct, &value_vec_ct, ADD_VALUES));
855  }
856  else
857  {
858  PetscScalar value_ct = alpha * _subchannel_mesh.getCrossflowSign(i_ch, counter) *
859  _Wij(i_gap, cross_index) / S_j / rho_j;
860  PetscInt row_ct = i_ch + _n_channels * iz_ind;
861  PetscInt col_ct = jj_ch + _n_channels * (iz_ind - 1);
862  LibmeshPetscCall(MatSetValues(
863  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
864  }
865  PetscScalar value_ct = (1.0 - alpha) *
866  _subchannel_mesh.getCrossflowSign(i_ch, counter) *
867  _Wij(i_gap, cross_index) / S_j / rho_j;
868  PetscInt row_ct = i_ch + _n_channels * iz_ind;
869  PetscInt col_ct = jj_ch + _n_channels * iz_ind;
870  LibmeshPetscCall(MatSetValues(
871  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
872  }
873 
874  if (iz == first_node)
875  {
876  PetscScalar value_vec_ct = -2.0 * alpha * (*_mdot_soln)(node_in)*_CT *
877  _WijPrime(i_gap, cross_index) / (rho_interp * S_interp);
878  value_vec_ct += alpha * (*_mdot_soln)(node_in_j)*_CT * _WijPrime(i_gap, cross_index) /
879  (rho_j * S_j);
880  value_vec_ct += alpha * (*_mdot_soln)(node_in_i)*_CT * _WijPrime(i_gap, cross_index) /
881  (rho_i * S_i);
882  PetscInt row_vec_ct = i_ch + _n_channels * iz_ind;
883  LibmeshPetscCall(
884  VecSetValues(_amc_cross_derivative_rhs, 1, &row_vec_ct, &value_vec_ct, ADD_VALUES));
885  }
886  else
887  {
888  PetscScalar value_center_ct =
889  2.0 * alpha * _CT * _WijPrime(i_gap, cross_index) / (rho_interp * S_interp);
890  PetscInt row_ct = i_ch + _n_channels * iz_ind;
891  PetscInt col_ct = i_ch + _n_channels * (iz_ind - 1);
892  LibmeshPetscCall(MatSetValues(
893  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_center_ct, ADD_VALUES));
894 
895  PetscScalar value_left_ct =
896  -1.0 * alpha * _CT * _WijPrime(i_gap, cross_index) / (rho_j * S_j);
897  row_ct = i_ch + _n_channels * iz_ind;
898  col_ct = jj_ch + _n_channels * (iz_ind - 1);
899  LibmeshPetscCall(MatSetValues(
900  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_left_ct, ADD_VALUES));
901 
902  PetscScalar value_right_ct =
903  -1.0 * alpha * _CT * _WijPrime(i_gap, cross_index) / (rho_i * S_i);
904  row_ct = i_ch + _n_channels * iz_ind;
905  col_ct = ii_ch + _n_channels * (iz_ind - 1);
906  LibmeshPetscCall(MatSetValues(
907  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_right_ct, ADD_VALUES));
908  }
909 
910  PetscScalar value_center_ct =
911  2.0 * (1.0 - alpha) * _CT * _WijPrime(i_gap, cross_index) / (rho_interp * S_interp);
912  PetscInt row_ct = i_ch + _n_channels * iz_ind;
913  PetscInt col_ct = i_ch + _n_channels * iz_ind;
914  LibmeshPetscCall(MatSetValues(
915  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_center_ct, ADD_VALUES));
916 
917  PetscScalar value_left_ct =
918  -1.0 * (1.0 - alpha) * _CT * _WijPrime(i_gap, cross_index) / (rho_j * S_j);
919  row_ct = i_ch + _n_channels * iz_ind;
920  col_ct = jj_ch + _n_channels * iz_ind;
921  LibmeshPetscCall(MatSetValues(
922  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_left_ct, ADD_VALUES));
923 
924  PetscScalar value_right_ct =
925  -1.0 * (1.0 - alpha) * _CT * _WijPrime(i_gap, cross_index) / (rho_i * S_i);
926  row_ct = i_ch + _n_channels * iz_ind;
927  col_ct = ii_ch + _n_channels * iz_ind;
928  LibmeshPetscCall(MatSetValues(
929  _amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_right_ct, ADD_VALUES));
930  counter++;
931  }
932 
934  PetscScalar mdot_interp =
935  computeInterpolatedValue((*_mdot_soln)(node_out), (*_mdot_soln)(node_in), Pe);
936  auto Re = ((mdot_interp / S_interp) * Dh_i / mu_interp);
937  _friction_args.Re = Re;
938  _friction_args.i_ch = i_ch;
939  _friction_args.S = S_interp;
940  _friction_args.w_perim = w_perim_interp;
943  auto ki = 0.0;
944  if ((*_mdot_soln)(node_out) >= 0)
945  ki = k_grid[i_ch][iz - 1];
946  else
947  ki = k_grid[i_ch][iz];
948  auto coef = (fi * dz / Dh_i + ki) * 0.5 * std::abs((*_mdot_soln)(node_out)) /
949  (S_interp * rho_interp);
950  if (iz == first_node)
951  {
952  PetscScalar value_vec = -1.0 * alpha * coef * (*_mdot_soln)(node_in);
953  PetscInt row_vec = i_ch + _n_channels * iz_ind;
954  LibmeshPetscCall(
955  VecSetValues(_amc_friction_force_rhs, 1, &row_vec, &value_vec, ADD_VALUES));
956  }
957  else
958  {
959  PetscInt row = i_ch + _n_channels * iz_ind;
960  PetscInt col = i_ch + _n_channels * (iz_ind - 1);
961  PetscScalar value = alpha * coef;
962  LibmeshPetscCall(
963  MatSetValues(_amc_friction_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
964  }
965 
966  // Adding diagonal elements
967  PetscInt row = i_ch + _n_channels * iz_ind;
968  PetscInt col = i_ch + _n_channels * iz_ind;
969  PetscScalar value = (1.0 - alpha) * coef;
970  LibmeshPetscCall(
971  MatSetValues(_amc_friction_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
972 
974  PetscScalar value_vec = _dir_grav * -1.0 * _g_grav * rho_interp * dz * S_interp;
975  PetscInt row_vec = i_ch + _n_channels * iz_ind;
976  LibmeshPetscCall(VecSetValues(_amc_gravity_rhs, 1, &row_vec, &value_vec, ADD_VALUES));
977  }
978  }
980  LibmeshPetscCall(MatZeroEntries(_amc_sys_mdot_mat));
981  LibmeshPetscCall(VecZeroEntries(_amc_sys_mdot_rhs));
982  LibmeshPetscCall(MatAssemblyBegin(_amc_time_derivative_mat, MAT_FINAL_ASSEMBLY));
983  LibmeshPetscCall(MatAssemblyEnd(_amc_time_derivative_mat, MAT_FINAL_ASSEMBLY));
984  LibmeshPetscCall(MatAssemblyBegin(_amc_advective_derivative_mat, MAT_FINAL_ASSEMBLY));
985  LibmeshPetscCall(MatAssemblyEnd(_amc_advective_derivative_mat, MAT_FINAL_ASSEMBLY));
986  LibmeshPetscCall(MatAssemblyBegin(_amc_cross_derivative_mat, MAT_FINAL_ASSEMBLY));
987  LibmeshPetscCall(MatAssemblyEnd(_amc_cross_derivative_mat, MAT_FINAL_ASSEMBLY));
988  LibmeshPetscCall(MatAssemblyBegin(_amc_friction_force_mat, MAT_FINAL_ASSEMBLY));
989  LibmeshPetscCall(MatAssemblyEnd(_amc_friction_force_mat, MAT_FINAL_ASSEMBLY));
990  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
991  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
992  // Matrix
993 #if !PETSC_VERSION_LESS_THAN(3, 15, 0)
994  LibmeshPetscCall(
995  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_time_derivative_mat, UNKNOWN_NONZERO_PATTERN));
996  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
997  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
998  LibmeshPetscCall(
999  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_advective_derivative_mat, UNKNOWN_NONZERO_PATTERN));
1000  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1001  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1002  LibmeshPetscCall(
1003  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_cross_derivative_mat, UNKNOWN_NONZERO_PATTERN));
1004  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1005  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1006  LibmeshPetscCall(
1007  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_friction_force_mat, UNKNOWN_NONZERO_PATTERN));
1008 #else
1009  LibmeshPetscCall(
1010  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_time_derivative_mat, DIFFERENT_NONZERO_PATTERN));
1011  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1012  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1013  LibmeshPetscCall(
1014  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_advective_derivative_mat, DIFFERENT_NONZERO_PATTERN));
1015  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1016  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1017  LibmeshPetscCall(
1018  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_cross_derivative_mat, DIFFERENT_NONZERO_PATTERN));
1019  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1020  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1021  LibmeshPetscCall(
1022  MatAXPY(_amc_sys_mdot_mat, 1.0, _amc_friction_force_mat, DIFFERENT_NONZERO_PATTERN));
1023 #endif
1024  LibmeshPetscCall(MatAssemblyBegin(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1025  LibmeshPetscCall(MatAssemblyEnd(_amc_sys_mdot_mat, MAT_FINAL_ASSEMBLY));
1026  // RHS
1027  LibmeshPetscCall(VecAXPY(_amc_sys_mdot_rhs, 1.0, _amc_time_derivative_rhs));
1028  LibmeshPetscCall(VecAXPY(_amc_sys_mdot_rhs, 1.0, _amc_advective_derivative_rhs));
1029  LibmeshPetscCall(VecAXPY(_amc_sys_mdot_rhs, 1.0, _amc_cross_derivative_rhs));
1030  LibmeshPetscCall(VecAXPY(_amc_sys_mdot_rhs, 1.0, _amc_friction_force_rhs));
1031  LibmeshPetscCall(VecAXPY(_amc_sys_mdot_rhs, 1.0, _amc_gravity_rhs));
1032  if (_segregated_bool)
1033  {
1034  // Assembly the matrix system
1035  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
1036  _prod, *_mdot_soln, first_node, last_node, _n_channels));
1037  Vec ls;
1038  LibmeshPetscCall(VecDuplicate(_amc_sys_mdot_rhs, &ls));
1039  LibmeshPetscCall(MatMult(_amc_sys_mdot_mat, _prod, ls));
1040  LibmeshPetscCall(VecAXPY(ls, -1.0, _amc_sys_mdot_rhs));
1041  PetscScalar * xx;
1042  LibmeshPetscCall(VecGetArray(ls, &xx));
1043  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1044  {
1045  auto iz_ind = iz - first_node;
1046  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1047  {
1048  // Setting nodes
1049  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1050  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
1051 
1052  // inlet, outlet, and interpolated axial surface area
1053  auto S_in = (*_S_flow_soln)(node_in);
1054  auto S_out = (*_S_flow_soln)(node_out);
1055  auto S_interp = computeInterpolatedValue(S_out, S_in, 0.5);
1056 
1057  // Setting solutions
1058  if (S_interp != 0)
1059  {
1060  auto DP = std::pow(S_interp, -1.0) * xx[iz_ind * _n_channels + i_ch];
1061  _DP_soln->set(node_out, DP);
1062  }
1063  else
1064  {
1065  auto DP = 0.0;
1066  _DP_soln->set(node_out, DP);
1067  }
1068  }
1069  }
1070  LibmeshPetscCall(VecDestroy(&ls));
1071  }
1072  }
1073 }
1074 
1075 void
1077 {
1078  unsigned int last_node = (iblock + 1) * _block_size;
1079  unsigned int first_node = iblock * _block_size + 1;
1080  if (!_implicit_bool)
1081  {
1083  {
1084  for (unsigned int iz = last_node; iz > first_node - 1; iz--)
1085  {
1086  // Calculate pressure in the inlet of the cell assuming known outlet
1087  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1088  {
1089  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
1090  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1091  // update Pressure solution
1092  _P_soln->set(node_in, (*_P_soln)(node_out) + (*_DP_soln)(node_out));
1093  }
1094  }
1095  }
1096  else
1097  {
1098  for (unsigned int iz = last_node; iz > first_node - 1; iz--)
1099  {
1100  // Calculate pressure in the inlet of the cell assuming known outlet
1101  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1102  {
1103  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
1104  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1105  // update Pressure solution
1106  // Note: assuming uniform axial discretization in the curren code
1107  // We will need to update this later if we allow non-uniform refinements in the axial
1108  // direction
1109  PetscScalar Pe = 0.5;
1111  if (iz == last_node)
1112  {
1113  _P_soln->set(node_in, (*_P_soln)(node_out) + (*_DP_soln)(node_out) / 2.0);
1114  }
1115  else
1116  {
1117  _P_soln->set(node_in,
1118  (*_P_soln)(node_out) + (1.0 - alpha) * (*_DP_soln)(node_out) +
1119  alpha * (*_DP_soln)(node_in));
1120  }
1121  }
1122  }
1123  }
1124  }
1125  else
1126  {
1128  {
1129  LibmeshPetscCall(VecZeroEntries(_amc_pressure_force_rhs));
1130  for (unsigned int iz = last_node; iz > first_node - 1; iz--)
1131  {
1132  auto iz_ind = iz - first_node;
1133  // Calculate pressure in the inlet of the cell assuming known outlet
1134  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1135  {
1136  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
1137  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1138 
1139  // inlet, outlet, and interpolated axial surface area
1140  auto S_in = (*_S_flow_soln)(node_in);
1141  auto S_out = (*_S_flow_soln)(node_out);
1142  auto S_interp = computeInterpolatedValue(S_out, S_in, 0.5);
1143 
1144  // Creating matrix of coefficients
1145  PetscInt row = i_ch + _n_channels * iz_ind;
1146  PetscInt col = i_ch + _n_channels * iz_ind;
1147  PetscScalar value = -1.0 * S_interp;
1148  LibmeshPetscCall(
1149  MatSetValues(_amc_pressure_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
1150 
1151  if (iz == last_node)
1152  {
1153  PetscScalar value = -1.0 * (*_P_soln)(node_out)*S_interp;
1154  PetscInt row = i_ch + _n_channels * iz_ind;
1155  LibmeshPetscCall(VecSetValues(_amc_pressure_force_rhs, 1, &row, &value, ADD_VALUES));
1156  }
1157  else
1158  {
1159  PetscInt row = i_ch + _n_channels * iz_ind;
1160  PetscInt col = i_ch + _n_channels * (iz_ind + 1);
1161  PetscScalar value = 1.0 * S_interp;
1162  LibmeshPetscCall(
1163  MatSetValues(_amc_pressure_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
1164  }
1165 
1166  if (_segregated_bool)
1167  {
1168  auto dp_out = (*_DP_soln)(node_out);
1169  PetscScalar value_v = -1.0 * dp_out * S_interp;
1170  PetscInt row_v = i_ch + _n_channels * iz_ind;
1171  LibmeshPetscCall(
1172  VecSetValues(_amc_pressure_force_rhs, 1, &row_v, &value_v, ADD_VALUES));
1173  }
1174  }
1175  }
1176  // Solving pressure problem
1177  LibmeshPetscCall(MatAssemblyBegin(_amc_pressure_force_mat, MAT_FINAL_ASSEMBLY));
1178  LibmeshPetscCall(MatAssemblyEnd(_amc_pressure_force_mat, MAT_FINAL_ASSEMBLY));
1179  if (_segregated_bool)
1180  {
1181  KSP ksploc;
1182  PC pc;
1183  Vec sol;
1184  LibmeshPetscCall(VecDuplicate(_amc_pressure_force_rhs, &sol));
1185  LibmeshPetscCall(KSPCreate(PETSC_COMM_WORLD, &ksploc));
1186  LibmeshPetscCall(KSPSetOperators(ksploc, _amc_pressure_force_mat, _amc_pressure_force_mat));
1187  LibmeshPetscCall(KSPGetPC(ksploc, &pc));
1188  LibmeshPetscCall(PCSetType(pc, PCJACOBI));
1189  LibmeshPetscCall(KSPSetTolerances(ksploc, _rtol, _atol, _dtol, _maxit));
1190  LibmeshPetscCall(KSPSetFromOptions(ksploc));
1191  LibmeshPetscCall(KSPSolve(ksploc, _amc_pressure_force_rhs, sol));
1192  PetscScalar * xx;
1193  LibmeshPetscCall(VecGetArray(sol, &xx));
1194  // update Pressure solution
1195  for (unsigned int iz = last_node; iz > first_node - 1; iz--)
1196  {
1197  auto iz_ind = iz - first_node;
1198  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1199  {
1200  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1201  PetscScalar value = xx[iz_ind * _n_channels + i_ch];
1202  _P_soln->set(node_in, value);
1203  }
1204  }
1205  LibmeshPetscCall(VecZeroEntries(_amc_pressure_force_rhs));
1206  LibmeshPetscCall(KSPDestroy(&ksploc));
1207  LibmeshPetscCall(VecDestroy(&sol));
1208  }
1209  }
1210  else
1211  {
1212  LibmeshPetscCall(VecZeroEntries(_amc_pressure_force_rhs));
1213  for (unsigned int iz = last_node; iz > first_node - 1; iz--)
1214  {
1215  auto iz_ind = iz - first_node;
1216  // Calculate pressure in the inlet of the cell assuming known outlet
1217  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1218  {
1219  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
1220  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1221 
1222  // inlet, outlet, and interpolated axial surface area
1223  auto S_in = (*_S_flow_soln)(node_in);
1224  auto S_out = (*_S_flow_soln)(node_out);
1225  auto S_interp = computeInterpolatedValue(S_out, S_in, 0.5);
1226 
1227  // Creating matrix of coefficients
1228  PetscInt row = i_ch + _n_channels * iz_ind;
1229  PetscInt col = i_ch + _n_channels * iz_ind;
1230  PetscScalar value = -1.0 * S_interp;
1231  LibmeshPetscCall(
1232  MatSetValues(_amc_pressure_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
1233 
1234  if (iz == last_node)
1235  {
1236  PetscScalar value = -1.0 * (*_P_soln)(node_out)*S_interp;
1237  PetscInt row = i_ch + _n_channels * iz_ind;
1238  LibmeshPetscCall(VecSetValues(_amc_pressure_force_rhs, 1, &row, &value, ADD_VALUES));
1239 
1240  auto dp_out = (*_DP_soln)(node_out);
1241  PetscScalar value_v = -1.0 * dp_out / 2.0 * S_interp;
1242  PetscInt row_v = i_ch + _n_channels * iz_ind;
1243  LibmeshPetscCall(
1244  VecSetValues(_amc_pressure_force_rhs, 1, &row_v, &value_v, ADD_VALUES));
1245  }
1246  else
1247  {
1248  PetscInt row = i_ch + _n_channels * iz_ind;
1249  PetscInt col = i_ch + _n_channels * (iz_ind + 1);
1250  PetscScalar value = 1.0 * S_interp;
1251  LibmeshPetscCall(
1252  MatSetValues(_amc_pressure_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
1253 
1254  if (_segregated_bool)
1255  {
1256  auto dp_in = (*_DP_soln)(node_in);
1257  auto dp_out = (*_DP_soln)(node_out);
1258  auto dp_interp = computeInterpolatedValue(dp_out, dp_in, 0.5);
1259  PetscScalar value_v = -1.0 * dp_interp * S_interp;
1260  PetscInt row_v = i_ch + _n_channels * iz_ind;
1261  LibmeshPetscCall(
1262  VecSetValues(_amc_pressure_force_rhs, 1, &row_v, &value_v, ADD_VALUES));
1263  }
1264  }
1265  }
1266  }
1267  // Solving pressure problem
1268  LibmeshPetscCall(MatAssemblyBegin(_amc_pressure_force_mat, MAT_FINAL_ASSEMBLY));
1269  LibmeshPetscCall(MatAssemblyEnd(_amc_pressure_force_mat, MAT_FINAL_ASSEMBLY));
1270  if (_verbose_subchannel)
1271  _console << "Block: " << iblock << " - Axial momentum pressure force matrix assembled"
1272  << std::endl;
1273 
1274  if (_segregated_bool)
1275  {
1276  KSP ksploc;
1277  PC pc;
1278  Vec sol;
1279  LibmeshPetscCall(VecDuplicate(_amc_pressure_force_rhs, &sol));
1280  LibmeshPetscCall(KSPCreate(PETSC_COMM_WORLD, &ksploc));
1281  LibmeshPetscCall(KSPSetOperators(ksploc, _amc_pressure_force_mat, _amc_pressure_force_mat));
1282  LibmeshPetscCall(KSPGetPC(ksploc, &pc));
1283  LibmeshPetscCall(PCSetType(pc, PCJACOBI));
1284  LibmeshPetscCall(KSPSetTolerances(ksploc, _rtol, _atol, _dtol, _maxit));
1285  LibmeshPetscCall(KSPSetFromOptions(ksploc));
1286  LibmeshPetscCall(KSPSolve(ksploc, _amc_pressure_force_rhs, sol));
1287  PetscScalar * xx;
1288  LibmeshPetscCall(VecGetArray(sol, &xx));
1289  // update Pressure solution
1290  for (unsigned int iz = last_node; iz > first_node - 1; iz--)
1291  {
1292  auto iz_ind = iz - first_node;
1293  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1294  {
1295  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1296  PetscScalar value = xx[iz_ind * _n_channels + i_ch];
1297  _P_soln->set(node_in, value);
1298  }
1299  }
1300  LibmeshPetscCall(VecZeroEntries(_amc_pressure_force_rhs));
1301  LibmeshPetscCall(KSPDestroy(&ksploc));
1302  LibmeshPetscCall(VecDestroy(&sol));
1303  }
1304  }
1305  }
1306 }
1307 
1308 void
1310 {
1311  unsigned int last_node = (iblock + 1) * _block_size;
1312  unsigned int first_node = iblock * _block_size + 1;
1313  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1314  {
1315  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1316  {
1317  auto * node = _subchannel_mesh.getChannelNode(i_ch, iz);
1318  _T_soln->set(node, _fp->T_from_p_h((*_P_soln)(node) + _P_out, (*_h_soln)(node)));
1319  }
1320  }
1321 }
1322 
1323 void
1325 {
1326  unsigned int last_node = (iblock + 1) * _block_size;
1327  unsigned int first_node = iblock * _block_size + 1;
1328  if (iblock == 0)
1329  {
1330  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1331  {
1332  auto * node = _subchannel_mesh.getChannelNode(i_ch, 0);
1333  _rho_soln->set(node, _fp->rho_from_p_T((*_P_soln)(node) + _P_out, (*_T_soln)(node)));
1334  }
1335  }
1336  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1337  {
1338  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1339  {
1340  auto * node = _subchannel_mesh.getChannelNode(i_ch, iz);
1341  _rho_soln->set(node, _fp->rho_from_p_T((*_P_soln)(node) + _P_out, (*_T_soln)(node)));
1342  }
1343  }
1344 }
1345 
1346 void
1348 {
1349  unsigned int last_node = (iblock + 1) * _block_size;
1350  unsigned int first_node = iblock * _block_size + 1;
1351  if (iblock == 0)
1352  {
1353  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1354  {
1355  auto * node = _subchannel_mesh.getChannelNode(i_ch, 0);
1356  _mu_soln->set(node, _fp->mu_from_p_T((*_P_soln)(node) + _P_out, (*_T_soln)(node)));
1357  }
1358  }
1359  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1360  {
1361  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1362  {
1363  auto * node = _subchannel_mesh.getChannelNode(i_ch, iz);
1364  _mu_soln->set(node, _fp->mu_from_p_T((*_P_soln)(node) + _P_out, (*_T_soln)(node)));
1365  }
1366  }
1367 }
1368 
1369 void
1371 {
1372  // Cross flow residual
1373  if (!_implicit_bool)
1374  {
1375  unsigned int last_node = (iblock + 1) * _block_size;
1376  unsigned int first_node = iblock * _block_size + 1;
1377  const Real & pitch = _subchannel_mesh.getPitch();
1378  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1379  {
1380  auto dz = _z_grid[iz] - _z_grid[iz - 1];
1381  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
1382  {
1383  auto chans = _subchannel_mesh.getGapChannels(i_gap);
1384  unsigned int i_ch = chans.first;
1385  unsigned int j_ch = chans.second;
1386  auto * node_in_i = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1387  auto * node_out_i = _subchannel_mesh.getChannelNode(i_ch, iz);
1388  auto * node_in_j = _subchannel_mesh.getChannelNode(j_ch, iz - 1);
1389  auto * node_out_j = _subchannel_mesh.getChannelNode(j_ch, iz);
1390  auto rho_i = (*_rho_soln)(node_in_i);
1391  auto rho_j = (*_rho_soln)(node_in_j);
1392  auto Si = (*_S_flow_soln)(node_in_i);
1393  auto Sj = (*_S_flow_soln)(node_in_j);
1394  auto Sij = dz * _subchannel_mesh.getGapWidth(iz, i_gap);
1395  auto Lij = pitch;
1396  // total local form loss in the ij direction
1397  auto friction_term = _kij * _Wij(i_gap, iz) * std::abs(_Wij(i_gap, iz));
1398  auto DPij = (*_P_soln)(node_in_i) - (*_P_soln)(node_in_j);
1399  // Figure out donor cell density
1400  auto rho_star = 0.0;
1401  if (_Wij(i_gap, iz) > 0.0)
1402  rho_star = rho_i;
1403  else if (_Wij(i_gap, iz) < 0.0)
1404  rho_star = rho_j;
1405  else
1406  rho_star = (rho_i + rho_j) / 2.0;
1407  auto mass_term_out =
1408  (*_mdot_soln)(node_out_i) / Si / rho_i + (*_mdot_soln)(node_out_j) / Sj / rho_j;
1409  auto mass_term_in =
1410  (*_mdot_soln)(node_in_i) / Si / rho_i + (*_mdot_soln)(node_in_j) / Sj / rho_j;
1411  auto term_out = Sij * rho_star * (Lij / dz) * mass_term_out * _Wij(i_gap, iz);
1412  auto term_in = Sij * rho_star * (Lij / dz) * mass_term_in * _Wij(i_gap, iz - 1);
1413  auto inertia_term = term_out - term_in;
1414  auto pressure_term = 2 * std::pow(Sij, 2.0) * DPij * rho_star;
1415  auto time_term =
1416  _TR * 2.0 * (_Wij(i_gap, iz) - _Wij_old(i_gap, iz)) * Lij * Sij * rho_star / _dt;
1417 
1418  _Wij_residual_matrix(i_gap, iz - 1 - iblock * _block_size) =
1419  time_term + friction_term + inertia_term - pressure_term;
1420  }
1421  }
1422  }
1423  else
1424  {
1425  // Initializing to zero the elements of the lateral momentum assembly
1426  LibmeshPetscCall(MatZeroEntries(_cmc_time_derivative_mat));
1427  LibmeshPetscCall(MatZeroEntries(_cmc_advective_derivative_mat));
1428  LibmeshPetscCall(MatZeroEntries(_cmc_friction_force_mat));
1429  LibmeshPetscCall(MatZeroEntries(_cmc_pressure_force_mat));
1430  LibmeshPetscCall(VecZeroEntries(_cmc_time_derivative_rhs));
1431  LibmeshPetscCall(VecZeroEntries(_cmc_advective_derivative_rhs));
1432  LibmeshPetscCall(VecZeroEntries(_cmc_friction_force_rhs));
1433  LibmeshPetscCall(VecZeroEntries(_cmc_pressure_force_rhs));
1434  LibmeshPetscCall(MatZeroEntries(_cmc_sys_Wij_mat));
1435  LibmeshPetscCall(VecZeroEntries(_cmc_sys_Wij_rhs));
1436  unsigned int last_node = (iblock + 1) * _block_size;
1437  unsigned int first_node = iblock * _block_size + 1;
1438  const Real & pitch = _subchannel_mesh.getPitch();
1439  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1440  {
1441  auto dz = _z_grid[iz] - _z_grid[iz - 1];
1442  auto iz_ind = iz - first_node;
1443  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
1444  {
1445  auto chans = _subchannel_mesh.getGapChannels(i_gap);
1446  unsigned int i_ch = chans.first;
1447  unsigned int j_ch = chans.second;
1448  auto * node_in_i = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1449  auto * node_out_i = _subchannel_mesh.getChannelNode(i_ch, iz);
1450  auto * node_in_j = _subchannel_mesh.getChannelNode(j_ch, iz - 1);
1451  auto * node_out_j = _subchannel_mesh.getChannelNode(j_ch, iz);
1452 
1453  // inlet, outlet, and interpolated densities
1454  auto rho_i_in = (*_rho_soln)(node_in_i);
1455  auto rho_i_out = (*_rho_soln)(node_out_i);
1456  auto rho_i_interp = computeInterpolatedValue(rho_i_out, rho_i_in, 0.5);
1457  auto rho_j_in = (*_rho_soln)(node_in_j);
1458  auto rho_j_out = (*_rho_soln)(node_out_j);
1459  auto rho_j_interp = computeInterpolatedValue(rho_j_out, rho_j_in, 0.5);
1460 
1461  // inlet, outlet, and interpolated areas
1462  auto S_i_in = (*_S_flow_soln)(node_in_i);
1463  auto S_i_out = (*_S_flow_soln)(node_out_i);
1464  auto S_j_in = (*_S_flow_soln)(node_in_j);
1465  auto S_j_out = (*_S_flow_soln)(node_out_j);
1466 
1467  // Cross-sectional gap area
1468  auto Sij = dz * _subchannel_mesh.getGapWidth(iz, i_gap);
1469  auto Lij = pitch;
1470 
1471  // Figure out donor cell density
1472  auto rho_star = 0.0;
1473  if (_Wij(i_gap, iz) > 0.0)
1474  rho_star = rho_i_interp;
1475  else if (_Wij(i_gap, iz) < 0.0)
1476  rho_star = rho_j_interp;
1477  else
1478  rho_star = (rho_i_interp + rho_j_interp) / 2.0;
1479 
1480  // Assembling time derivative
1481  PetscScalar time_factor = _TR * Lij * Sij * rho_star / _dt;
1482  PetscInt row_td = i_gap + _n_gaps * iz_ind;
1483  PetscInt col_td = i_gap + _n_gaps * iz_ind;
1484  PetscScalar value_td = time_factor;
1485  LibmeshPetscCall(MatSetValues(
1486  _cmc_time_derivative_mat, 1, &row_td, 1, &col_td, &value_td, INSERT_VALUES));
1487  PetscScalar value_td_rhs = time_factor * _Wij_old(i_gap, iz);
1488  LibmeshPetscCall(
1489  VecSetValues(_cmc_time_derivative_rhs, 1, &row_td, &value_td_rhs, INSERT_VALUES));
1490 
1491  // Assembling inertial term
1492  PetscScalar Pe = 0.5;
1494  auto mass_term_out = (*_mdot_soln)(node_out_i) / S_i_out / rho_i_out +
1495  (*_mdot_soln)(node_out_j) / S_j_out / rho_j_out;
1496  auto mass_term_in = (*_mdot_soln)(node_in_i) / S_i_in / rho_i_in +
1497  (*_mdot_soln)(node_in_j) / S_j_in / rho_j_in;
1498  auto term_out = Sij * rho_star * (Lij / dz) * mass_term_out / 2.0;
1499  auto term_in = Sij * rho_star * (Lij / dz) * mass_term_in / 2.0;
1500  if (iz == first_node)
1501  {
1502  PetscInt row_ad = i_gap + _n_gaps * iz_ind;
1503  PetscScalar value_ad = term_in * alpha * _Wij(i_gap, iz - 1);
1504  LibmeshPetscCall(
1505  VecSetValues(_cmc_advective_derivative_rhs, 1, &row_ad, &value_ad, ADD_VALUES));
1506 
1507  PetscInt col_ad = i_gap + _n_gaps * iz_ind;
1508  value_ad = -1.0 * term_in * (1.0 - alpha) + term_out * alpha;
1509  LibmeshPetscCall(MatSetValues(
1510  _cmc_advective_derivative_mat, 1, &row_ad, 1, &col_ad, &value_ad, INSERT_VALUES));
1511 
1512  col_ad = i_gap + _n_gaps * (iz_ind + 1);
1513  value_ad = term_out * (1.0 - alpha);
1514  LibmeshPetscCall(MatSetValues(
1515  _cmc_advective_derivative_mat, 1, &row_ad, 1, &col_ad, &value_ad, INSERT_VALUES));
1516  }
1517  else if (iz == last_node)
1518  {
1519  PetscInt row_ad = i_gap + _n_gaps * iz_ind;
1520  PetscInt col_ad = i_gap + _n_gaps * (iz_ind - 1);
1521  PetscScalar value_ad = -1.0 * term_in * alpha;
1522  LibmeshPetscCall(MatSetValues(
1523  _cmc_advective_derivative_mat, 1, &row_ad, 1, &col_ad, &value_ad, INSERT_VALUES));
1524 
1525  col_ad = i_gap + _n_gaps * iz_ind;
1526  value_ad = -1.0 * term_in * (1.0 - alpha) + term_out * alpha;
1527  LibmeshPetscCall(MatSetValues(
1528  _cmc_advective_derivative_mat, 1, &row_ad, 1, &col_ad, &value_ad, INSERT_VALUES));
1529 
1530  value_ad = -1.0 * term_out * (1.0 - alpha) * _Wij(i_gap, iz);
1531  LibmeshPetscCall(
1532  VecSetValues(_cmc_advective_derivative_rhs, 1, &row_ad, &value_ad, ADD_VALUES));
1533  }
1534  else
1535  {
1536  PetscInt row_ad = i_gap + _n_gaps * iz_ind;
1537  PetscInt col_ad = i_gap + _n_gaps * (iz_ind - 1);
1538  PetscScalar value_ad = -1.0 * term_in * alpha;
1539  LibmeshPetscCall(MatSetValues(
1540  _cmc_advective_derivative_mat, 1, &row_ad, 1, &col_ad, &value_ad, INSERT_VALUES));
1541 
1542  col_ad = i_gap + _n_gaps * iz_ind;
1543  value_ad = -1.0 * term_in * (1.0 - alpha) + term_out * alpha;
1544  LibmeshPetscCall(MatSetValues(
1545  _cmc_advective_derivative_mat, 1, &row_ad, 1, &col_ad, &value_ad, INSERT_VALUES));
1546 
1547  col_ad = i_gap + _n_gaps * (iz_ind + 1);
1548  value_ad = term_out * (1.0 - alpha);
1549  LibmeshPetscCall(MatSetValues(
1550  _cmc_advective_derivative_mat, 1, &row_ad, 1, &col_ad, &value_ad, INSERT_VALUES));
1551  }
1552  // Assembling friction force
1553  PetscInt row_ff = i_gap + _n_gaps * iz_ind;
1554  PetscInt col_ff = i_gap + _n_gaps * iz_ind;
1555  PetscScalar value_ff = _kij * std::abs(_Wij(i_gap, iz)) / 2.0;
1556  LibmeshPetscCall(MatSetValues(
1557  _cmc_friction_force_mat, 1, &row_ff, 1, &col_ff, &value_ff, INSERT_VALUES));
1558 
1559  // Assembling pressure force
1561 
1563  {
1564  PetscScalar pressure_factor = std::pow(Sij, 2.0) * rho_star;
1565  PetscInt row_pf = i_gap + _n_gaps * iz_ind;
1566  PetscInt col_pf = i_ch + _n_channels * iz_ind;
1567  PetscScalar value_pf = -1.0 * alpha * pressure_factor;
1568  LibmeshPetscCall(
1569  MatSetValues(_cmc_pressure_force_mat, 1, &row_pf, 1, &col_pf, &value_pf, ADD_VALUES));
1570  col_pf = j_ch + _n_channels * iz_ind;
1571  value_pf = alpha * pressure_factor;
1572  LibmeshPetscCall(
1573  MatSetValues(_cmc_pressure_force_mat, 1, &row_pf, 1, &col_pf, &value_pf, ADD_VALUES));
1574 
1575  if (iz == last_node)
1576  {
1577  PetscInt row_pf = i_gap + _n_gaps * iz_ind;
1578  PetscScalar value_pf = (1.0 - alpha) * pressure_factor * (*_P_soln)(node_out_i);
1579  LibmeshPetscCall(
1580  VecSetValues(_cmc_pressure_force_rhs, 1, &row_pf, &value_pf, ADD_VALUES));
1581  value_pf = -1.0 * (1.0 - alpha) * pressure_factor * (*_P_soln)(node_out_j);
1582  LibmeshPetscCall(
1583  VecSetValues(_cmc_pressure_force_rhs, 1, &row_pf, &value_pf, ADD_VALUES));
1584  }
1585  else
1586  {
1587  row_pf = i_gap + _n_gaps * iz_ind;
1588  col_pf = i_ch + _n_channels * (iz_ind + 1);
1589  value_pf = -1.0 * (1.0 - alpha) * pressure_factor;
1590  LibmeshPetscCall(MatSetValues(
1591  _cmc_pressure_force_mat, 1, &row_pf, 1, &col_pf, &value_pf, ADD_VALUES));
1592  col_pf = j_ch + _n_channels * (iz_ind + 1);
1593  value_pf = (1.0 - alpha) * pressure_factor;
1594  LibmeshPetscCall(MatSetValues(
1595  _cmc_pressure_force_mat, 1, &row_pf, 1, &col_pf, &value_pf, ADD_VALUES));
1596  }
1597  }
1598  else
1599  {
1600  PetscScalar pressure_factor = std::pow(Sij, 2.0) * rho_star;
1601  PetscInt row_pf = i_gap + _n_gaps * iz_ind;
1602  PetscInt col_pf = i_ch + _n_channels * iz_ind;
1603  PetscScalar value_pf = -1.0 * pressure_factor;
1604  LibmeshPetscCall(
1605  MatSetValues(_cmc_pressure_force_mat, 1, &row_pf, 1, &col_pf, &value_pf, ADD_VALUES));
1606  col_pf = j_ch + _n_channels * iz_ind;
1607  value_pf = pressure_factor;
1608  LibmeshPetscCall(
1609  MatSetValues(_cmc_pressure_force_mat, 1, &row_pf, 1, &col_pf, &value_pf, ADD_VALUES));
1610  }
1611  }
1612  }
1614  LibmeshPetscCall(MatZeroEntries(_cmc_sys_Wij_mat));
1615  LibmeshPetscCall(VecZeroEntries(_cmc_sys_Wij_rhs));
1616  LibmeshPetscCall(MatAssemblyBegin(_cmc_time_derivative_mat, MAT_FINAL_ASSEMBLY));
1617  LibmeshPetscCall(MatAssemblyEnd(_cmc_time_derivative_mat, MAT_FINAL_ASSEMBLY));
1618  LibmeshPetscCall(MatAssemblyBegin(_cmc_advective_derivative_mat, MAT_FINAL_ASSEMBLY));
1619  LibmeshPetscCall(MatAssemblyEnd(_cmc_advective_derivative_mat, MAT_FINAL_ASSEMBLY));
1620  LibmeshPetscCall(MatAssemblyBegin(_cmc_friction_force_mat, MAT_FINAL_ASSEMBLY));
1621  LibmeshPetscCall(MatAssemblyEnd(_cmc_friction_force_mat, MAT_FINAL_ASSEMBLY));
1622  LibmeshPetscCall(MatAssemblyBegin(_cmc_pressure_force_mat, MAT_FINAL_ASSEMBLY));
1623  LibmeshPetscCall(MatAssemblyEnd(_cmc_pressure_force_mat, MAT_FINAL_ASSEMBLY));
1624  LibmeshPetscCall(MatAssemblyBegin(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1625  LibmeshPetscCall(MatAssemblyEnd(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1626  // Matrix
1627 #if !PETSC_VERSION_LESS_THAN(3, 15, 0)
1628  LibmeshPetscCall(
1629  MatAXPY(_cmc_sys_Wij_mat, 1.0, _cmc_time_derivative_mat, UNKNOWN_NONZERO_PATTERN));
1630  LibmeshPetscCall(MatAssemblyBegin(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1631  LibmeshPetscCall(MatAssemblyEnd(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1632  LibmeshPetscCall(
1633  MatAXPY(_cmc_sys_Wij_mat, 1.0, _cmc_advective_derivative_mat, UNKNOWN_NONZERO_PATTERN));
1634  LibmeshPetscCall(MatAssemblyBegin(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1635  LibmeshPetscCall(MatAssemblyEnd(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1636  LibmeshPetscCall(
1637  MatAXPY(_cmc_sys_Wij_mat, 1.0, _cmc_friction_force_mat, UNKNOWN_NONZERO_PATTERN));
1638 #else
1639  LibmeshPetscCall(
1640  MatAXPY(_cmc_sys_Wij_mat, 1.0, _cmc_time_derivative_mat, DIFFERENT_NONZERO_PATTERN));
1641  LibmeshPetscCall(MatAssemblyBegin(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1642  LibmeshPetscCall(MatAssemblyEnd(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1643  LibmeshPetscCall(
1644  MatAXPY(_cmc_sys_Wij_mat, 1.0, _cmc_advective_derivative_mat, DIFFERENT_NONZERO_PATTERN));
1645  LibmeshPetscCall(MatAssemblyBegin(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1646  LibmeshPetscCall(MatAssemblyEnd(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1647  LibmeshPetscCall(
1648  MatAXPY(_cmc_sys_Wij_mat, 1.0, _cmc_friction_force_mat, DIFFERENT_NONZERO_PATTERN));
1649 #endif
1650  LibmeshPetscCall(MatAssemblyBegin(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1651  LibmeshPetscCall(MatAssemblyEnd(_cmc_sys_Wij_mat, MAT_FINAL_ASSEMBLY));
1652  // RHS
1653  LibmeshPetscCall(VecAXPY(_cmc_sys_Wij_rhs, 1.0, _cmc_time_derivative_rhs));
1654  LibmeshPetscCall(VecAXPY(_cmc_sys_Wij_rhs, 1.0, _cmc_advective_derivative_rhs));
1655  LibmeshPetscCall(VecAXPY(_cmc_sys_Wij_rhs, 1.0, _cmc_friction_force_rhs));
1656 
1657  if (_segregated_bool)
1658  {
1659  // Assembly the matrix system
1660  Vec sol_holder_P;
1661  LibmeshPetscCall(createPetscVector(sol_holder_P, _block_size * _n_gaps));
1662  Vec sol_holder_W;
1663  LibmeshPetscCall(createPetscVector(sol_holder_W, _block_size * _n_gaps));
1664  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
1665  _prodp, *_P_soln, iblock * _block_size, (iblock + 1) * _block_size - 1, _n_channels));
1666  LibmeshPetscCall(MatMult(_cmc_sys_Wij_mat, _Wij_vec, sol_holder_W));
1667  LibmeshPetscCall(VecAXPY(sol_holder_W, -1.0, _cmc_sys_Wij_rhs));
1668  LibmeshPetscCall(MatMult(_cmc_pressure_force_mat, _prodp, sol_holder_P));
1669  LibmeshPetscCall(VecAXPY(sol_holder_P, -1.0, _cmc_pressure_force_rhs));
1670  LibmeshPetscCall(VecAXPY(sol_holder_W, 1.0, sol_holder_P));
1671  PetscScalar * xx;
1672  LibmeshPetscCall(VecGetArray(sol_holder_W, &xx));
1673  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1674  {
1675  auto iz_ind = iz - first_node;
1676  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
1677  {
1678  _Wij_residual_matrix(i_gap, iz - 1 - iblock * _block_size) = xx[iz_ind * _n_gaps + i_gap];
1679  }
1680  }
1681  LibmeshPetscCall(VecDestroy(&sol_holder_P));
1682  LibmeshPetscCall(VecDestroy(&sol_holder_W));
1683  }
1684  }
1685 }
1686 
1687 void
1689 {
1690  unsigned int last_node = (iblock + 1) * _block_size;
1691  unsigned int first_node = iblock * _block_size + 1;
1692  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1693  {
1694  auto dz = _z_grid[iz] - _z_grid[iz - 1];
1695  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
1696  {
1697  auto chans = _subchannel_mesh.getGapChannels(i_gap);
1698  unsigned int i_ch = chans.first;
1699  unsigned int j_ch = chans.second;
1700  auto * node_in_i = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1701  auto * node_out_i = _subchannel_mesh.getChannelNode(i_ch, iz);
1702  auto * node_in_j = _subchannel_mesh.getChannelNode(j_ch, iz - 1);
1703  auto * node_out_j = _subchannel_mesh.getChannelNode(j_ch, iz);
1704  auto Si_in = (*_S_flow_soln)(node_in_i);
1705  auto Sj_in = (*_S_flow_soln)(node_in_j);
1706  auto Si_out = (*_S_flow_soln)(node_out_i);
1707  auto Sj_out = (*_S_flow_soln)(node_out_j);
1708  auto gap = _subchannel_mesh.getGapWidth(iz, i_gap);
1709  auto Sij = dz * gap;
1710  auto avg_massflux =
1711  0.5 * (((*_mdot_soln)(node_in_i) + (*_mdot_soln)(node_in_j)) / (Si_in + Sj_in) +
1712  ((*_mdot_soln)(node_out_i) + (*_mdot_soln)(node_out_j)) / (Si_out + Sj_out));
1713  auto beta = computeBeta(i_gap, iz, /*enthalpy=*/false);
1714 
1715  if (!_implicit_bool)
1716  {
1717  _WijPrime(i_gap, iz) = beta * avg_massflux * Sij;
1718  }
1719  else
1720  {
1721  auto iz_ind = iz - first_node;
1722  PetscScalar base_value = beta * 0.5 * Sij;
1723 
1724  // Bottom values
1725  if (iz == first_node)
1726  {
1727  PetscScalar value_tl = -1.0 * base_value / (Si_in + Sj_in) *
1728  ((*_mdot_soln)(node_in_i) + (*_mdot_soln)(node_in_j));
1729  PetscInt row = i_gap + _n_gaps * iz_ind;
1730  LibmeshPetscCall(
1731  VecSetValues(_amc_turbulent_cross_flows_rhs, 1, &row, &value_tl, INSERT_VALUES));
1732  }
1733  else
1734  {
1735  PetscScalar value_tl = base_value / (Si_in + Sj_in);
1736  PetscInt row = i_gap + _n_gaps * iz_ind;
1737 
1738  PetscInt col_ich = i_ch + _n_channels * (iz_ind - 1);
1739  LibmeshPetscCall(MatSetValues(
1740  _amc_turbulent_cross_flows_mat, 1, &row, 1, &col_ich, &value_tl, INSERT_VALUES));
1741 
1742  PetscInt col_jch = j_ch + _n_channels * (iz_ind - 1);
1743  LibmeshPetscCall(MatSetValues(
1744  _amc_turbulent_cross_flows_mat, 1, &row, 1, &col_jch, &value_tl, INSERT_VALUES));
1745  }
1746 
1747  // Top values
1748  PetscScalar value_bl = base_value / (Si_out + Sj_out);
1749  PetscInt row = i_gap + _n_gaps * iz_ind;
1750 
1751  PetscInt col_ich = i_ch + _n_channels * iz_ind;
1752  LibmeshPetscCall(MatSetValues(
1753  _amc_turbulent_cross_flows_mat, 1, &row, 1, &col_ich, &value_bl, INSERT_VALUES));
1754 
1755  PetscInt col_jch = j_ch + _n_channels * iz_ind;
1756  LibmeshPetscCall(MatSetValues(
1757  _amc_turbulent_cross_flows_mat, 1, &row, 1, &col_jch, &value_bl, INSERT_VALUES));
1758  }
1759  }
1760  }
1761 
1762  if (_implicit_bool)
1763  {
1764  LibmeshPetscCall(MatAssemblyBegin(_amc_turbulent_cross_flows_mat, MAT_FINAL_ASSEMBLY));
1765  LibmeshPetscCall(MatAssemblyEnd(_amc_turbulent_cross_flows_mat, MAT_FINAL_ASSEMBLY));
1766 
1768  Vec loc_prod;
1769  Vec loc_Wij;
1770  LibmeshPetscCall(VecDuplicate(_amc_sys_mdot_rhs, &loc_prod));
1771  LibmeshPetscCall(VecDuplicate(_Wij_vec, &loc_Wij));
1772  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
1773  loc_prod, *_mdot_soln, first_node, last_node, _n_channels));
1774  LibmeshPetscCall(MatMult(_amc_turbulent_cross_flows_mat, loc_prod, loc_Wij));
1775  LibmeshPetscCall(VecAXPY(loc_Wij, -1.0, _amc_turbulent_cross_flows_rhs));
1777  loc_Wij, _WijPrime, first_node, last_node, _n_gaps));
1778  LibmeshPetscCall(VecDestroy(&loc_prod));
1779  LibmeshPetscCall(VecDestroy(&loc_Wij));
1780  }
1781 }
1782 
1785 {
1786  unsigned int last_node = (iblock + 1) * _block_size;
1787  unsigned int first_node = iblock * _block_size + 1;
1788  libMesh::DenseVector<Real> Wij_residual_vector(_n_gaps * _block_size, 0.0);
1789  // Assign the solution to the cross-flow matrix
1790  int i = 0;
1791  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1792  {
1793  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
1794  {
1795  _Wij(i_gap, iz) = solution(i);
1796  i++;
1797  }
1798  }
1799 
1800  // Calculating sum of crossflows
1801  computeSumWij(iblock);
1802  // Solving axial flux
1803  computeMdot(iblock);
1804  // Calculation of turbulent Crossflow
1805  computeWijPrime(iblock);
1806  // Solving for Pressure Drop
1807  computeDP(iblock);
1808  // Solving for pressure
1809  computeP(iblock);
1810  // Populating lateral crossflow residual matrix
1811  computeWijResidual(iblock);
1812 
1813  // Turn the residual matrix into a residual vector
1814  for (unsigned int iz = 0; iz < _block_size; iz++)
1815  {
1816  for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
1817  {
1818  int i = _n_gaps * iz + i_gap; // column wise transfer
1819  Wij_residual_vector(i) = _Wij_residual_matrix(i_gap, iz);
1820  }
1821  }
1822  return Wij_residual_vector;
1823 }
1824 
1825 PetscErrorCode
1827  const libMesh::DenseVector<Real> & solution,
1829 {
1830  SNES snes;
1831  KSP ksp;
1832  PC pc;
1833  Vec x, r;
1834  PetscMPIInt size;
1835  PetscScalar * xx;
1836 
1838  PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1839  if (size > 1)
1840  SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Example is only for sequential runs");
1841  LibmeshPetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
1842  LibmeshPetscCall(VecCreate(PETSC_COMM_WORLD, &x));
1843  LibmeshPetscCall(VecSetSizes(x, PETSC_DECIDE, _block_size * _n_gaps));
1844  LibmeshPetscCall(VecSetFromOptions(x));
1845  LibmeshPetscCall(VecDuplicate(x, &r));
1846 
1847 #if PETSC_VERSION_LESS_THAN(3, 13, 0)
1848  LibmeshPetscCall(PetscOptionsSetValue(PETSC_NULL, "-snes_mf", PETSC_NULL));
1849 #else
1850  LibmeshPetscCall(SNESSetUseMatrixFree(snes, PETSC_FALSE, PETSC_TRUE));
1851 #endif
1852  Ctx ctx;
1853  ctx.iblock = iblock;
1854  ctx.schp = this;
1855  LibmeshPetscCall(SNESSetFunction(snes, r, formFunction, &ctx));
1856  LibmeshPetscCall(SNESGetKSP(snes, &ksp));
1857  LibmeshPetscCall(KSPGetPC(ksp, &pc));
1858  LibmeshPetscCall(PCSetType(pc, PCNONE));
1859  LibmeshPetscCall(KSPSetTolerances(ksp, _rtol, _atol, _dtol, _maxit));
1860  LibmeshPetscCall(SNESSetFromOptions(snes));
1861  LibmeshPetscCall(VecGetArray(x, &xx));
1862  for (unsigned int i = 0; i < _block_size * _n_gaps; i++)
1863  {
1864  xx[i] = solution(i);
1865  }
1866  LibmeshPetscCall(VecRestoreArray(x, &xx));
1867 
1868  LibmeshPetscCall(SNESSolve(snes, NULL, x));
1869  LibmeshPetscCall(VecGetArray(x, &xx));
1870  for (unsigned int i = 0; i < _block_size * _n_gaps; i++)
1871  root(i) = xx[i];
1872 
1873  LibmeshPetscCall(VecRestoreArray(x, &xx));
1874  LibmeshPetscCall(VecDestroy(&x));
1875  LibmeshPetscCall(VecDestroy(&r));
1876  LibmeshPetscCall(SNESDestroy(&snes));
1877  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
1878 }
1879 
1880 Real
1881 SubChannel1PhaseProblem::computeAddedHeatDuct(unsigned int i_ch, unsigned int iz)
1882 {
1883  mooseAssert(iz > 0, "Trapezoidal rule requires starting at index 1 at least");
1884  if (_duct_mesh_exist)
1885  {
1886  auto subch_type = _subchannel_mesh.getSubchannelType(i_ch);
1887  if (subch_type == EChannelType::EDGE || subch_type == EChannelType::CORNER)
1888  {
1889  auto dz = _z_grid[iz] - _z_grid[iz - 1];
1890  auto * node_in_chan = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
1891  auto * node_out_chan = _subchannel_mesh.getChannelNode(i_ch, iz);
1892  auto * node_in_duct = _subchannel_mesh.getDuctNodeFromChannel(node_in_chan);
1893  auto * node_out_duct = _subchannel_mesh.getDuctNodeFromChannel(node_out_chan);
1894  auto heat_rate_in = (*_duct_heat_flux_soln)(node_in_duct);
1895  auto heat_rate_out = (*_duct_heat_flux_soln)(node_out_duct);
1896  auto width = getSubChannelPeripheralDuctWidth(i_ch);
1897  return 0.5 * (heat_rate_in + heat_rate_out) * dz * width;
1898  }
1899  else
1900  {
1901  return 0.0;
1902  }
1903  }
1904  else
1905  {
1906  return 0.0;
1907  }
1908 }
1909 
1910 PetscErrorCode
1912 {
1913  bool lag_block_thermal_solve = true;
1914  Vec b_nest, x_nest; /* approx solution, RHS, exact solution */
1915  Mat A_nest; /* linear system matrix */
1916  KSP ksp; /* linear solver context */
1917  PC pc; /* preconditioner context */
1918 
1920  PetscInt Q = _monolithic_thermal_bool ? 4 : 3;
1921  std::vector<Mat> mat_array(Q * Q);
1922  std::vector<Vec> vec_array(Q);
1923 
1925  bool _axial_mass_flow_tight_coupling = true;
1926  bool _pressure_axial_momentum_tight_coupling = true;
1927  bool _pressure_cross_momentum_tight_coupling = true;
1928  unsigned int first_node = iblock * _block_size + 1;
1929  unsigned int last_node = (iblock + 1) * _block_size;
1930 
1932  // Computing sum of crossflows with previous iteration
1933  computeSumWij(iblock);
1934  // Assembling axial flux matrix
1935  computeMdot(iblock);
1936  // Computing turbulent crossflow with previous step axial mass flows
1937  computeWijPrime(iblock);
1938  // Assembling for Pressure Drop matrix
1939  computeDP(iblock);
1940  // Assembling pressure matrix
1941  computeP(iblock);
1942  // Assembling cross fluxes matrix
1943  computeWijResidual(iblock);
1944  // If monolithic solve - Assembling enthalpy matrix
1946  computeh(iblock);
1947 
1948  if (_verbose_subchannel)
1949  {
1950  _console << "Starting nested system." << std::endl;
1951  _console << "Number of simultaneous variables: " << Q << std::endl;
1952  }
1953  // Mass conservation
1954  PetscInt field_num = 0;
1955  LibmeshPetscCall(
1956  MatDuplicate(_mc_axial_convection_mat, MAT_COPY_VALUES, &mat_array[Q * field_num + 0]));
1957  LibmeshPetscCall(MatAssemblyBegin(mat_array[Q * field_num + 0], MAT_FINAL_ASSEMBLY));
1958  LibmeshPetscCall(MatAssemblyEnd(mat_array[Q * field_num + 0], MAT_FINAL_ASSEMBLY));
1959  mat_array[Q * field_num + 1] = NULL;
1960  if (_axial_mass_flow_tight_coupling)
1961  {
1962  LibmeshPetscCall(MatDuplicate(_mc_sumWij_mat, MAT_COPY_VALUES, &mat_array[Q * field_num + 2]));
1963  LibmeshPetscCall(MatAssemblyBegin(mat_array[Q * field_num + 2], MAT_FINAL_ASSEMBLY));
1964  LibmeshPetscCall(MatAssemblyEnd(mat_array[Q * field_num + 2], MAT_FINAL_ASSEMBLY));
1965  }
1966  else
1967  {
1968  mat_array[Q * field_num + 2] = NULL;
1969  }
1971  {
1972  mat_array[Q * field_num + 3] = NULL;
1973  }
1974  LibmeshPetscCall(VecDuplicate(_mc_axial_convection_rhs, &vec_array[field_num]));
1975  LibmeshPetscCall(VecCopy(_mc_axial_convection_rhs, vec_array[field_num]));
1976  if (!_axial_mass_flow_tight_coupling)
1977  {
1978  Vec sumWij_loc;
1979  LibmeshPetscCall(VecDuplicate(_mc_axial_convection_rhs, &sumWij_loc));
1980  LibmeshPetscCall(VecSet(sumWij_loc, 0.0));
1981  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
1982  {
1983  auto iz_ind = iz - first_node;
1984  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
1985  {
1986  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
1987 
1988  PetscScalar value_vec_2 = -1.0 * (*_SumWij_soln)(node_out);
1989  PetscInt row_vec_2 = i_ch + _n_channels * iz_ind;
1990  LibmeshPetscCall(VecSetValues(sumWij_loc, 1, &row_vec_2, &value_vec_2, ADD_VALUES));
1991  }
1992  }
1993  LibmeshPetscCall(VecAXPY(vec_array[field_num], 1.0, sumWij_loc));
1994  LibmeshPetscCall(VecDestroy(&sumWij_loc));
1995  }
1996 
1997  if (_verbose_subchannel)
1998  _console << "Mass ok." << std::endl;
1999  // Axial momentum conservation
2000  field_num = 1;
2001  if (_pressure_axial_momentum_tight_coupling)
2002  {
2003  LibmeshPetscCall(
2004  MatDuplicate(_amc_sys_mdot_mat, MAT_COPY_VALUES, &mat_array[Q * field_num + 0]));
2005  LibmeshPetscCall(MatAssemblyBegin(mat_array[Q * field_num + 0], MAT_FINAL_ASSEMBLY));
2006  LibmeshPetscCall(MatAssemblyEnd(mat_array[Q * field_num + 0], MAT_FINAL_ASSEMBLY));
2007  }
2008  else
2009  {
2010  mat_array[Q * field_num + 0] = NULL;
2011  }
2012  LibmeshPetscCall(
2013  MatDuplicate(_amc_pressure_force_mat, MAT_COPY_VALUES, &mat_array[Q * field_num + 1]));
2014  LibmeshPetscCall(MatAssemblyBegin(mat_array[Q * field_num + 1], MAT_FINAL_ASSEMBLY));
2015  LibmeshPetscCall(MatAssemblyEnd(mat_array[Q * field_num + 1], MAT_FINAL_ASSEMBLY));
2016  mat_array[Q * field_num + 2] = NULL;
2018  {
2019  mat_array[Q * field_num + 3] = NULL;
2020  }
2021  LibmeshPetscCall(VecDuplicate(_amc_pressure_force_rhs, &vec_array[field_num]));
2022  LibmeshPetscCall(VecCopy(_amc_pressure_force_rhs, vec_array[field_num]));
2023  if (_pressure_axial_momentum_tight_coupling)
2024  {
2025  LibmeshPetscCall(VecAXPY(vec_array[field_num], 1.0, _amc_sys_mdot_rhs));
2026  }
2027  else
2028  {
2029  unsigned int last_node = (iblock + 1) * _block_size;
2030  unsigned int first_node = iblock * _block_size + 1;
2031  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
2032  _prod, *_mdot_soln, first_node, last_node, _n_channels));
2033  Vec ls;
2034  LibmeshPetscCall(VecDuplicate(_amc_sys_mdot_rhs, &ls));
2035  LibmeshPetscCall(MatMult(_amc_sys_mdot_mat, _prod, ls));
2036  LibmeshPetscCall(VecAXPY(ls, -1.0, _amc_sys_mdot_rhs));
2037  LibmeshPetscCall(VecAXPY(vec_array[field_num], -1.0, ls));
2038  LibmeshPetscCall(VecDestroy(&ls));
2039  }
2040 
2041  if (_verbose_subchannel)
2042  _console << "Lin mom OK." << std::endl;
2043 
2044  // Cross momentum conservation
2045  field_num = 2;
2046  mat_array[Q * field_num + 0] = NULL;
2047  if (_pressure_cross_momentum_tight_coupling)
2048  {
2049  LibmeshPetscCall(
2050  MatDuplicate(_cmc_pressure_force_mat, MAT_COPY_VALUES, &mat_array[Q * field_num + 1]));
2051  LibmeshPetscCall(MatAssemblyBegin(mat_array[Q * field_num + 1], MAT_FINAL_ASSEMBLY));
2052  LibmeshPetscCall(MatAssemblyEnd(mat_array[Q * field_num + 1], MAT_FINAL_ASSEMBLY));
2053  }
2054  else
2055  {
2056  mat_array[Q * field_num + 1] = NULL;
2057  }
2058 
2059  LibmeshPetscCall(MatDuplicate(_cmc_sys_Wij_mat, MAT_COPY_VALUES, &mat_array[Q * field_num + 2]));
2060  LibmeshPetscCall(MatAssemblyBegin(mat_array[Q * field_num + 2], MAT_FINAL_ASSEMBLY));
2061  LibmeshPetscCall(MatAssemblyEnd(mat_array[Q * field_num + 2], MAT_FINAL_ASSEMBLY));
2063  {
2064  mat_array[Q * field_num + 3] = NULL;
2065  }
2066 
2067  LibmeshPetscCall(VecDuplicate(_cmc_sys_Wij_rhs, &vec_array[field_num]));
2068  LibmeshPetscCall(VecCopy(_cmc_sys_Wij_rhs, vec_array[field_num]));
2069  if (_pressure_cross_momentum_tight_coupling)
2070  {
2071  LibmeshPetscCall(VecAXPY(vec_array[field_num], 1.0, _cmc_pressure_force_rhs));
2072  }
2073  else
2074  {
2075  Vec sol_holder_P;
2076  LibmeshPetscCall(createPetscVector(sol_holder_P, _block_size * _n_gaps));
2077  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
2078  _prodp, *_P_soln, iblock * _block_size, (iblock + 1) * _block_size - 1, _n_channels));
2079 
2080  LibmeshPetscCall(MatMult(_cmc_pressure_force_mat, _prodp, sol_holder_P));
2081  LibmeshPetscCall(VecAXPY(sol_holder_P, -1.0, _cmc_pressure_force_rhs));
2082  LibmeshPetscCall(VecScale(sol_holder_P, 1.0));
2083  LibmeshPetscCall(VecAXPY(vec_array[field_num], 1.0, sol_holder_P));
2084  LibmeshPetscCall(VecDestroy(&sol_holder_P));
2085  }
2086 
2087  if (_verbose_subchannel)
2088  _console << "Cross mom ok." << std::endl;
2089 
2090  // Energy conservation
2092  {
2093  field_num = 3;
2094  mat_array[Q * field_num + 0] = NULL;
2095  mat_array[Q * field_num + 1] = NULL;
2096  mat_array[Q * field_num + 2] = NULL;
2097  LibmeshPetscCall(MatDuplicate(_hc_sys_h_mat, MAT_COPY_VALUES, &mat_array[Q * field_num + 3]));
2098  if (lag_block_thermal_solve)
2099  {
2100  LibmeshPetscCall(MatZeroEntries(mat_array[Q * field_num + 3]));
2101  LibmeshPetscCall(MatShift(mat_array[Q * field_num + 3], 1.0));
2102  }
2103  LibmeshPetscCall(MatAssemblyBegin(mat_array[Q * field_num + 3], MAT_FINAL_ASSEMBLY));
2104  LibmeshPetscCall(MatAssemblyEnd(mat_array[Q * field_num + 3], MAT_FINAL_ASSEMBLY));
2105  LibmeshPetscCall(VecDuplicate(_hc_sys_h_rhs, &vec_array[field_num]));
2106  LibmeshPetscCall(VecCopy(_hc_sys_h_rhs, vec_array[field_num]));
2107  if (lag_block_thermal_solve)
2108  {
2109  LibmeshPetscCall(VecZeroEntries(vec_array[field_num]));
2110  LibmeshPetscCall(VecShift(vec_array[field_num], 1.0));
2111  }
2112 
2113  if (_verbose_subchannel)
2114  _console << "Energy ok." << std::endl;
2115  }
2116 
2117  // Relaxing linear system
2118  // Weaker relaxation
2119  if (true)
2120  {
2121  // Estimating cross-flow resistances to achieve realizable solves
2122  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
2123  _prod, *_mdot_soln, first_node, last_node, _n_channels));
2124  Vec mdot_estimate;
2125  LibmeshPetscCall(createPetscVector(mdot_estimate, _block_size * _n_channels));
2126  Vec pmat_diag;
2127  LibmeshPetscCall(createPetscVector(pmat_diag, _block_size * _n_channels));
2128  Vec p_estimate;
2129  LibmeshPetscCall(createPetscVector(p_estimate, _block_size * _n_channels));
2130  Vec unity_vec;
2131  LibmeshPetscCall(createPetscVector(unity_vec, _block_size * _n_channels));
2132  LibmeshPetscCall(VecSet(unity_vec, 1.0));
2133  Vec sol_holder_P;
2134  LibmeshPetscCall(createPetscVector(sol_holder_P, _block_size * _n_gaps));
2135  Vec diag_Wij_loc;
2136  LibmeshPetscCall(createPetscVector(diag_Wij_loc, _block_size * _n_gaps));
2137  Vec Wij_estimate;
2138  LibmeshPetscCall(createPetscVector(Wij_estimate, _block_size * _n_gaps));
2139  Vec unity_vec_Wij;
2140  LibmeshPetscCall(createPetscVector(unity_vec_Wij, _block_size * _n_gaps));
2141  LibmeshPetscCall(VecSet(unity_vec_Wij, 1.0));
2142  Vec _Wij_loc_vec;
2143  LibmeshPetscCall(createPetscVector(_Wij_loc_vec, _block_size * _n_gaps));
2144  Vec _Wij_old_loc_vec;
2145  LibmeshPetscCall(createPetscVector(_Wij_old_loc_vec, _block_size * _n_gaps));
2146  LibmeshPetscCall(MatMult(mat_array[Q], _prod, mdot_estimate));
2147  LibmeshPetscCall(MatGetDiagonal(mat_array[Q + 1], pmat_diag));
2148  LibmeshPetscCall(VecAXPY(pmat_diag, 1e-10, unity_vec));
2149  LibmeshPetscCall(VecPointwiseDivide(p_estimate, mdot_estimate, pmat_diag));
2150  LibmeshPetscCall(MatMult(mat_array[2 * Q + 1], p_estimate, sol_holder_P));
2151  LibmeshPetscCall(VecAXPY(sol_holder_P, -1.0, _cmc_pressure_force_rhs));
2152  LibmeshPetscCall(MatGetDiagonal(mat_array[2 * Q + 2], diag_Wij_loc));
2153  LibmeshPetscCall(VecAXPY(diag_Wij_loc, 1e-10, unity_vec_Wij));
2154  LibmeshPetscCall(VecPointwiseDivide(Wij_estimate, sol_holder_P, diag_Wij_loc));
2155  Vec sumWij_loc;
2156  LibmeshPetscCall(createPetscVector(sumWij_loc, _block_size * _n_channels));
2157  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
2158  {
2159  auto iz_ind = iz - first_node;
2160  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
2161  {
2162  PetscScalar sumWij = 0.0;
2163  unsigned int counter = 0;
2164  for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
2165  {
2166  auto chans = _subchannel_mesh.getGapChannels(i_gap);
2167  unsigned int i_ch_loc = chans.first;
2168  PetscInt row_vec = i_ch_loc + _n_channels * iz_ind;
2169  PetscScalar loc_Wij_value;
2170  LibmeshPetscCall(VecGetValues(sol_holder_P, 1, &row_vec, &loc_Wij_value));
2171  sumWij += _subchannel_mesh.getCrossflowSign(i_ch, counter) * loc_Wij_value;
2172  counter++;
2173  }
2174  PetscInt row_vec = i_ch + _n_channels * iz_ind;
2175  LibmeshPetscCall(VecSetValues(sumWij_loc, 1, &row_vec, &sumWij, INSERT_VALUES));
2176  }
2177  }
2178 
2179  PetscScalar min_mdot;
2180  LibmeshPetscCall(VecAbs(_prod));
2181  LibmeshPetscCall(VecMin(_prod, NULL, &min_mdot));
2182  if (_verbose_subchannel)
2183  _console << "Minimum estimated mdot: " << min_mdot << std::endl;
2184 
2185  LibmeshPetscCall(VecAbs(sumWij_loc));
2186  LibmeshPetscCall(VecMax(sumWij_loc, NULL, &_max_sumWij));
2187  _max_sumWij = std::max(1e-10, _max_sumWij);
2188  if (_verbose_subchannel)
2189  _console << "Maximum estimated Wij: " << _max_sumWij << std::endl;
2190 
2192  _Wij_loc_vec, _Wij, first_node, last_node, _n_gaps));
2193  LibmeshPetscCall(VecAbs(_Wij_loc_vec));
2195  _Wij_old_loc_vec, _Wij_old, first_node, last_node, _n_gaps));
2196  LibmeshPetscCall(VecAbs(_Wij_old_loc_vec));
2197  LibmeshPetscCall(VecAXPY(_Wij_loc_vec, -1.0, _Wij_old_loc_vec));
2198  PetscScalar relax_factor;
2199  LibmeshPetscCall(VecAbs(_Wij_loc_vec));
2200 #if !PETSC_VERSION_LESS_THAN(3, 16, 0)
2201  LibmeshPetscCall(VecMean(_Wij_loc_vec, &relax_factor));
2202 #else
2203  VecSum(_Wij_loc_vec, &relax_factor);
2204  relax_factor /= _block_size * _n_gaps;
2205 #endif
2206  relax_factor = relax_factor / _max_sumWij + 0.5;
2207  if (_verbose_subchannel)
2208  _console << "Relax base value: " << relax_factor << std::endl;
2209 
2210  PetscScalar resistance_relaxation = 0.9;
2211  _added_K = _max_sumWij / min_mdot;
2212  if (_verbose_subchannel)
2213  _console << "New cross resistance: " << _added_K << std::endl;
2214  _added_K = (_added_K * resistance_relaxation + (1.0 - resistance_relaxation) * _added_K_old) *
2215  relax_factor;
2216  if (_verbose_subchannel)
2217  _console << "Relaxed cross resistance: " << _added_K << std::endl;
2218  if (_added_K < 10 && _added_K >= 1.0)
2219  _added_K = 1.0; //(1.0 - resistance_relaxation);
2220  if (_added_K < 1.0 && _added_K >= 0.1)
2221  _added_K = 0.5;
2222  if (_added_K < 0.1 && _added_K >= 0.01)
2223  _added_K = 1. / 3.;
2224  if (_added_K < 1e-2 && _added_K >= 1e-3)
2225  _added_K = 0.1;
2226  if (_added_K < 1e-3)
2227  _added_K = 1.0 * _added_K;
2228  if (_verbose_subchannel)
2229  _console << "Actual added cross resistance: " << _added_K << std::endl;
2230  LibmeshPetscCall(VecScale(unity_vec_Wij, _added_K));
2232 
2233  // Adding cross resistances
2234  LibmeshPetscCall(MatDiagonalSet(mat_array[2 * Q + 2], unity_vec_Wij, ADD_VALUES));
2235  LibmeshPetscCall(VecDestroy(&mdot_estimate));
2236  LibmeshPetscCall(VecDestroy(&pmat_diag));
2237  LibmeshPetscCall(VecDestroy(&unity_vec));
2238  LibmeshPetscCall(VecDestroy(&p_estimate));
2239  LibmeshPetscCall(VecDestroy(&sol_holder_P));
2240  LibmeshPetscCall(VecDestroy(&diag_Wij_loc));
2241  LibmeshPetscCall(VecDestroy(&unity_vec_Wij));
2242  LibmeshPetscCall(VecDestroy(&Wij_estimate));
2243  LibmeshPetscCall(VecDestroy(&sumWij_loc));
2244  LibmeshPetscCall(VecDestroy(&_Wij_loc_vec));
2245  LibmeshPetscCall(VecDestroy(&_Wij_old_loc_vec));
2246 
2247  // Auto-computing relaxation factors
2248  PetscScalar relaxation_factor_mdot, relaxation_factor_P, relaxation_factor_Wij;
2249  relaxation_factor_mdot = 1.0;
2250  relaxation_factor_P = 1.0; // std::exp(-5.0);
2251  relaxation_factor_Wij = 0.1;
2252 
2253  if (_verbose_subchannel)
2254  {
2255  _console << "Relax mdot: " << relaxation_factor_mdot << std::endl;
2256  _console << "Relax P: " << relaxation_factor_P << std::endl;
2257  _console << "Relax Wij: " << relaxation_factor_Wij << std::endl;
2258  }
2259 
2260  PetscInt field_num = 0;
2261  Vec diag_mdot;
2262  LibmeshPetscCall(VecDuplicate(vec_array[field_num], &diag_mdot));
2263  LibmeshPetscCall(MatGetDiagonal(mat_array[Q * field_num + field_num], diag_mdot));
2264  LibmeshPetscCall(VecScale(diag_mdot, 1.0 / relaxation_factor_mdot));
2265  LibmeshPetscCall(
2266  MatDiagonalSet(mat_array[Q * field_num + field_num], diag_mdot, INSERT_VALUES));
2267  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
2268  _prod, *_mdot_soln, first_node, last_node, _n_channels));
2269  LibmeshPetscCall(VecScale(diag_mdot, (1.0 - relaxation_factor_mdot)));
2270  LibmeshPetscCall(VecPointwiseMult(_prod, _prod, diag_mdot));
2271  LibmeshPetscCall(VecAXPY(vec_array[field_num], 1.0, _prod));
2272  LibmeshPetscCall(VecDestroy(&diag_mdot));
2273 
2274  if (_verbose_subchannel)
2275  _console << "mdot relaxed" << std::endl;
2276 
2277  field_num = 1;
2278  Vec diag_P;
2279  LibmeshPetscCall(VecDuplicate(vec_array[field_num], &diag_P));
2280  LibmeshPetscCall(MatGetDiagonal(mat_array[Q * field_num + field_num], diag_P));
2281  LibmeshPetscCall(VecScale(diag_P, 1.0 / relaxation_factor_P));
2282  LibmeshPetscCall(MatDiagonalSet(mat_array[Q * field_num + field_num], diag_P, INSERT_VALUES));
2283  if (_verbose_subchannel)
2284  _console << "Mat assembled" << std::endl;
2285  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
2286  _prod, *_P_soln, first_node, last_node, _n_channels));
2287  LibmeshPetscCall(VecScale(diag_P, (1.0 - relaxation_factor_P)));
2288  LibmeshPetscCall(VecPointwiseMult(_prod, _prod, diag_P));
2289  LibmeshPetscCall(VecAXPY(vec_array[field_num], 1.0, _prod));
2290  LibmeshPetscCall(VecDestroy(&diag_P));
2291 
2292  if (_verbose_subchannel)
2293  _console << "P relaxed" << std::endl;
2294 
2295  field_num = 2;
2296  Vec diag_Wij;
2297  LibmeshPetscCall(VecDuplicate(vec_array[field_num], &diag_Wij));
2298  LibmeshPetscCall(MatGetDiagonal(mat_array[Q * field_num + field_num], diag_Wij));
2299  LibmeshPetscCall(VecScale(diag_Wij, 1.0 / relaxation_factor_Wij));
2300  LibmeshPetscCall(MatDiagonalSet(mat_array[Q * field_num + field_num], diag_Wij, INSERT_VALUES));
2302  _Wij_vec, _Wij, first_node, last_node, _n_gaps));
2303  LibmeshPetscCall(VecScale(diag_Wij, (1.0 - relaxation_factor_Wij)));
2304  LibmeshPetscCall(VecPointwiseMult(_Wij_vec, _Wij_vec, diag_Wij));
2305  LibmeshPetscCall(VecAXPY(vec_array[field_num], 1.0, _Wij_vec));
2306  LibmeshPetscCall(VecDestroy(&diag_Wij));
2307 
2308  if (_verbose_subchannel)
2309  _console << "Wij relaxed" << std::endl;
2310  }
2311  if (_verbose_subchannel)
2312  _console << "Linear solver relaxed." << std::endl;
2313 
2314  // Creating nested matrices
2315  LibmeshPetscCall(MatCreateNest(PETSC_COMM_WORLD, Q, NULL, Q, NULL, mat_array.data(), &A_nest));
2316  LibmeshPetscCall(VecCreateNest(PETSC_COMM_WORLD, Q, NULL, vec_array.data(), &b_nest));
2317  if (_verbose_subchannel)
2318  _console << "Nested system created." << std::endl;
2319 
2321  // Creating linear solver
2322  LibmeshPetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
2323  LibmeshPetscCall(KSPSetType(ksp, KSPFGMRES));
2324  // Setting KSP operators
2325  LibmeshPetscCall(KSPSetOperators(ksp, A_nest, A_nest));
2326  // Set KSP and PC options
2327  LibmeshPetscCall(KSPGetPC(ksp, &pc));
2328  LibmeshPetscCall(PCSetType(pc, PCFIELDSPLIT));
2329  LibmeshPetscCall(KSPSetTolerances(ksp, _rtol, _atol, _dtol, _maxit));
2330  // Splitting fields
2331  std::vector<IS> rows(Q);
2332  // IS rows[Q];
2333  PetscInt M = 0;
2334  LibmeshPetscCall(MatNestGetISs(A_nest, rows.data(), NULL));
2335  for (PetscInt j = 0; j < Q; ++j)
2336  {
2337  IS expand1;
2338  LibmeshPetscCall(ISDuplicate(rows[M], &expand1));
2339  M += 1;
2340  LibmeshPetscCall(PCFieldSplitSetIS(pc, NULL, expand1));
2341  LibmeshPetscCall(ISDestroy(&expand1));
2342  }
2343  if (_verbose_subchannel)
2344  _console << "Linear solver assembled." << std::endl;
2345 
2347  LibmeshPetscCall(VecDuplicate(b_nest, &x_nest));
2348  LibmeshPetscCall(VecSet(x_nest, 0.0));
2349  LibmeshPetscCall(KSPSolve(ksp, b_nest, x_nest));
2350 
2352  LibmeshPetscCall(VecDestroy(&b_nest));
2353  LibmeshPetscCall(MatDestroy(&A_nest));
2354  LibmeshPetscCall(KSPDestroy(&ksp));
2355  for (PetscInt i = 0; i < Q * Q; i++)
2356  {
2357  LibmeshPetscCall(MatDestroy(&mat_array[i]));
2358  }
2359  for (PetscInt i = 0; i < Q; i++)
2360  {
2361  LibmeshPetscCall(VecDestroy(&vec_array[i]));
2362  }
2363  if (_verbose_subchannel)
2364  _console << "Solver elements destroyed." << std::endl;
2365 
2367  Vec sol_mdot, sol_p, sol_Wij;
2368  if (_verbose_subchannel)
2369  _console << "Vectors created." << std::endl;
2370  PetscInt num_vecs;
2371  Vec * loc_vecs;
2372  LibmeshPetscCall(VecNestGetSubVecs(x_nest, &num_vecs, &loc_vecs));
2373  if (_verbose_subchannel)
2374  _console << "Starting extraction." << std::endl;
2375  LibmeshPetscCall(VecDuplicate(_mc_axial_convection_rhs, &sol_mdot));
2376  LibmeshPetscCall(VecCopy(loc_vecs[0], sol_mdot));
2377  LibmeshPetscCall(VecDuplicate(_amc_sys_mdot_rhs, &sol_p));
2378  LibmeshPetscCall(VecCopy(loc_vecs[1], sol_p));
2379  LibmeshPetscCall(VecDuplicate(_cmc_sys_Wij_rhs, &sol_Wij));
2380  LibmeshPetscCall(VecCopy(loc_vecs[2], sol_Wij));
2381  if (_verbose_subchannel)
2382  _console << "Getting individual field solutions from coupled solver." << std::endl;
2383 
2385  PetscScalar * sol_p_array;
2386  LibmeshPetscCall(VecGetArray(sol_p, &sol_p_array));
2387  PetscScalar * sol_Wij_array;
2388  LibmeshPetscCall(VecGetArray(sol_Wij, &sol_Wij_array));
2389 
2391  LibmeshPetscCall(populateSolutionChan<SolutionHandle>(
2392  sol_mdot, *_mdot_soln, first_node, last_node, _n_channels));
2393 
2395  for (unsigned int iz = last_node; iz > first_node - 1; iz--)
2396  {
2397  auto iz_ind = iz - first_node;
2398  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
2399  {
2400  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
2401  PetscScalar value = sol_p_array[iz_ind * _n_channels + i_ch];
2402  _P_soln->set(node_in, value);
2403  }
2404  }
2405 
2408  sol_Wij, _Wij, first_node, last_node - 1, _n_gaps));
2409 
2412  {
2413  if (lag_block_thermal_solve)
2414  {
2415  KSP ksploc;
2416  PC pc;
2417  Vec sol;
2418  LibmeshPetscCall(VecDuplicate(_hc_sys_h_rhs, &sol));
2419  LibmeshPetscCall(KSPCreate(PETSC_COMM_WORLD, &ksploc));
2420  LibmeshPetscCall(KSPSetOperators(ksploc, _hc_sys_h_mat, _hc_sys_h_mat));
2421  LibmeshPetscCall(KSPGetPC(ksploc, &pc));
2422  LibmeshPetscCall(PCSetType(pc, PCJACOBI));
2423  LibmeshPetscCall(KSPSetTolerances(ksploc, _rtol, _atol, _dtol, _maxit));
2424  LibmeshPetscCall(KSPSetFromOptions(ksploc));
2425  LibmeshPetscCall(KSPSolve(ksploc, _hc_sys_h_rhs, sol));
2426  PetscScalar * xx;
2427  LibmeshPetscCall(VecGetArray(sol, &xx));
2428  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
2429  {
2430  auto iz_ind = iz - first_node;
2431  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
2432  {
2433  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
2434  auto h_out = xx[iz_ind * _n_channels + i_ch];
2435  if (h_out < 0)
2436  {
2437  mooseError(name(),
2438  " : Calculation of negative Enthalpy h_out = : ",
2439  h_out,
2440  " Axial Level= : ",
2441  iz);
2442  }
2443  _h_soln->set(node_out, h_out);
2444  }
2445  }
2446  LibmeshPetscCall(KSPDestroy(&ksploc));
2447  LibmeshPetscCall(VecDestroy(&sol));
2448  }
2449  else
2450  {
2451  Vec sol_h;
2452  LibmeshPetscCall(VecDuplicate(_hc_sys_h_rhs, &sol_h));
2453  LibmeshPetscCall(VecCopy(loc_vecs[3], sol_h));
2454  PetscScalar * sol_h_array;
2455  LibmeshPetscCall(VecGetArray(sol_h, &sol_h_array));
2456  for (unsigned int iz = first_node; iz < last_node + 1; iz++)
2457  {
2458  auto iz_ind = iz - first_node;
2459  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
2460  {
2461  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
2462  auto h_out = sol_h_array[iz_ind * _n_channels + i_ch];
2463  if (h_out < 0)
2464  {
2465  mooseError(name(),
2466  " : Calculation of negative Enthalpy h_out = : ",
2467  h_out,
2468  " Axial Level= : ",
2469  iz);
2470  }
2471  _h_soln->set(node_out, h_out);
2472  }
2473  }
2474  LibmeshPetscCall(VecDestroy(&sol_h));
2475  }
2476  }
2477 
2479  LibmeshPetscCall(MatMult(_mc_sumWij_mat, sol_Wij, _prod));
2480  LibmeshPetscCall(populateSolutionChan<SolutionHandle>(
2481  _prod, *_SumWij_soln, first_node, last_node, _n_channels));
2482 
2483  Vec sumWij_loc;
2484  LibmeshPetscCall(createPetscVector(sumWij_loc, _block_size * _n_channels));
2485  LibmeshPetscCall(populateVectorFromHandle<SolutionHandle>(
2486  _prod, *_SumWij_soln, first_node, last_node, _n_channels));
2487 
2488  LibmeshPetscCall(VecAbs(_prod));
2489  LibmeshPetscCall(VecMax(_prod, NULL, &_max_sumWij_new));
2490  if (_verbose_subchannel)
2491  _console << "Maximum estimated Wij new: " << _max_sumWij_new << std::endl;
2493  if (_verbose_subchannel)
2494  _console << "Correction factor: " << _correction_factor << std::endl;
2495  if (_verbose_subchannel)
2496  _console << "Solutions assigned to MOOSE variables." << std::endl;
2497 
2499  LibmeshPetscCall(VecDestroy(&x_nest));
2500  LibmeshPetscCall(VecDestroy(&sol_mdot));
2501  LibmeshPetscCall(VecDestroy(&sol_p));
2502  LibmeshPetscCall(VecDestroy(&sol_Wij));
2503  LibmeshPetscCall(VecDestroy(&sumWij_loc));
2504  if (_verbose_subchannel)
2505  _console << "Solutions destroyed." << std::endl;
2506 
2507  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
2508 }
2509 
2510 void
2512 {
2513  _console << "Executing subchannel solver\n";
2514  _dt = (isTransient() ? dt() : _one);
2515  _TR = isTransient();
2517  if (_verbose_subchannel)
2518  _console << "Solution initialized" << std::endl;
2519  Real P_error = 1.0;
2520  unsigned int P_it = 0;
2521  unsigned int P_it_max;
2522 
2523  if (_segregated_bool)
2524  P_it_max = 20 * _n_blocks;
2525  else
2526  P_it_max = 100;
2527 
2528  if ((_n_blocks == 1) && (_segregated_bool))
2529  P_it_max = 5;
2530 
2531  while ((P_error > _P_tol && P_it < P_it_max))
2532  {
2533  P_it += 1;
2534  if (P_it == P_it_max && _n_blocks != 1)
2535  {
2536  _console << "Reached maximum number of axial pressure iterations" << std::endl;
2537  _converged = false;
2538  }
2539  _console << "Solving Outer Iteration : " << P_it << std::endl;
2540  auto P_L2norm_old_axial = _P_soln->L2norm();
2541  for (unsigned int iblock = 0; iblock < _n_blocks; iblock++)
2542  {
2543  int last_level = (iblock + 1) * _block_size;
2544  int first_level = iblock * _block_size + 1;
2545  Real T_block_error = 1.0;
2546  auto T_it = 0;
2547  _console << "Solving Block: " << iblock << " From first level: " << first_level
2548  << " to last level: " << last_level << std::endl;
2549  while (T_block_error > _T_tol && T_it < _T_maxit)
2550  {
2551  T_it += 1;
2552  if (T_it == _T_maxit)
2553  {
2554  _console << "Reached maximum number of temperature iterations for block: " << iblock
2555  << std::endl;
2556  _converged = false;
2557  }
2558  auto T_L2norm_old_block = _T_soln->L2norm();
2559 
2560  if (_segregated_bool)
2561  {
2562  computeWijFromSolve(iblock);
2563 
2564  if (_compute_power)
2565  {
2566  computeh(iblock);
2567  computeT(iblock);
2568  }
2569  }
2570  else
2571  {
2572  LibmeshPetscCall(implicitPetscSolve(iblock));
2573  computeWijPrime(iblock);
2574  if (_verbose_subchannel)
2575  _console << "Done with main solve." << std::endl;
2577  {
2578  // Enthalpy is already solved from the monolithic solve
2579  computeT(iblock);
2580  }
2581  else
2582  {
2583  if (_verbose_subchannel)
2584  _console << "Starting thermal solve." << std::endl;
2585  if (_compute_power)
2586  {
2587  computeh(iblock);
2588  computeT(iblock);
2589  }
2590  if (_verbose_subchannel)
2591  _console << "Done with thermal solve." << std::endl;
2592  }
2593  }
2594 
2595  if (_verbose_subchannel)
2596  _console << "Start updating thermophysical properties." << std::endl;
2597 
2598  if (_compute_density)
2599  computeRho(iblock);
2600 
2601  if (_compute_viscosity)
2602  computeMu(iblock);
2603 
2604  if (_verbose_subchannel)
2605  _console << "Done updating thermophysical properties." << std::endl;
2606 
2607  // We must do a global assembly to make sure data is parallel consistent before we do things
2608  // like compute L2 norms
2609  _aux->solution().close();
2610 
2611  auto T_L2norm_new = _T_soln->L2norm();
2612  T_block_error =
2613  std::abs((T_L2norm_new - T_L2norm_old_block) / (T_L2norm_old_block + 1E-14));
2614  _console << "T_block_error: " << T_block_error << std::endl;
2615  }
2616  }
2617  auto P_L2norm_new_axial = _P_soln->L2norm();
2618  P_error =
2619  std::abs((P_L2norm_new_axial - P_L2norm_old_axial) / (P_L2norm_old_axial + _P_out + 1E-14));
2620  _console << "P_error :" << P_error << std::endl;
2621  if (_verbose_subchannel)
2622  {
2623  _console << "Iteration: " << P_it << std::endl;
2624  _console << "Maximum iterations: " << P_it_max << std::endl;
2625  }
2626  }
2627  // update old crossflow matrix
2628  _Wij_old = _Wij;
2629  _console << "Finished executing subchannel solver\n";
2630 
2632  if (_pin_mesh_exist)
2633  {
2634  _console << "Commencing calculation of Pin surface temperature \n";
2635  for (unsigned int i_pin = 0; i_pin < _n_pins; i_pin++)
2636  {
2637  for (unsigned int iz = 0; iz < _n_cells + 1; ++iz)
2638  {
2639  const auto * pin_node = _subchannel_mesh.getPinNode(i_pin, iz);
2640  Real sumTemp = 0.0;
2641  Real rod_counter = 0.0;
2642  // Calculate sum of pin surface temperatures that the channels around the pin see
2643  for (auto i_ch : _subchannel_mesh.getPinChannels(i_pin))
2644  {
2645  const auto * node = _subchannel_mesh.getChannelNode(i_ch, iz);
2646  auto mu = (*_mu_soln)(node);
2647  auto S = (*_S_flow_soln)(node);
2648  auto w_perim = (*_w_perim_soln)(node);
2649  auto Dh_i = 4.0 * S / w_perim;
2650  auto Re = (((*_mdot_soln)(node) / S) * Dh_i / mu);
2651  auto k = _fp->k_from_p_T((*_P_soln)(node) + _P_out, (*_T_soln)(node));
2652  auto cp = _fp->cp_from_p_T((*_P_soln)(node) + _P_out, (*_T_soln)(node));
2653  auto Pr = (*_mu_soln)(node)*cp / k;
2654  auto Nu = 0.023 * std::pow(Re, 0.8) * std::pow(Pr, 0.4);
2655  auto hw = Nu * k / Dh_i;
2656  if ((*_Dpin_soln)(pin_node) <= 0)
2657  mooseError("Dpin should not be null or negative when computing pin powers: ",
2658  (*_Dpin_soln)(pin_node));
2659  sumTemp +=
2660  (*_q_prime_soln)(pin_node) / ((*_Dpin_soln)(pin_node)*M_PI * hw) + (*_T_soln)(node);
2661  rod_counter += 1.0;
2662  }
2663  if (rod_counter > 0)
2664  _Tpin_soln->set(pin_node, sumTemp / rod_counter);
2665  else
2666  mooseError("Pin was not found for pin index: " + std::to_string(i_pin));
2667  }
2668  }
2669  }
2670 
2672  if (_duct_mesh_exist)
2673  {
2674  _console << "Commencing calculation of duct surface temperature " << std::endl;
2675  auto duct_nodes = _subchannel_mesh.getDuctNodes();
2676  for (Node * dn : duct_nodes)
2677  {
2678  auto * node_chan = _subchannel_mesh.getChannelNodeFromDuct(dn);
2679  auto mu = (*_mu_soln)(node_chan);
2680  auto S = (*_S_flow_soln)(node_chan);
2681  auto w_perim = (*_w_perim_soln)(node_chan);
2682  auto Dh_i = 4.0 * S / w_perim;
2683  auto Re = (((*_mdot_soln)(node_chan) / S) * Dh_i / mu);
2684  auto k = _fp->k_from_p_T((*_P_soln)(node_chan) + _P_out, (*_T_soln)(node_chan));
2685  auto cp = _fp->cp_from_p_T((*_P_soln)(node_chan) + _P_out, (*_T_soln)(node_chan));
2686  auto Pr = (*_mu_soln)(node_chan)*cp / k;
2688  auto Nu = 0.023 * std::pow(Re, 0.8) * std::pow(Pr, 0.4);
2689  auto hw = Nu * k / Dh_i;
2690  auto T_chan = (*_duct_heat_flux_soln)(dn) / hw + (*_T_soln)(node_chan);
2691  _Tduct_soln->set(dn, T_chan);
2692  }
2693  }
2694  _aux->solution().close();
2695  _aux->update();
2696 
2697  Real power_in = 0.0;
2698  Real power_out = 0.0;
2699  Real Total_surface_area = 0.0;
2700  Real Total_wetted_perimeter = 0.0;
2701  Real mass_flow_in = 0.0;
2702  Real mass_flow_out = 0.0;
2703  for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
2704  {
2705  auto * node_in = _subchannel_mesh.getChannelNode(i_ch, 0);
2706  auto * node_out = _subchannel_mesh.getChannelNode(i_ch, _n_cells);
2707  Total_surface_area += (*_S_flow_soln)(node_in);
2708  Total_wetted_perimeter += (*_w_perim_soln)(node_in);
2709  power_in += (*_mdot_soln)(node_in) * (*_h_soln)(node_in);
2710  power_out += (*_mdot_soln)(node_out) * (*_h_soln)(node_out);
2711  mass_flow_in += (*_mdot_soln)(node_in);
2712  mass_flow_out += (*_mdot_soln)(node_out);
2713  }
2714  auto h_bulk_out = power_out / mass_flow_out;
2715  auto T_bulk_out = _fp->T_from_p_h(_P_out, h_bulk_out);
2716 
2717  Real bulk_Dh = 4.0 * Total_surface_area / Total_wetted_perimeter;
2718  Real inlet_mu = (*_mu_soln)(_subchannel_mesh.getChannelNode(0, 0));
2719  Real bulk_Re = mass_flow_in * bulk_Dh / (inlet_mu * Total_surface_area);
2720  if (_verbose_subchannel)
2721  {
2722  _console << " ======================================= " << std::endl;
2723  _console << " ======== Subchannel Print Outs ======== " << std::endl;
2724  _console << " ======================================= " << std::endl;
2725  _console << "Total flow area :" << Total_surface_area << " m^2" << std::endl;
2726  _console << "Assembly hydraulic diameter :" << bulk_Dh << " m" << std::endl;
2727  _console << "Assembly Re number :" << bulk_Re << " [-]" << std::endl;
2728  _console << "Bulk coolant temperature at outlet :" << T_bulk_out << " K" << std::endl;
2729  _console << "Power added to coolant is : " << power_out - power_in << " Watt" << std::endl;
2730  _console << "Mass flow rate in is : " << mass_flow_in << " kg/sec" << std::endl;
2731  _console << "Mass balance is : " << mass_flow_out - mass_flow_in << " kg/sec" << std::endl;
2732  _console << "User defined outlet pressure is : " << _P_out << " Pa" << std::endl;
2733  _console << " ======================================= " << std::endl;
2734  }
2735 
2736  if (MooseUtils::absoluteFuzzyLessEqual((power_out - power_in), -1.0))
2737  mooseWarning(
2738  "Energy conservation equation might not be solved correctly, Power added to coolant: " +
2739  std::to_string(power_out - power_in) + " Watt ");
2740 }
2741 
2742 void
2744 {
2745 }
Mat _amc_sys_mdot_mat
Axial momentum system matrix.
static const std::string PRESSURE_DROP
pressure drop
Definition: SubChannelApp.h:37
static const std::string MASS_FLOW_RATE
mass flow rate
Definition: SubChannelApp.h:29
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.
virtual const unsigned int & getNumOfChannels() const =0
Return the number of channels per layer.
std::unique_ptr< SolutionHandle > _Tduct_soln
std::unique_ptr< SolutionHandle > _duct_heat_flux_soln
static const std::string DENSITY
density
Definition: SubChannelApp.h:47
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.
static InputParameters validParams()
std::unique_ptr< SolutionHandle > _T_soln
std::unique_ptr< SolutionHandle > _h_soln
const PetscReal & _dtol
The divergence tolerance for the ksp linear solver.
const bool _monolithic_thermal_bool
Thermal monolithic bool.
virtual const std::vector< unsigned int > & getPinChannels(unsigned int i_pin) const =0
Return a vector of channel indices for a given Pin index.
virtual Node * getChannelNodeFromDuct(Node *channel_node)=0
Function that gets the channel node from the duct node.
T & getMesh(MooseMesh &mesh)
function to cast mesh
Definition: SCM.h:35
Mat _cmc_friction_force_mat
Cross momentum conservation - friction force.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
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.
virtual void zero() override final
virtual const std::vector< Real > & getZGrid() const
Get axial location of layers.
void computeMdot(int iblock)
Computes mass flow per channel for block iblock.
virtual const Real & getPitch() const
Return the pitch between 2 subchannels.
virtual EChannelType getSubchannelType(unsigned int index) const =0
Return the type of the subchannel for given subchannel index.
std::unique_ptr< SolutionHandle > _S_flow_soln
const bool _compute_density
Flag that activates or deactivates the calculation of density.
Mat _amc_friction_force_mat
Axial momentum conservation - friction force.
static const std::string PIN_DIAMETER
pin diameter
Definition: SubChannelApp.h:45
PetscScalar _added_K
Added resistances for monolithic convergence.
PetscErrorCode createPetscVector(Vec &v, PetscInt n)
Petsc Functions.
PetscFunctionBegin
static InputParameters validParams()
const bool _duct_mesh_exist
Flag that informs if there is a pin mesh or not.
static InputParameters validParams()
const Real & _T_tol
Convergence tolerance for the temperature loop in external solve.
void computeMu(int iblock)
Computes Viscosity per channel for block iblock.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
virtual Node * getPinNode(unsigned int i_pin, unsigned iz) const =0
Get the pin mesh node for a given pin index and elevation index.
PetscScalar computeInterpolatedValue(PetscScalar topValue, PetscScalar botValue, PetscScalar Peclet=0.0)
bool isRestarting() const
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...
virtual const std::vector< unsigned int > & getChannelGaps(unsigned int i_chan) const =0
Return a vector of gap indices for a given channel index.
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.
void addRequiredParam(const std::string &name, const std::string &doc_string)
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 override
Real computeAddedHeatDuct(unsigned int i_ch, unsigned int iz)
Function that computes the heat added by the duct, for channel i_ch and cell iz.
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)
static const std::string DUCT_HEAT_FLUX
duct linear heat rate
Definition: SubChannelApp.h:55
static const std::string DUCT_TEMPERATURE
duct temperature
Definition: SubChannelApp.h:57
SubChannel1PhaseProblem * schp
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.
static const std::string WETTED_PERIMETER
wetted perimeter
Definition: SubChannelApp.h:51
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.
static const std::string cp
Definition: NS.h:121
Mat _amc_advective_derivative_mat
Axial momentum conservation - advective (Eulerian) derivative.
static const std::string VISCOSITY
viscosity
Definition: SubChannelApp.h:49
Vec _hc_added_heat_rhs
Enthalpy conservation - source and sink.
processor_id_type n_processors() const
std::vector< Real > _z_grid
axial location of nodes
const std::string & name() const
static const std::string PIN_TEMPERATURE
pin temperature
Definition: SubChannelApp.h:43
const int & _T_maxit
Maximum iterations for the inner temperature loop.
Mat _amc_pressure_force_mat
Axial momentum conservation - pressure force.
static const std::string LINEAR_HEAT_RATE
linear heat rate
Definition: SubChannelApp.h:53
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
void computeWijResidual(int iblock)
Computes Residual Matrix based on the lateral momentum conservation equation for block iblock...
const std::vector< double > x
static const std::string S
Definition: NS.h:163
Real f(Real x)
Test function for Brents method.
std::unique_ptr< SolutionHandle > _q_prime_soln
virtual void syncSolutions(Direction direction) override
static const std::string mu
Definition: NS.h:123
static const std::string pitch
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.
static const std::string ENTHALPY
enthalpy
Definition: SubChannelApp.h:39
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
std::shared_ptr< AuxiliarySystem > _aux
virtual const unsigned int & getNumOfPins() const =0
Return the number of pins.
Mat _hc_cross_derivative_mat
Enthalpy conservation - cross flux derivative.
Mat _mc_axial_convection_mat
Mass conservation - axial convection.
virtual Node * getDuctNodeFromChannel(Node *channel_node)=0
Function that gets the duct node from the channel node.
void initialSetup() override
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.
virtual const std::vector< std::vector< Real > > & getKGrid() const
Get axial cell location and value of loss coefficient.
const MooseEnum _interpolation_scheme
The interpolation method used in constructing the systems.
std::unique_ptr< SolutionHandle > _Dpin_soln
Real volume(const MeshBase &mesh, unsigned int dim=libMesh::invalid_uint)
virtual Real getSubChannelPeripheralDuctWidth(unsigned int i_ch)=0
Function that computes the width of the duct cell that the peripheral subchannel i_ch sees...
Base class for the 1-phase steady-state/transient subchannel solver.
libMesh::DenseMatrix< Real > _Wij_residual_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)
static const std::string SUM_CROSSFLOW
sum of diversion crossflow
Definition: SubChannelApp.h:33
virtual const unsigned int & getNumOfAxialCells() const
Return the number of axial cells.
virtual Real getGapWidth(unsigned int axial_index, unsigned int gap_index) const =0
Return gap width for a given gap index.
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
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
PetscErrorCode formFunction(SNES, Vec x, Vec f, void *ctx)
virtual const std::vector< Node * > getDuctNodes() const =0
Function that return the vector with the maps to the nodes if they exist.
std::unique_ptr< SolutionHandle > _DP_soln
const bool _compute_viscosity
Flag that activates or deactivates the calculation of viscosity.
static const std::string PRESSURE
pressure
Definition: SubChannelApp.h:35
MooseApp & _app
libMesh::DenseMatrix< Real > & _Wij
virtual void externalSolve() override
void computeWijPrime(int iblock)
Computes turbulent crossflow per gap for block iblock.
static const std::string alpha
Definition: NS.h:134
virtual const Real & getCrossflowSign(unsigned int i_chan, unsigned int i_local) const =0
Return a sign for the crossflow given a subchannel index and local neighbor index.
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...
static const std::string SURFACE_AREA
surface area
Definition: SubChannelApp.h:31
bool absoluteFuzzyLessEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
void mooseWarning(Args &&... args) const
void resize(const unsigned int new_m, const unsigned int new_n)
const bool _verbose_subchannel
Boolean to printout information related to subchannel solve.
Mat _amc_cross_derivative_mat
Axial momentum conservation - cross flux derivative.
void mooseError(Args &&... args) const
const Real & _P_tol
Convergence tolerance for the pressure loop in external solve.
virtual Real computeBeta(unsigned int i_gap, unsigned int iz, bool enthalpy)=0
Computes turbulent mixing coefficient.
void addClassDescription(const std::string &doc_string)
virtual bool solverSystemConverged(const unsigned int) override
virtual const unsigned int & getNumOfGapsPerLayer() const =0
Return the number of gaps per layer.
struct SubChannel1PhaseProblem::FrictionStruct _friction_args
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
Mat _hc_advective_derivative_mat
Enthalpy conservation - advective (Eulerian) derivative;.
Base class for subchannel meshes.
void * ctx
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.
const ConsoleStream _console
static const std::string DISPLACEMENT
subchannel displacement
Definition: SubChannelApp.h:59
Mat _hc_time_derivative_mat
Enthalpy Enthalpy conservation - time derivative.
virtual bool isTransient() const override
PetscErrorCode implicitPetscSolve(int iblock)
Computes implicit solve using PetSc.
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)
bool isRecovering() const
static const std::string TEMPERATURE
temperature
Definition: SubChannelApp.h:41
Mat _cmc_pressure_force_mat
Cross momentum conservation - pressure force.
virtual Real & dt() const
MooseUnits pow(const MooseUnits &, int)
std::unique_ptr< SolutionHandle > _P_soln
static const std::string k
Definition: NS.h:130
std::unique_ptr< SolutionHandle > _w_perim_soln
void ErrorVector unsigned int
Definition: SCM.h:16
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
virtual const std::pair< unsigned int, unsigned int > & getGapChannels(unsigned int i_gap) const =0
Return a pair of subchannel indices for a given gap index.
std::unique_ptr< SolutionHandle > _Tpin_soln