-
Notifications
You must be signed in to change notification settings - Fork 1.2k
Expand file tree
/
Copy pathSubChannel1PhaseProblem.C
More file actions
2501 lines (2309 loc) · 104 KB
/
SubChannel1PhaseProblem.C
File metadata and controls
2501 lines (2309 loc) · 104 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html
#include "SubChannel1PhaseProblem.h"
#include "SystemBase.h"
#include "libmesh/petsc_vector.h"
#include "libmesh/dense_matrix.h"
#include "libmesh/dense_vector.h"
#include <iostream>
#include <cmath>
#include "AuxiliarySystem.h"
#include "SCM.h"
#include "SinglePhaseFluidProperties.h"
#include "SCMFrictionClosureBase.h"
#include "SCMHTCClosureBase.h"
struct Ctx
{
int iblock;
SubChannel1PhaseProblem * schp;
};
PetscErrorCode
formFunction(SNES, Vec x, Vec f, void * ctx)
{
const PetscScalar * xx;
PetscScalar * ff;
PetscInt size;
PetscFunctionBegin;
Ctx * cc = static_cast<Ctx *>(ctx);
LibmeshPetscCallQ(VecGetSize(x, &size));
libMesh::DenseVector<Real> solution_seed(size, 0.0);
LibmeshPetscCallQ(VecGetArrayRead(x, &xx));
for (PetscInt i = 0; i < size; i++)
solution_seed(i) = xx[i];
LibmeshPetscCallQ(VecRestoreArrayRead(x, &xx));
libMesh::DenseVector<Real> Wij_residual_vector =
cc->schp->residualFunction(cc->iblock, solution_seed);
LibmeshPetscCallQ(VecGetArray(f, &ff));
for (int i = 0; i < size; i++)
ff[i] = Wij_residual_vector(i);
LibmeshPetscCallQ(VecRestoreArray(f, &ff));
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
}
InputParameters
SubChannel1PhaseProblem::validParams()
{
// Enumerations
MooseEnum schemes("upwind downwind central_difference exponential", "central_difference");
MooseEnum gravity_direction("counter_flow co_flow none", "counter_flow");
// Inputs
InputParameters params = ExternalProblem::validParams();
params += PostprocessorInterface::validParams();
params.addClassDescription("Base class of the subchannel solvers");
params.addRequiredParam<unsigned int>("n_blocks", "The number of blocks in the axial direction");
params.addRequiredParam<Real>("CT", "Turbulent modeling parameter");
params.addParam<Real>("P_tol", 1e-6, "Pressure tolerance");
params.addParam<Real>("T_tol", 1e-6, "Temperature tolerance");
params.addParam<int>("T_maxit", 100, "Maximum number of iterations for inner temperature loop");
params.addParam<PetscReal>("rtol", 1e-6, "Relative tolerance for ksp solver");
params.addParam<PetscReal>("atol", 1e-6, "Absolute tolerance for ksp solver");
params.addParam<PetscReal>("dtol", 1e5, "Divergence tolerance or ksp solver");
params.addParam<PetscInt>("maxit", 1e4, "Maximum number of iterations for ksp solver");
params.addParam<MooseEnum>(
"interpolation_scheme",
schemes,
"Interpolation scheme used for the method. Default is central_difference");
params.addParam<MooseEnum>(
"gravity", gravity_direction, "Direction of gravity. Default is counter_flow");
params.addParam<bool>(
"implicit", false, "Boolean to define the use of explicit or implicit solution.");
params.addParam<bool>(
"segregated", true, "Boolean to define whether to use a segregated solution.");
params.addParam<bool>(
"verbose_subchannel", false, "Boolean to print out information related to subchannel solve.");
params.addParam<bool>(
"deformation",
false,
"Boolean that activates the deformation effect based on values for: displacement, Dpin");
params.addRequiredParam<bool>("compute_density", "Flag that enables the calculation of density");
params.addRequiredParam<bool>("compute_viscosity",
"Flag that enables the calculation of viscosity");
params.addRequiredParam<bool>(
"compute_power",
"Flag that informs whether we solve the Enthalpy/Temperature equations or not");
params.addRequiredParam<PostprocessorName>(
"P_out", "The postprocessor (or scalar) that provides the value of outlet pressure [Pa]");
params.addRequiredParam<UserObjectName>("fp", "Fluid properties user object name");
params.addRequiredParam<UserObjectName>("friction_closure",
"Closure computing the friction factor");
// Make these OPTIONAL here; enforce them conditionally
params.addParam<UserObjectName>(
"pin_HTC_closure", "Closure computing HTC on fuel pin (required if pin mesh exists).");
params.addParam<UserObjectName>("duct_HTC_closure",
"Closure computing HTC on duct (required if duct mesh exists).");
return params;
}
SubChannel1PhaseProblem::SubChannel1PhaseProblem(const InputParameters & params)
: ExternalProblem(params),
PostprocessorInterface(this),
_friction_args(/*i_ch=*/0, /*Re=*/1.0, /*S=*/0.0, /*w_perim=*/0.0),
_nusselt_args(
/*Re=*/1.0, /*Pr=*/1.0, std::numeric_limits<unsigned int>::max(), /*iz=*/0, /*i_ch=*/0),
_subchannel_mesh(SCM::getMesh<SubChannelMesh>(_mesh)),
_n_blocks(getParam<unsigned int>("n_blocks")),
_Wij(declareRestartableData<libMesh::DenseMatrix<Real>>("Wij")),
_g_grav(9.81),
_kij(_subchannel_mesh.getKij()),
_one(1.0),
_compute_density(getParam<bool>("compute_density")),
_compute_viscosity(getParam<bool>("compute_viscosity")),
_compute_power(getParam<bool>("compute_power")),
_pin_mesh_exist(_subchannel_mesh.pinMeshExist()),
_duct_mesh_exist(_subchannel_mesh.ductMeshExist()),
_P_out(getPostprocessorValue("P_out")),
_CT(getParam<Real>("CT")),
_P_tol(getParam<Real>("P_tol")),
_T_tol(getParam<Real>("T_tol")),
_T_maxit(getParam<int>("T_maxit")),
_rtol(getParam<PetscReal>("rtol")),
_atol(getParam<PetscReal>("atol")),
_dtol(getParam<PetscReal>("dtol")),
_maxit(getParam<PetscInt>("maxit")),
_interpolation_scheme(getParam<MooseEnum>("interpolation_scheme")),
_gravity_direction(getParam<MooseEnum>("gravity")),
_dir_grav(computeGravityDir(_gravity_direction)),
_implicit_bool(getParam<bool>("implicit")),
_segregated_bool(getParam<bool>("segregated")),
_verbose_subchannel(getParam<bool>("verbose_subchannel")),
_deformation(getParam<bool>("deformation")),
_fp(nullptr),
_Tpin_soln(nullptr),
_duct_heat_flux_soln(nullptr),
_Tduct_soln(nullptr)
{
if (_pin_mesh_exist && !isParamValid("pin_HTC_closure"))
paramError("pin_HTC_closure", "required when a pin mesh exists.");
if (_duct_mesh_exist && !isParamValid("duct_HTC_closure"))
paramError("duct_HTC_closure", "required when a duct mesh exists.");
// NOTE: The four quantities below are 0 for processor_id != 0
_n_cells = _subchannel_mesh.getNumOfAxialCells();
_n_gaps = _subchannel_mesh.getNumOfGapsPerLayer();
_n_pins = _subchannel_mesh.getNumOfPins();
_n_channels = _subchannel_mesh.getNumOfChannels();
// NOTE: The four quantities above are 0 for processor_id != 0
_z_grid = _subchannel_mesh.getZGrid();
_block_size = _n_cells / _n_blocks;
// Turbulent crossflow (stuff that live on the gaps)
if (!_app.isRestarting() && !_app.isRecovering())
{
_Wij.resize(_n_gaps, _n_cells + 1);
_Wij.zero();
}
_Wij_old.resize(_n_gaps, _n_cells + 1);
_Wij_old.zero();
_WijPrime.resize(_n_gaps, _n_cells + 1);
_WijPrime.zero();
_Wij_residual_matrix.resize(_n_gaps, _block_size);
_Wij_residual_matrix.zero();
_converged = true;
// Mass conservation components
LibmeshPetscCall(
createPetscMatrix(_mc_sumWij_mat, _block_size * _n_channels, _block_size * _n_gaps));
LibmeshPetscCall(createPetscVector(_Wij_vec, _block_size * _n_gaps));
LibmeshPetscCall(createPetscVector(_prod, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_prodp, _block_size * _n_channels));
LibmeshPetscCall(createPetscMatrix(
_mc_axial_convection_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_mc_axial_convection_rhs, _block_size * _n_channels));
// Axial momentum conservation components
LibmeshPetscCall(createPetscMatrix(
_amc_turbulent_cross_flows_mat, _block_size * _n_gaps, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_turbulent_cross_flows_rhs, _block_size * _n_gaps));
LibmeshPetscCall(createPetscMatrix(
_amc_time_derivative_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_time_derivative_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscMatrix(
_amc_advective_derivative_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_advective_derivative_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscMatrix(
_amc_cross_derivative_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_cross_derivative_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscMatrix(
_amc_friction_force_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_friction_force_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_gravity_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscMatrix(
_amc_pressure_force_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_pressure_force_rhs, _block_size * _n_channels));
LibmeshPetscCall(
createPetscMatrix(_amc_sys_mdot_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_amc_sys_mdot_rhs, _block_size * _n_channels));
// Lateral momentum conservation components
LibmeshPetscCall(
createPetscMatrix(_cmc_time_derivative_mat, _block_size * _n_gaps, _block_size * _n_gaps));
LibmeshPetscCall(createPetscVector(_cmc_time_derivative_rhs, _block_size * _n_gaps));
LibmeshPetscCall(createPetscMatrix(
_cmc_advective_derivative_mat, _block_size * _n_gaps, _block_size * _n_gaps));
LibmeshPetscCall(createPetscVector(_cmc_advective_derivative_rhs, _block_size * _n_gaps));
LibmeshPetscCall(
createPetscMatrix(_cmc_friction_force_mat, _block_size * _n_gaps, _block_size * _n_gaps));
LibmeshPetscCall(createPetscVector(_cmc_friction_force_rhs, _block_size * _n_gaps));
LibmeshPetscCall(
createPetscMatrix(_cmc_pressure_force_mat, _block_size * _n_gaps, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_cmc_pressure_force_rhs, _block_size * _n_gaps));
LibmeshPetscCall(
createPetscMatrix(_cmc_sys_Wij_mat, _block_size * _n_gaps, _block_size * _n_gaps));
LibmeshPetscCall(createPetscVector(_cmc_sys_Wij_rhs, _block_size * _n_gaps));
// Energy conservation components
LibmeshPetscCall(createPetscMatrix(
_hc_time_derivative_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_hc_time_derivative_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscMatrix(
_hc_advective_derivative_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_hc_advective_derivative_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscMatrix(
_hc_cross_derivative_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_hc_cross_derivative_rhs, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_hc_added_heat_rhs, _block_size * _n_channels));
LibmeshPetscCall(
createPetscMatrix(_hc_sys_h_mat, _block_size * _n_channels, _block_size * _n_channels));
LibmeshPetscCall(createPetscVector(_hc_sys_h_rhs, _block_size * _n_channels));
if ((_n_blocks == _n_cells) && _implicit_bool)
{
mooseError(name(),
": When implicit number of blocks can't be equal to number of cells. This will "
"cause problems with the subchannel interpolation scheme.");
}
}
void
SubChannel1PhaseProblem::initialSetup()
{
ExternalProblem::initialSetup();
_fp = &getUserObject<SinglePhaseFluidProperties>(getParam<UserObjectName>("fp"));
_friction_closure =
&getUserObject<SCMFrictionClosureBase>(getParam<UserObjectName>("friction_closure"));
// Create variables for output and storage
_mdot_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::MASS_FLOW_RATE));
_SumWij_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::SUM_CROSSFLOW));
_P_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PRESSURE));
_DP_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PRESSURE_DROP));
_h_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::ENTHALPY));
_T_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::TEMPERATURE));
if (_pin_mesh_exist)
{
_Tpin_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PIN_TEMPERATURE));
_Dpin_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::PIN_DIAMETER));
_pin_HTC_closure =
&getUserObject<SCMHTCClosureBase>(getParam<UserObjectName>("pin_HTC_closure"));
}
_rho_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DENSITY));
_mu_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::VISCOSITY));
_S_flow_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::SURFACE_AREA));
_w_perim_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::WETTED_PERIMETER));
_q_prime_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::LINEAR_HEAT_RATE));
_displacement_soln =
std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DISPLACEMENT));
_ff_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::FRICTION_FACTOR));
if (_duct_mesh_exist)
{
_duct_heat_flux_soln =
std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DUCT_HEAT_FLUX));
_Tduct_soln = std::make_unique<SolutionHandle>(getVariable(0, SubChannelApp::DUCT_TEMPERATURE));
_duct_HTC_closure =
&getUserObject<SCMHTCClosureBase>(getParam<UserObjectName>("duct_HTC_closure"));
}
}
SubChannel1PhaseProblem::~SubChannel1PhaseProblem()
{
PetscErrorCode ierr = cleanUp();
if (ierr)
mooseError(name(), ": Error in memory cleanup");
}
PetscErrorCode
SubChannel1PhaseProblem::cleanUp()
{
PetscFunctionBegin;
// We need to clean up the petsc matrices/vectors
// Mass conservation components
LibmeshPetscCall(MatDestroy(&_mc_sumWij_mat));
LibmeshPetscCall(VecDestroy(&_Wij_vec));
LibmeshPetscCall(VecDestroy(&_prod));
LibmeshPetscCall(VecDestroy(&_prodp));
LibmeshPetscCall(MatDestroy(&_mc_axial_convection_mat));
LibmeshPetscCall(VecDestroy(&_mc_axial_convection_rhs));
// Axial momentum conservation components
LibmeshPetscCall(MatDestroy(&_amc_turbulent_cross_flows_mat));
LibmeshPetscCall(VecDestroy(&_amc_turbulent_cross_flows_rhs));
LibmeshPetscCall(MatDestroy(&_amc_time_derivative_mat));
LibmeshPetscCall(VecDestroy(&_amc_time_derivative_rhs));
LibmeshPetscCall(MatDestroy(&_amc_advective_derivative_mat));
LibmeshPetscCall(VecDestroy(&_amc_advective_derivative_rhs));
LibmeshPetscCall(MatDestroy(&_amc_cross_derivative_mat));
LibmeshPetscCall(VecDestroy(&_amc_cross_derivative_rhs));
LibmeshPetscCall(MatDestroy(&_amc_friction_force_mat));
LibmeshPetscCall(VecDestroy(&_amc_friction_force_rhs));
LibmeshPetscCall(VecDestroy(&_amc_gravity_rhs));
LibmeshPetscCall(MatDestroy(&_amc_pressure_force_mat));
LibmeshPetscCall(VecDestroy(&_amc_pressure_force_rhs));
LibmeshPetscCall(MatDestroy(&_amc_sys_mdot_mat));
LibmeshPetscCall(VecDestroy(&_amc_sys_mdot_rhs));
// Lateral momentum conservation components
LibmeshPetscCall(MatDestroy(&_cmc_time_derivative_mat));
LibmeshPetscCall(VecDestroy(&_cmc_time_derivative_rhs));
LibmeshPetscCall(MatDestroy(&_cmc_advective_derivative_mat));
LibmeshPetscCall(VecDestroy(&_cmc_advective_derivative_rhs));
LibmeshPetscCall(MatDestroy(&_cmc_friction_force_mat));
LibmeshPetscCall(VecDestroy(&_cmc_friction_force_rhs));
LibmeshPetscCall(MatDestroy(&_cmc_pressure_force_mat));
LibmeshPetscCall(VecDestroy(&_cmc_pressure_force_rhs));
LibmeshPetscCall(MatDestroy(&_cmc_sys_Wij_mat));
LibmeshPetscCall(VecDestroy(&_cmc_sys_Wij_rhs));
// Energy conservation components
LibmeshPetscCall(MatDestroy(&_hc_time_derivative_mat));
LibmeshPetscCall(VecDestroy(&_hc_time_derivative_rhs));
LibmeshPetscCall(MatDestroy(&_hc_advective_derivative_mat));
LibmeshPetscCall(VecDestroy(&_hc_advective_derivative_rhs));
LibmeshPetscCall(MatDestroy(&_hc_cross_derivative_mat));
LibmeshPetscCall(VecDestroy(&_hc_cross_derivative_rhs));
LibmeshPetscCall(VecDestroy(&_hc_added_heat_rhs));
LibmeshPetscCall(MatDestroy(&_hc_sys_h_mat));
LibmeshPetscCall(VecDestroy(&_hc_sys_h_rhs));
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
}
bool
SubChannel1PhaseProblem::solverSystemConverged(const unsigned int)
{
return _converged;
}
PetscScalar
SubChannel1PhaseProblem::computeInterpolationCoefficients(PetscScalar Peclet)
{
switch (_interpolation_scheme)
{
case 0: // upwind interpolation
return 1.0;
case 1: // downwind interpolation
return 0.0;
case 2: // central_difference interpolation
return 0.5;
case 3: // exponential interpolation (Peclet limited)
return ((Peclet - 1.0) * std::exp(Peclet) + 1) / (Peclet * (std::exp(Peclet) - 1.) + 1e-10);
default:
mooseError(name(),
": Interpolation scheme should be a string: upwind, downwind, central_difference, "
"exponential");
}
}
PetscScalar
SubChannel1PhaseProblem::computeInterpolatedValue(PetscScalar topValue,
PetscScalar botValue,
PetscScalar Peclet)
{
PetscScalar alpha = computeInterpolationCoefficients(Peclet);
return alpha * botValue + (1.0 - alpha) * topValue;
}
void
SubChannel1PhaseProblem::computeWijFromSolve(int iblock)
{
const unsigned int last_node = (iblock + 1) * _block_size;
const unsigned int first_node = iblock * _block_size + 1;
// Initial guess, port crossflow of block (iblock) into a vector that will act as my initial guess
libMesh::DenseVector<Real> solution_seed(_n_gaps * _block_size, 0.0);
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
{
int i = _n_gaps * (iz - first_node) + i_gap; // column wise transfer
solution_seed(i) = _Wij(i_gap, iz);
}
}
// Solving the combined lateral momentum equation for Wij using a PETSc solver and update vector
// root
libMesh::DenseVector<Real> root(_n_gaps * _block_size, 0.0);
LibmeshPetscCall(petscSnesSolver(iblock, solution_seed, root));
// Assign the solution to the cross-flow matrix
int i = 0;
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
for (unsigned int i_gap = 0; i_gap < _n_gaps; i_gap++)
{
_Wij(i_gap, iz) = root(i);
i++;
}
}
}
void
SubChannel1PhaseProblem::computeSumWij(int iblock)
{
const unsigned int last_node = (iblock + 1) * _block_size;
const unsigned int first_node = iblock * _block_size + 1;
// Add to solution vector if explicit
if (!_implicit_bool)
{
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
{
auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
Real sumWij = 0.0;
// Calculate sum of crossflow into channel i from channels j around i
unsigned int counter = 0;
for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
{
sumWij += _subchannel_mesh.getCrossflowSign(i_ch, counter) * _Wij(i_gap, iz);
counter++;
}
// The net crossflow coming out of cell i [kg/sec]
_SumWij_soln->set(node_out, sumWij);
}
}
}
// Add to matrix if implicit
else
{
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
unsigned int iz_ind = iz - first_node;
for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
{
// Calculate sum of crossflow into channel i from channels j around i
unsigned int counter = 0;
for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
{
PetscInt row = i_ch + _n_channels * iz_ind;
PetscInt col = i_gap + _n_gaps * iz_ind;
PetscScalar value = _subchannel_mesh.getCrossflowSign(i_ch, counter);
LibmeshPetscCall(MatSetValues(_mc_sumWij_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
counter++;
}
}
}
LibmeshPetscCall(MatAssemblyBegin(_mc_sumWij_mat, MAT_FINAL_ASSEMBLY));
LibmeshPetscCall(MatAssemblyEnd(_mc_sumWij_mat, MAT_FINAL_ASSEMBLY));
if (_segregated_bool)
{
Vec loc_prod;
Vec loc_Wij;
LibmeshPetscCall(VecDuplicate(_amc_sys_mdot_rhs, &loc_prod));
LibmeshPetscCall(VecDuplicate(_Wij_vec, &loc_Wij));
LibmeshPetscCall(populateVectorFromDense<libMesh::DenseMatrix<Real>>(
loc_Wij, _Wij, first_node, last_node, _n_gaps));
LibmeshPetscCall(MatMult(_mc_sumWij_mat, loc_Wij, loc_prod));
LibmeshPetscCall(populateSolutionChan<SolutionHandle>(
loc_prod, *_SumWij_soln, first_node, last_node, _n_channels));
LibmeshPetscCall(VecDestroy(&loc_prod));
LibmeshPetscCall(VecDestroy(&loc_Wij));
}
}
}
void
SubChannel1PhaseProblem::computeMdot(int iblock)
{
const unsigned int last_node = (iblock + 1) * _block_size;
const unsigned int first_node = iblock * _block_size + 1;
if (!_implicit_bool)
{
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
auto dz = _z_grid[iz] - _z_grid[iz - 1];
for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
{
auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
auto volume = dz * (*_S_flow_soln)(node_in);
auto time_term = _TR * ((*_rho_soln)(node_out)-_rho_soln->old(node_out)) * volume / _dt;
// Wij positive out of i into j;
auto mdot_out = (*_mdot_soln)(node_in) - (*_SumWij_soln)(node_out)-time_term;
if (mdot_out < 0)
{
_console << "Wij = : " << _Wij << "\n";
mooseError(name(),
" : Calculation of negative mass flow mdot_out = : ",
mdot_out,
" Axial Level= : ",
iz,
" - Implicit solves are required for recirculating flow.");
}
_mdot_soln->set(node_out, mdot_out); // kg/sec
}
}
}
else
{
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
auto dz = _z_grid[iz] - _z_grid[iz - 1];
auto iz_ind = iz - first_node;
for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
{
auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
auto volume = dz * (*_S_flow_soln)(node_in);
// Adding time derivative to the RHS
auto time_term = _TR * ((*_rho_soln)(node_out)-_rho_soln->old(node_out)) * volume / _dt;
PetscInt row_vec = i_ch + _n_channels * iz_ind;
PetscScalar value_vec = -1.0 * time_term;
LibmeshPetscCall(
VecSetValues(_mc_axial_convection_rhs, 1, &row_vec, &value_vec, INSERT_VALUES));
// Imposing bottom boundary condition or adding of diagonal elements
if (iz == first_node)
{
PetscScalar value_vec = (*_mdot_soln)(node_in);
PetscInt row_vec = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(
VecSetValues(_mc_axial_convection_rhs, 1, &row_vec, &value_vec, ADD_VALUES));
}
else
{
PetscInt row = i_ch + _n_channels * iz_ind;
PetscInt col = i_ch + _n_channels * (iz_ind - 1);
PetscScalar value = -1.0;
LibmeshPetscCall(
MatSetValues(_mc_axial_convection_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
}
// Adding diagonal elements
PetscInt row = i_ch + _n_channels * iz_ind;
PetscInt col = i_ch + _n_channels * iz_ind;
PetscScalar value = 1.0;
LibmeshPetscCall(
MatSetValues(_mc_axial_convection_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
// Adding cross flows RHS
if (_segregated_bool)
{
PetscScalar value_vec_2 = -1.0 * (*_SumWij_soln)(node_out);
PetscInt row_vec_2 = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(
VecSetValues(_mc_axial_convection_rhs, 1, &row_vec_2, &value_vec_2, ADD_VALUES));
}
}
}
LibmeshPetscCall(MatAssemblyBegin(_mc_axial_convection_mat, MAT_FINAL_ASSEMBLY));
LibmeshPetscCall(MatAssemblyEnd(_mc_axial_convection_mat, MAT_FINAL_ASSEMBLY));
if (_segregated_bool)
{
KSP ksploc;
PC pc;
Vec sol;
LibmeshPetscCall(VecDuplicate(_mc_axial_convection_rhs, &sol));
LibmeshPetscCall(KSPCreate(PETSC_COMM_SELF, &ksploc));
LibmeshPetscCall(KSPSetOperators(ksploc, _mc_axial_convection_mat, _mc_axial_convection_mat));
LibmeshPetscCall(KSPGetPC(ksploc, &pc));
LibmeshPetscCall(PCSetType(pc, PCJACOBI));
LibmeshPetscCall(KSPSetTolerances(ksploc, _rtol, _atol, _dtol, _maxit));
LibmeshPetscCall(KSPSetFromOptions(ksploc));
LibmeshPetscCall(KSPSolve(ksploc, _mc_axial_convection_rhs, sol));
LibmeshPetscCall(populateSolutionChan<SolutionHandle>(
sol, *_mdot_soln, first_node, last_node, _n_channels));
LibmeshPetscCall(VecZeroEntries(_mc_axial_convection_rhs));
LibmeshPetscCall(KSPDestroy(&ksploc));
LibmeshPetscCall(VecDestroy(&sol));
}
}
}
void
SubChannel1PhaseProblem::computeDP(int iblock)
{
const unsigned int last_node = (iblock + 1) * _block_size;
const unsigned int first_node = iblock * _block_size + 1;
if (!_implicit_bool)
{
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
auto k_grid = _subchannel_mesh.getKGrid();
auto dz = _z_grid[iz] - _z_grid[iz - 1];
for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
{
auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
auto rho_in = (*_rho_soln)(node_in);
auto rho_out = (*_rho_soln)(node_out);
auto mu_in = (*_mu_soln)(node_in);
auto S = (*_S_flow_soln)(node_in);
auto w_perim = (*_w_perim_soln)(node_in);
// hydraulic diameter in the i direction
auto Dh_i = 4.0 * S / w_perim;
auto time_term = _TR * ((*_mdot_soln)(node_out)-_mdot_soln->old(node_out)) * dz / _dt -
dz * 2.0 * (*_mdot_soln)(node_out) * (rho_out - _rho_soln->old(node_out)) /
rho_in / _dt;
auto mass_term1 =
Utility::pow<2>((*_mdot_soln)(node_out)) * (1.0 / S / rho_out - 1.0 / S / rho_in);
auto mass_term2 = -2.0 * (*_mdot_soln)(node_out) * (*_SumWij_soln)(node_out) / S / rho_in;
auto crossflow_term = 0.0;
auto turbulent_term = 0.0;
unsigned int counter = 0;
for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
{
auto chans = _subchannel_mesh.getGapChannels(i_gap);
unsigned int ii_ch = chans.first;
unsigned int jj_ch = chans.second;
auto * node_in_i = _subchannel_mesh.getChannelNode(ii_ch, iz - 1);
auto * node_in_j = _subchannel_mesh.getChannelNode(jj_ch, iz - 1);
auto * node_out_i = _subchannel_mesh.getChannelNode(ii_ch, iz);
auto * node_out_j = _subchannel_mesh.getChannelNode(jj_ch, iz);
auto rho_i = (*_rho_soln)(node_in_i);
auto rho_j = (*_rho_soln)(node_in_j);
auto Si = (*_S_flow_soln)(node_in_i);
auto Sj = (*_S_flow_soln)(node_in_j);
Real u_star = 0.0;
// figure out donor axial velocity
if (_Wij(i_gap, iz) > 0.0)
u_star = (*_mdot_soln)(node_out_i) / Si / rho_i;
else
u_star = (*_mdot_soln)(node_out_j) / Sj / rho_j;
crossflow_term +=
_subchannel_mesh.getCrossflowSign(i_ch, counter) * _Wij(i_gap, iz) * u_star;
turbulent_term += _WijPrime(i_gap, iz) * (2 * (*_mdot_soln)(node_out) / rho_in / S -
(*_mdot_soln)(node_out_j) / Sj / rho_j -
(*_mdot_soln)(node_out_i) / Si / rho_i);
counter++;
}
turbulent_term *= _CT;
auto Re = (((*_mdot_soln)(node_in) / S) * Dh_i / mu_in);
_friction_args = FrictionStruct(i_ch, Re, S, w_perim);
Real ff = _friction_closure->computeFrictionFactor(_friction_args);
_ff_soln->set(node_out, ff);
/// Upwind local form loss
auto ki = 0.0;
if ((*_mdot_soln)(node_out) >= 0)
ki = k_grid[i_ch][iz - 1];
else
ki = k_grid[i_ch][iz];
auto friction_term = (ff * dz / Dh_i + ki) * 0.5 *
(*_mdot_soln)(node_out)*std::abs((*_mdot_soln)(node_out)) /
(S * (*_rho_soln)(node_out));
auto gravity_term = _dir_grav * _g_grav * (*_rho_soln)(node_out)*dz * S;
auto DP = (1 / S) * (time_term + mass_term1 + mass_term2 + crossflow_term + turbulent_term +
friction_term + gravity_term); // Pa
_DP_soln->set(node_out, DP);
}
}
}
else
{
LibmeshPetscCall(MatZeroEntries(_amc_time_derivative_mat));
LibmeshPetscCall(MatZeroEntries(_amc_advective_derivative_mat));
LibmeshPetscCall(MatZeroEntries(_amc_cross_derivative_mat));
LibmeshPetscCall(MatZeroEntries(_amc_friction_force_mat));
LibmeshPetscCall(VecZeroEntries(_amc_time_derivative_rhs));
LibmeshPetscCall(VecZeroEntries(_amc_advective_derivative_rhs));
LibmeshPetscCall(VecZeroEntries(_amc_cross_derivative_rhs));
LibmeshPetscCall(VecZeroEntries(_amc_friction_force_rhs));
LibmeshPetscCall(VecZeroEntries(_amc_gravity_rhs));
LibmeshPetscCall(MatZeroEntries(_amc_sys_mdot_mat));
LibmeshPetscCall(VecZeroEntries(_amc_sys_mdot_rhs));
for (unsigned int iz = first_node; iz < last_node + 1; iz++)
{
auto k_grid = _subchannel_mesh.getKGrid();
auto dz = _z_grid[iz] - _z_grid[iz - 1];
auto iz_ind = iz - first_node;
for (unsigned int i_ch = 0; i_ch < _n_channels; i_ch++)
{
// inlet and outlet nodes
auto * node_in = _subchannel_mesh.getChannelNode(i_ch, iz - 1);
auto * node_out = _subchannel_mesh.getChannelNode(i_ch, iz);
// interpolation weight coefficient
PetscScalar Pe = 0.5;
if (_interpolation_scheme == 3)
{
// Compute the Peclet number
auto S_in = (*_S_flow_soln)(node_in);
auto S_out = (*_S_flow_soln)(node_out);
auto S_interp = computeInterpolatedValue(S_out, S_in, 0.5);
auto w_perim_in = (*_w_perim_soln)(node_in);
auto w_perim_out = (*_w_perim_soln)(node_out);
auto w_perim_interp = this->computeInterpolatedValue(w_perim_out, w_perim_in, 0.5);
auto mdot_loc =
this->computeInterpolatedValue((*_mdot_soln)(node_out), (*_mdot_soln)(node_in), 0.5);
auto mu_in = (*_mu_soln)(node_in);
auto mu_out = (*_mu_soln)(node_out);
auto mu_interp = this->computeInterpolatedValue(mu_out, mu_in, 0.5);
auto Dh_i = 4.0 * S_interp / w_perim_interp;
// Compute friction factor
auto Re = ((mdot_loc / S_interp) * Dh_i / mu_interp);
_friction_args = FrictionStruct(i_ch, Re, S_interp, w_perim_interp);
Real ff = _friction_closure->computeFrictionFactor(_friction_args);
_ff_soln->set(node_out, ff);
/// Upwind local form loss
auto ki = 0.0;
if ((*_mdot_soln)(node_out) >= 0)
ki = k_grid[i_ch][iz - 1];
else
ki = k_grid[i_ch][iz];
Pe = 1.0 / ((ff * dz / Dh_i + ki) * 0.5) * mdot_loc / std::abs(mdot_loc);
}
auto alpha = computeInterpolationCoefficients(Pe);
// inlet, outlet, and interpolated density
auto rho_in = (*_rho_soln)(node_in);
auto rho_out = (*_rho_soln)(node_out);
auto rho_interp = computeInterpolatedValue(rho_out, rho_in, Pe);
// inlet, outlet, and interpolated viscosity
auto mu_in = (*_mu_soln)(node_in);
auto mu_out = (*_mu_soln)(node_out);
auto mu_interp = computeInterpolatedValue(mu_out, mu_in, Pe);
// inlet, outlet, and interpolated axial surface area
auto S_in = (*_S_flow_soln)(node_in);
auto S_out = (*_S_flow_soln)(node_out);
auto S_interp = computeInterpolatedValue(S_out, S_in, Pe);
// inlet, outlet, and interpolated wetted perimeter
auto w_perim_in = (*_w_perim_soln)(node_in);
auto w_perim_out = (*_w_perim_soln)(node_out);
auto w_perim_interp = computeInterpolatedValue(w_perim_out, w_perim_in, Pe);
// hydraulic diameter in the i direction
auto Dh_i = 4.0 * S_interp / w_perim_interp;
/// Time derivative term
if (iz == first_node)
{
PetscScalar value_vec_tt = -1.0 * _TR * alpha * (*_mdot_soln)(node_in)*dz / _dt;
PetscInt row_vec_tt = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(
VecSetValues(_amc_time_derivative_rhs, 1, &row_vec_tt, &value_vec_tt, ADD_VALUES));
}
else
{
PetscInt row_tt = i_ch + _n_channels * iz_ind;
PetscInt col_tt = i_ch + _n_channels * (iz_ind - 1);
PetscScalar value_tt = _TR * alpha * dz / _dt;
LibmeshPetscCall(MatSetValues(
_amc_time_derivative_mat, 1, &row_tt, 1, &col_tt, &value_tt, INSERT_VALUES));
}
// Adding diagonal elements
PetscInt row_tt = i_ch + _n_channels * iz_ind;
PetscInt col_tt = i_ch + _n_channels * iz_ind;
PetscScalar value_tt = _TR * (1.0 - alpha) * dz / _dt;
LibmeshPetscCall(MatSetValues(
_amc_time_derivative_mat, 1, &row_tt, 1, &col_tt, &value_tt, INSERT_VALUES));
// Adding RHS elements
PetscScalar mdot_old_interp =
computeInterpolatedValue(_mdot_soln->old(node_out), _mdot_soln->old(node_in), Pe);
PetscScalar value_vec_tt = _TR * mdot_old_interp * dz / _dt;
PetscInt row_vec_tt = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(
VecSetValues(_amc_time_derivative_rhs, 1, &row_vec_tt, &value_vec_tt, ADD_VALUES));
/// Advective derivative term
if (iz == first_node)
{
PetscScalar value_vec_at = Utility::pow<2>((*_mdot_soln)(node_in)) / (S_in * rho_in);
PetscInt row_vec_at = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(VecSetValues(
_amc_advective_derivative_rhs, 1, &row_vec_at, &value_vec_at, ADD_VALUES));
}
else
{
PetscInt row_at = i_ch + _n_channels * iz_ind;
PetscInt col_at = i_ch + _n_channels * (iz_ind - 1);
PetscScalar value_at = -1.0 * std::abs((*_mdot_soln)(node_in)) / (S_in * rho_in);
LibmeshPetscCall(MatSetValues(
_amc_advective_derivative_mat, 1, &row_at, 1, &col_at, &value_at, INSERT_VALUES));
}
// Adding diagonal elements
PetscInt row_at = i_ch + _n_channels * iz_ind;
PetscInt col_at = i_ch + _n_channels * iz_ind;
PetscScalar value_at = std::abs((*_mdot_soln)(node_out)) / (S_out * rho_out);
LibmeshPetscCall(MatSetValues(
_amc_advective_derivative_mat, 1, &row_at, 1, &col_at, &value_at, INSERT_VALUES));
/// Cross derivative term
unsigned int counter = 0;
unsigned int cross_index = iz; // iz-1;
for (auto i_gap : _subchannel_mesh.getChannelGaps(i_ch))
{
auto chans = _subchannel_mesh.getGapChannels(i_gap);
unsigned int ii_ch = chans.first;
unsigned int jj_ch = chans.second;
auto * node_in_i = _subchannel_mesh.getChannelNode(ii_ch, iz - 1);
auto * node_in_j = _subchannel_mesh.getChannelNode(jj_ch, iz - 1);
auto * node_out_i = _subchannel_mesh.getChannelNode(ii_ch, iz);
auto * node_out_j = _subchannel_mesh.getChannelNode(jj_ch, iz);
auto rho_i =
computeInterpolatedValue((*_rho_soln)(node_out_i), (*_rho_soln)(node_in_i), Pe);
auto rho_j =
computeInterpolatedValue((*_rho_soln)(node_out_j), (*_rho_soln)(node_in_j), Pe);
auto S_i =
computeInterpolatedValue((*_S_flow_soln)(node_out_i), (*_S_flow_soln)(node_in_i), Pe);
auto S_j =
computeInterpolatedValue((*_S_flow_soln)(node_out_j), (*_S_flow_soln)(node_in_j), Pe);
auto u_star = 0.0;
// figure out donor axial velocity
if (_Wij(i_gap, cross_index) > 0.0)
{
if (iz == first_node)
{
u_star = (*_mdot_soln)(node_in_i) / S_i / rho_i;
PetscScalar value_vec_ct = -1.0 * alpha *
_subchannel_mesh.getCrossflowSign(i_ch, counter) *
_Wij(i_gap, cross_index) * u_star;
PetscInt row_vec_ct = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(VecSetValues(
_amc_cross_derivative_rhs, 1, &row_vec_ct, &value_vec_ct, ADD_VALUES));
}
else
{
PetscScalar value_ct = alpha * _subchannel_mesh.getCrossflowSign(i_ch, counter) *
_Wij(i_gap, cross_index) / S_i / rho_i;
PetscInt row_ct = i_ch + _n_channels * iz_ind;
PetscInt col_ct = ii_ch + _n_channels * (iz_ind - 1);
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
}
PetscScalar value_ct = (1.0 - alpha) *
_subchannel_mesh.getCrossflowSign(i_ch, counter) *
_Wij(i_gap, cross_index) / S_i / rho_i;
PetscInt row_ct = i_ch + _n_channels * iz_ind;
PetscInt col_ct = ii_ch + _n_channels * iz_ind;
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
}
else if (_Wij(i_gap, cross_index) < 0.0) // _Wij=0 operations not necessary
{
if (iz == first_node)
{
u_star = (*_mdot_soln)(node_in_j) / S_j / rho_j;
PetscScalar value_vec_ct = -1.0 * alpha *
_subchannel_mesh.getCrossflowSign(i_ch, counter) *
_Wij(i_gap, cross_index) * u_star;
PetscInt row_vec_ct = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(VecSetValues(
_amc_cross_derivative_rhs, 1, &row_vec_ct, &value_vec_ct, ADD_VALUES));
}
else
{
PetscScalar value_ct = alpha * _subchannel_mesh.getCrossflowSign(i_ch, counter) *
_Wij(i_gap, cross_index) / S_j / rho_j;
PetscInt row_ct = i_ch + _n_channels * iz_ind;
PetscInt col_ct = jj_ch + _n_channels * (iz_ind - 1);
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
}
PetscScalar value_ct = (1.0 - alpha) *
_subchannel_mesh.getCrossflowSign(i_ch, counter) *
_Wij(i_gap, cross_index) / S_j / rho_j;
PetscInt row_ct = i_ch + _n_channels * iz_ind;
PetscInt col_ct = jj_ch + _n_channels * iz_ind;
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_ct, ADD_VALUES));
}
if (iz == first_node)
{
PetscScalar value_vec_ct = -2.0 * alpha * (*_mdot_soln)(node_in)*_CT *
_WijPrime(i_gap, cross_index) / (rho_interp * S_interp);
value_vec_ct += alpha * (*_mdot_soln)(node_in_j)*_CT * _WijPrime(i_gap, cross_index) /
(rho_j * S_j);
value_vec_ct += alpha * (*_mdot_soln)(node_in_i)*_CT * _WijPrime(i_gap, cross_index) /
(rho_i * S_i);
PetscInt row_vec_ct = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(
VecSetValues(_amc_cross_derivative_rhs, 1, &row_vec_ct, &value_vec_ct, ADD_VALUES));
}
else
{
PetscScalar value_center_ct =
2.0 * alpha * _CT * _WijPrime(i_gap, cross_index) / (rho_interp * S_interp);
PetscInt row_ct = i_ch + _n_channels * iz_ind;
PetscInt col_ct = i_ch + _n_channels * (iz_ind - 1);
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_center_ct, ADD_VALUES));
PetscScalar value_left_ct =
-1.0 * alpha * _CT * _WijPrime(i_gap, cross_index) / (rho_j * S_j);
row_ct = i_ch + _n_channels * iz_ind;
col_ct = jj_ch + _n_channels * (iz_ind - 1);
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_left_ct, ADD_VALUES));
PetscScalar value_right_ct =
-1.0 * alpha * _CT * _WijPrime(i_gap, cross_index) / (rho_i * S_i);
row_ct = i_ch + _n_channels * iz_ind;
col_ct = ii_ch + _n_channels * (iz_ind - 1);
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_right_ct, ADD_VALUES));
}
PetscScalar value_center_ct =
2.0 * (1.0 - alpha) * _CT * _WijPrime(i_gap, cross_index) / (rho_interp * S_interp);
PetscInt row_ct = i_ch + _n_channels * iz_ind;
PetscInt col_ct = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_center_ct, ADD_VALUES));
PetscScalar value_left_ct =
-1.0 * (1.0 - alpha) * _CT * _WijPrime(i_gap, cross_index) / (rho_j * S_j);
row_ct = i_ch + _n_channels * iz_ind;
col_ct = jj_ch + _n_channels * iz_ind;
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_left_ct, ADD_VALUES));
PetscScalar value_right_ct =
-1.0 * (1.0 - alpha) * _CT * _WijPrime(i_gap, cross_index) / (rho_i * S_i);
row_ct = i_ch + _n_channels * iz_ind;
col_ct = ii_ch + _n_channels * iz_ind;
LibmeshPetscCall(MatSetValues(
_amc_cross_derivative_mat, 1, &row_ct, 1, &col_ct, &value_right_ct, ADD_VALUES));
counter++;
}
/// Friction term
PetscScalar mdot_interp =
computeInterpolatedValue((*_mdot_soln)(node_out), (*_mdot_soln)(node_in), Pe);
auto Re = ((mdot_interp / S_interp) * Dh_i / mu_interp);
_friction_args = FrictionStruct(i_ch, Re, S_interp, w_perim_interp);
Real ff = _friction_closure->computeFrictionFactor(_friction_args);
_ff_soln->set(node_out, ff);
/// Upwind local form loss
auto ki = 0.0;
if ((*_mdot_soln)(node_out) >= 0)
ki = k_grid[i_ch][iz - 1];
else
ki = k_grid[i_ch][iz];
auto coef = (ff * dz / Dh_i + ki) * 0.5 * std::abs((*_mdot_soln)(node_out)) /
(S_interp * rho_interp);
if (iz == first_node)
{
PetscScalar value_vec = -1.0 * alpha * coef * (*_mdot_soln)(node_in);
PetscInt row_vec = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(
VecSetValues(_amc_friction_force_rhs, 1, &row_vec, &value_vec, ADD_VALUES));
}
else
{
PetscInt row = i_ch + _n_channels * iz_ind;
PetscInt col = i_ch + _n_channels * (iz_ind - 1);
PetscScalar value = alpha * coef;
LibmeshPetscCall(
MatSetValues(_amc_friction_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
}
// Adding diagonal elements
PetscInt row = i_ch + _n_channels * iz_ind;
PetscInt col = i_ch + _n_channels * iz_ind;
PetscScalar value = (1.0 - alpha) * coef;
LibmeshPetscCall(
MatSetValues(_amc_friction_force_mat, 1, &row, 1, &col, &value, INSERT_VALUES));
/// Gravity force
PetscScalar value_vec = _dir_grav * -1.0 * _g_grav * rho_interp * dz * S_interp;
PetscInt row_vec = i_ch + _n_channels * iz_ind;
LibmeshPetscCall(VecSetValues(_amc_gravity_rhs, 1, &row_vec, &value_vec, ADD_VALUES));
}
}
/// Assembling system
LibmeshPetscCall(MatZeroEntries(_amc_sys_mdot_mat));
LibmeshPetscCall(VecZeroEntries(_amc_sys_mdot_rhs));
LibmeshPetscCall(MatAssemblyBegin(_amc_time_derivative_mat, MAT_FINAL_ASSEMBLY));