Loading...
Searching...
No Matches
fvMatrix.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::fvMatrix
29
30Description
31 A special matrix type and solver, designed for finite volume
32 solutions of scalar equations.
33 Face addressing is used to make all matrix assembly
34 and solution loops vectorise.
35
36SourceFiles
37 fvMatrix.C
38 fvMatrixSolve.C
39 fvScalarMatrix.C
40
41\*---------------------------------------------------------------------------*/
42
43#ifndef Foam_fvMatrix_H
44#define Foam_fvMatrix_H
45
46#include "volFields.H"
47#include "surfaceFields.H"
48#include "lduMatrix.H"
49#include "tmp.H"
50#include "autoPtr.H"
51#include "dimensionedTypes.H"
52#include "className.H"
54#include "lduMesh.H"
55
56#include "fvMatrixExpression.H"
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60namespace Foam
61{
62
63// Forward Declarations
64template<class Type> class fvMatrix;
65template<class T> class UIndirectList;
66
67template<class Type>
69
70
71template<class Type>
73(
76);
77
78template<class Type>
80(
83);
84
85template<class Type>
87(
90);
91
92template<class Type>
94(
97);
98
99template<class Type>
101(
104);
105
106template<class Type>
108(
109 const tmp<fvMatrix<Type>>&,
111);
112
113
114/*---------------------------------------------------------------------------*\
115 Class fvMatrix Declaration
116\*---------------------------------------------------------------------------*/
117
118template<class Type>
119class fvMatrix
120:
121 public refCount,
122 public lduMatrix
123{
124public:
125
126 // Public Types
127
128 //- The geometric field type for psi
129 typedef
132
133 //- Field type for face flux (for non-orthogonal correction)
134 typedef
137
140
141
142private:
143
144 // Private Data
145
146 //- Const reference to field
147 // Converted into a non-const reference at the point of solution.
148 const psiFieldType& psi_;
149
150 //- Originating fvMatrices when assembling matrices. Empty if not used.
151 PtrList<fvMatrix<Type>> subMatrices_;
152
153 //- Is the fvMatrix using implicit formulation
154 bool useImplicit_;
155
156 //- Name of the lduAssembly
157 word lduAssemblyName_;
158
159 //- Number of fvMatrices added to this
160 label nMatrix_;
161
162 //- Dimension set
163 dimensionSet dimensions_;
164
165 //- Source term
166 Field<Type> source_;
167
168 //- Boundary scalar field containing pseudo-matrix coeffs
169 //- for internal cells
170 FieldField<Field, Type> internalCoeffs_;
171
172 //- Boundary scalar field containing pseudo-matrix coeffs
173 //- for boundary cells
174 FieldField<Field, Type> boundaryCoeffs_;
175
176 //- Face flux field for non-orthogonal correction
177 mutable std::unique_ptr<faceFluxFieldType> faceFluxCorrectionPtr_;
178
179
180protected:
181
182 //- Declare friendship with the fvSolver class
183 friend class fvSolver;
184
185
186 // Protected Member Functions
187
188 //- Add patch contribution to internal field
189 template<class Type2>
191 (
192 const labelUList& addr,
193 const Field<Type2>& pf,
194 Field<Type2>& intf
195 ) const;
196
197 template<class Type2>
199 (
200 const labelUList& addr,
201 const tmp<Field<Type2>>& tpf,
202 Field<Type2>& intf
203 ) const;
204
205 //- Subtract patch contribution from internal field
206 template<class Type2>
208 (
209 const labelUList& addr,
210 const Field<Type2>& pf,
211 Field<Type2>& intf
212 ) const;
213
214 template<class Type2>
216 (
217 const labelUList& addr,
218 const tmp<Field<Type2>>& tpf,
219 Field<Type2>& intf
220 ) const;
221
223 // Implicit helper functions
224
225 //- Name the implicit assembly addressing
226 // \return true if assembly is used
227 bool checkImplicit(const label fieldi = 0);
228
229
230 // Matrix completion functionality
231
232 void addBoundaryDiag
235 const direction cmpt
236 ) const;
237
239
243 const bool couples=true
244 ) const;
245
246
247 // Matrix manipulation functionality
248
249 //- Set solution in given cells to the specified values
250 template<template<class> class ListType>
252 (
253 const labelUList& cellLabels,
254 const ListType<Type>& values
255 );
257
258public:
259
260 //- Solver class returned by the solver function
261 //- used for systems in which it is useful to cache the solver for reuse.
262 // E.g. if the solver is potentially expensive to construct (AMG) and can
263 // be used several times (PISO)
264 class fvSolver
265 {
266 fvMatrix<Type>& fvMat_;
270 public:
271
272 // Constructors
273
275 :
276 fvMat_(fvMat),
277 solver_(std::move(sol))
278 {}
279
280
281 // Member Functions
283 //- Solve returning the solution statistics.
284 // Solver controls read from dictionary
285 SolverPerformance<Type> solve(const dictionary& solverControls);
286
287 //- Solve returning the solution statistics.
288 // Solver controls read from fvSolution
290 };
291
292
293 // Runtime information
294 ClassName("fvMatrix");
295
296
297 // Constructors
299 //- Construct given a field to solve for
301 (
303 const dimensionSet& ds
304 );
305
306 //- Copy construct
307 fvMatrix(const fvMatrix<Type>&);
309 //- Copy/move construct from tmp<fvMatrix<Type>>
310 fvMatrix(const tmp<fvMatrix<Type>>&);
311
312 //- Deprecated(2022-05) - construct with dimensionSet instead
313 // \deprecated(2022-05) - construct with dimensionSet instead
314 FOAM_DEPRECATED_FOR(2022-05, "Construct with dimensionSet")
316 (
317 const GeometricField<Type, fvPatchField, volMesh>& psi,
318 Istream& is
319 )
320 :
321 fvMatrix<Type>(psi, dimensionSet(is))
322 {}
323
324 //- Construct and return a clone
326 {
327 return tmp<fvMatrix<Type>>::New(*this);
328 }
329
330
331 //- Destructor
332 virtual ~fvMatrix();
333
335 // Member Functions
336
337 // Access
338
339 // Coupling
340
341 label nMatrices() const
343 return (nMatrix_ == 0 ? 1 : nMatrix_);
344 }
345
346 const fvMatrix<Type>& matrix(const label i) const
347 {
348 return (nMatrix_ == 0 ? *this : subMatrices_[i]);
349 }
350
351 fvMatrix<Type>& matrix(const label i)
352 {
353 return (nMatrix_ == 0 ? *this : subMatrices_[i]);
354 }
355
357 (
358 const label fieldi,
359 const label patchi
360 ) const
361 {
362 if (!lduMeshPtr())
363 {
364 return patchi;
365 }
366 else
367 {
368 return lduMeshPtr()->patchMap()[fieldi][patchi];
369 }
370 }
371
372 //- Transfer lower, upper, diag and source to this fvMatrix
374
375 //- Create or update ldu assembly
377
378 //- Access to lduPrimitiveMeshAssembly
380
381 //- Const Access to lduPrimitiveMeshAssembly
383
384 //- Manipulate matrix
386
387 //- Manipulate boundary/internal coeffs for coupling
389
390 //- Set interfaces
391 void setInterfaces
392 (
395 );
396
397 //- Add internal and boundary contribution to local patches
400 label fieldi,
401 const FieldField<Field, Type>& fluxContrib,
403 bool internal
404 ) const;
405
406 //- Return optional lduAdressing
408 {
409 return *lduMeshPtr();
410 }
411
412 //- Return psi
414 (
415 const label i = 0
416 ) const
417 {
418 return
419 (
420 (i == 0 && nMatrix_ == 0) ? psi_ : matrix(i).psi()
421 );
422 }
423
424
425 GeometricField<Type, fvPatchField, volMesh>& psi
426 (
427 const label i = 0
429 {
430 return
431 (
432 (i == 0 && nMatrix_ == 0)
433 ? const_cast
434 <
436 >(psi_)
437 : const_cast
440 >
441 (
442 matrix(i).psi()
444 );
445 }
446
447 //- Clear multiple fvMatrices
448 void clear()
449 {
450 subMatrices_.clear();
451 nMatrix_ = 0;
452 }
454
455 const dimensionSet& dimensions() const noexcept
456 {
457 return dimensions_;
459
461 {
462 return source_;
463 }
464
465 const Field<Type>& source() const noexcept
466 {
467 return source_;
468 }
469
470 //- fvBoundary scalar field containing pseudo-matrix coeffs
471 //- for internal cells
473 {
474 return internalCoeffs_;
475 }
476
477 //- fvBoundary scalar field containing pseudo-matrix coeffs
478 //- for internal cells
480 {
481 return internalCoeffs_;
482 }
483
484 //- fvBoundary scalar field containing pseudo-matrix coeffs
485 //- for boundary cells
487 {
488 return boundaryCoeffs_;
489 }
490
491 //- fvBoundary scalar field containing pseudo-matrix coeffs
492 //- for boundary cells
494 {
495 return boundaryCoeffs_;
496 }
497
498 //- Declare return type of the faceFluxCorrectionPtr() function
499 typedef std::unique_ptr<faceFluxFieldType> faceFluxFieldPtrType;
500
501 //- Return pointer to face-flux non-orthogonal correction field
503 {
504 return faceFluxCorrectionPtr_;
505 }
506
507 //- Set pointer to face-flux non-orthogonal correction field
509 {
510 faceFluxCorrectionPtr_.reset(flux);
511 }
512
513 //- True if face-flux non-orthogonal correction field exists
515 {
516 return bool(faceFluxCorrectionPtr_);
517 }
518
519
520 // Operations
521
522 //- Set solution in given cells to the specified value
523 //- and eliminate the corresponding equations from the matrix.
524 void setValues
525 (
526 const labelUList& cellLabels,
527 const Type& value
528 );
529
530 //- Set solution in given cells to the specified values
531 //- and eliminate the corresponding equations from the matrix.
532 void setValues
533 (
534 const labelUList& cellLabels,
535 const UList<Type>& values
536 );
537
538 //- Set solution in given cells to the specified values
539 //- and eliminate the corresponding equations from the matrix.
541 (
542 const labelUList& cellLabels,
543 const UIndirectList<Type>& values
544 );
545
546 //- Set reference level for solution
547 void setReference
548 (
549 const label celli,
550 const Type& value,
551 const bool forceReference = false
552 );
553
554 //- Set reference level for solution
555 void setReferences
556 (
557 const labelUList& cellLabels,
558 const Type& value,
559 const bool forceReference = false
560 );
561
562 //- Set reference level for solution
563 void setReferences
564 (
565 const labelUList& cellLabels,
566 const UList<Type>& values,
567 const bool forceReference = false
568 );
569
570 //- Set reference level for a component of the solution
571 //- on a given patch face
573 (
574 const label patchi,
575 const label facei,
576 const direction cmpt,
577 const scalar value
578 );
579
580 //- Add fvMatrix
582
583 //- Relax matrix (for steady-state solution).
584 // alpha = 1 : diagonally equal
585 // alpha < 1 : diagonally dominant
586 // alpha = 0 : do nothing
587 // Note: Requires positive diagonal.
588 void relax(const scalar alpha);
590 //- Relax matrix (for steady-state solution).
591 // alpha is read from controlDict
592 void relax();
593
594 //- Manipulate based on a boundary field
596 (
598 Boundary& values
599 );
600
601 //- Construct and return the solver
602 // Use the given solver controls
604
605 //- Construct and return the solver
606 // Solver controls read from fvSolution
608
609 //- Solve segregated or coupled returning the solution statistics.
610 // Use the given solver controls
612
613 //- Solve segregated or coupled returning the solution statistics.
614 // Solver controls read from fvSolution
616
617 //- Solve segregated returning the solution statistics.
618 // Use the given solver controls
620
621 //- Solve coupled returning the solution statistics.
622 // Use the given solver controls
624
625 //- Solve returning the solution statistics.
626 // Use the given solver controls
628
629 //- Solve returning the solution statistics.
630 // Uses \p name solver controls from fvSolution
632
633 //- Solve returning the solution statistics.
634 // Solver controls read from fvSolution
636
637 //- Return the matrix residual
638 tmp<Field<Type>> residual() const;
639
640 //- Return the matrix scalar diagonal
641 tmp<scalarField> D() const;
642
643 //- Return the matrix Type diagonal
644 tmp<Field<Type>> DD() const;
645
646 //- Return the central coefficient
647 tmp<volScalarField> A() const;
648
649 //- Return the H operation source
651
652 //- Return H(1)
653 tmp<volScalarField> H1() const;
654
655 //- Return the face-flux field from the matrix
657 flux() const;
658
659 //- Return the solver dictionary (from fvSolution) for \p name
660 const dictionary& solverDict(const word& name) const;
661
662 //- Return the solver dictionary for psi,
663 //- taking into account finalIteration
664 const dictionary& solverDict() const;
665
667 // Member Operators
668
669 void operator=(const fvMatrix<Type>&);
670 void operator=(const tmp<fvMatrix<Type>>&);
671
672 //- Inplace negate
673 void negate();
674
675 void operator+=(const fvMatrix<Type>&);
676 void operator+=(const tmp<fvMatrix<Type>>&);
678 void operator-=(const fvMatrix<Type>&);
679 void operator-=(const tmp<fvMatrix<Type>>&);
680
683 void operator+=
684 (
686 );
687
690 void operator-=
691 (
693 );
694
695 void operator+=(const dimensioned<Type>&);
696 void operator-=(const dimensioned<Type>&);
697
699 void operator-=(Foam::zero) {}
700
702 void operator*=(const tmp<volScalarField::Internal>&);
703 void operator*=(const tmp<volScalarField>&);
704
706
707
708 // Friend Operators
709
711 operator& <Type>
712 (
713 const fvMatrix<Type>&,
715 );
716
718 operator& <Type>
719 (
720 const fvMatrix<Type>&,
722 );
723
725 operator& <Type>
726 (
727 const tmp<fvMatrix<Type>>&,
729 );
730
732 operator& <Type>
733 (
734 const tmp<fvMatrix<Type>>&,
736 );
737
738
739 // Ostream Operator
740
741 friend Ostream& operator<< <Type>
743 Ostream&,
744 const fvMatrix<Type>&
745 );
746
747
748 // Expression templates
750 //- Construct given a field to solve for and expressions for
751 //- all the components
752 template<typename E>
754 (
757 <
758 E,
759 typename E::DiagExpr,
760 typename E::UpperExpr,
761 typename E::LowerExpr,
762 typename E::FaceFluxExpr,
763 typename E::SourceExpr
764 >& expr
765 );
766
767 //- Wrap value as expression
769 <
771 > expr() const;
772
773 //- Assign values from expression
774 template<typename E>
775 void operator=
776 (
778 <
779 E,
780 typename E::DiagExpr,
781 typename E::UpperExpr,
782 typename E::LowerExpr,
783 typename E::FaceFluxExpr,
784 typename E::SourceExpr
785 >& expr
786 );
788
789
790// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
791
792template<class Type>
793void checkMethod
794(
795 const fvMatrix<Type>&,
796 const fvMatrix<Type>&,
797 const char*
798);
799
800template<class Type>
801void checkMethod
803 const fvMatrix<Type>&,
805 const char*
806);
808template<class Type>
809void checkMethod
810(
811 const fvMatrix<Type>&,
812 const dimensioned<Type>&,
813 const char*
814);
815
816
817//- Solve returning the solution statistics given convergence tolerance
818// Use the given solver controls
819template<class Type>
821(
823 const dictionary& solverControls
825
826//- Solve returning the solution statistics given convergence tolerance,
827//- deleting temporary matrix after solution.
828// Use the given solver controls
829template<class Type>
831(
832 const tmp<fvMatrix<Type>>&,
833 const dictionary& solverControls
834);
836
837//- Solve returning the solution statistics given convergence tolerance
838// Uses \p name solver controls from fvSolution
839template<class Type>
842//- Solve returning the solution statistics given convergence tolerance,
843//- deleting temporary matrix after solution.
844// Uses \p name solver controls from fvSolution
845template<class Type>
847
848
849//- Solve returning the solution statistics given convergence tolerance
850// Uses solver controls from fvSolution
851template<class Type>
853
854//- Solve returning the solution statistics given convergence tolerance,
855//- deleting temporary matrix after solution.
856// Uses solver controls from fvSolution
857template<class Type>
859
861//- Return the correction form of the given matrix
862//- by subtracting the matrix multiplied by the current field
863template<class Type>
866
867//- Return the correction form of the given temporary matrix
868//- by subtracting the matrix multiplied by the current field
869template<class Type>
871
872
873// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
874
875template<class Type>
876tmp<fvMatrix<Type>> operator==
877(
878 const fvMatrix<Type>&,
879 const fvMatrix<Type>&
881
882template<class Type>
883tmp<fvMatrix<Type>> operator==
884(
885 const tmp<fvMatrix<Type>>&,
886 const fvMatrix<Type>&
888
889template<class Type>
890tmp<fvMatrix<Type>> operator==
891(
892 const fvMatrix<Type>&,
893 const tmp<fvMatrix<Type>>&
895
896template<class Type>
897tmp<fvMatrix<Type>> operator==
898(
899 const tmp<fvMatrix<Type>>&,
900 const tmp<fvMatrix<Type>>&
901);
902
904template<class Type>
905tmp<fvMatrix<Type>> operator==
906(
907 const fvMatrix<Type>&,
909);
910
911template<class Type>
912tmp<fvMatrix<Type>> operator==
913(
914 const fvMatrix<Type>&,
916);
918template<class Type>
919tmp<fvMatrix<Type>> operator==
920(
921 const fvMatrix<Type>&,
923);
924
925template<class Type>
926tmp<fvMatrix<Type>> operator==
927(
928 const tmp<fvMatrix<Type>>&,
930);
931
932template<class Type>
933tmp<fvMatrix<Type>> operator==
934(
935 const tmp<fvMatrix<Type>>&,
938
939template<class Type>
940tmp<fvMatrix<Type>> operator==
941(
942 const tmp<fvMatrix<Type>>&,
944);
945
946template<class Type>
947tmp<fvMatrix<Type>> operator==
948(
949 const fvMatrix<Type>&,
950 const dimensioned<Type>&
951);
952
953template<class Type>
954tmp<fvMatrix<Type>> operator==
955(
956 const tmp<fvMatrix<Type>>&,
957 const dimensioned<Type>&
958);
959
960
961template<class Type>
962tmp<fvMatrix<Type>> operator==
963(
964 const fvMatrix<Type>&,
965 const Foam::zero
966);
967
968template<class Type>
969tmp<fvMatrix<Type>> operator==
970(
971 const tmp<fvMatrix<Type>>&,
972 const Foam::zero
973);
974
975
976//- Unary negation
977template<class Type>
978tmp<fvMatrix<Type>> operator-
979(
980 const fvMatrix<Type>&
981);
982
983//- Unary negation
984template<class Type>
985tmp<fvMatrix<Type>> operator-
986(
987 const tmp<fvMatrix<Type>>&
988);
989
990
991template<class Type>
992tmp<fvMatrix<Type>> operator+
993(
994 const fvMatrix<Type>&,
995 const fvMatrix<Type>&
996);
997
998template<class Type>
999tmp<fvMatrix<Type>> operator+
1000(
1001 const tmp<fvMatrix<Type>>&,
1002 const fvMatrix<Type>&
1003);
1005template<class Type>
1006tmp<fvMatrix<Type>> operator+
1007(
1008 const fvMatrix<Type>&,
1009 const tmp<fvMatrix<Type>>&
1010);
1011
1012template<class Type>
1013tmp<fvMatrix<Type>> operator+
1014(
1015 const tmp<fvMatrix<Type>>&,
1016 const tmp<fvMatrix<Type>>&
1018
1019
1020template<class Type>
1021tmp<fvMatrix<Type>> operator+
1022(
1023 const fvMatrix<Type>&,
1025);
1027template<class Type>
1028tmp<fvMatrix<Type>> operator+
1029(
1030 const fvMatrix<Type>&,
1032);
1033
1034template<class Type>
1036(
1037 const fvMatrix<Type>&,
1039);
1040
1041template<class Type>
1042tmp<fvMatrix<Type>> operator+
1043(
1046);
1047
1048template<class Type>
1049tmp<fvMatrix<Type>> operator+
1050(
1051 const tmp<fvMatrix<Type>>&,
1053);
1054
1055template<class Type>
1056tmp<fvMatrix<Type>> operator+
1057(
1058 const tmp<fvMatrix<Type>>&,
1061
1062template<class Type>
1063tmp<fvMatrix<Type>> operator+
1064(
1067);
1068
1069template<class Type>
1070tmp<fvMatrix<Type>> operator+
1071(
1074);
1075
1076template<class Type>
1077tmp<fvMatrix<Type>> operator+
1078(
1081);
1082
1083template<class Type>
1084tmp<fvMatrix<Type>> operator+
1085(
1088);
1089
1090template<class Type>
1091tmp<fvMatrix<Type>> operator+
1092(
1094 const tmp<fvMatrix<Type>>&
1096
1097template<class Type>
1098tmp<fvMatrix<Type>> operator+
1099(
1101 const tmp<fvMatrix<Type>>&
1103
1104
1105template<class Type>
1106tmp<fvMatrix<Type>> operator+
1107(
1108 const fvMatrix<Type>&,
1110);
1111
1112template<class Type>
1113tmp<fvMatrix<Type>> operator+
1114(
1115 const tmp<fvMatrix<Type>>&,
1117);
1118
1119template<class Type>
1120tmp<fvMatrix<Type>> operator+
1121(
1122 const dimensioned<Type>&,
1124);
1125
1126template<class Type>
1127tmp<fvMatrix<Type>> operator+
1128(
1129 const dimensioned<Type>&,
1131);
1132
1133
1134template<class Type>
1135tmp<fvMatrix<Type>> operator-
1136(
1138 const fvMatrix<Type>&
1139);
1140
1141template<class Type>
1142tmp<fvMatrix<Type>> operator-
1143(
1145 const fvMatrix<Type>&
1146);
1147
1148template<class Type>
1149tmp<fvMatrix<Type>> operator-
1150(
1151 const fvMatrix<Type>&,
1153);
1154
1155template<class Type>
1156tmp<fvMatrix<Type>> operator-
1157(
1158 const tmp<fvMatrix<Type>>&,
1160);
1161
1162
1163template<class Type>
1164tmp<fvMatrix<Type>> operator-
1165(
1166 const fvMatrix<Type>&,
1168);
1169
1170template<class Type>
1171tmp<fvMatrix<Type>> operator-
1172(
1173 const fvMatrix<Type>&,
1175);
1176
1177template<class Type>
1178tmp<fvMatrix<Type>> operator-
1180 const fvMatrix<Type>&,
1182);
1183
1184template<class Type>
1185tmp<fvMatrix<Type>> operator-
1187 const tmp<fvMatrix<Type>>&,
1189);
1190
1191template<class Type>
1192tmp<fvMatrix<Type>> operator-
1194 const tmp<fvMatrix<Type>>&,
1196);
1197
1198template<class Type>
1199tmp<fvMatrix<Type>> operator-
1201 const tmp<fvMatrix<Type>>&,
1203);
1204
1205template<class Type>
1206tmp<fvMatrix<Type>> operator-
1209 const fvMatrix<Type>&
1210);
1211
1212template<class Type>
1213tmp<fvMatrix<Type>> operator-
1214(
1216 const fvMatrix<Type>&
1217);
1218
1219template<class Type>
1220tmp<fvMatrix<Type>> operator-
1221(
1223 const fvMatrix<Type>&
1224);
1225
1226template<class Type>
1227tmp<fvMatrix<Type>> operator-
1228(
1230 const tmp<fvMatrix<Type>>&
1231);
1232
1233template<class Type>
1234tmp<fvMatrix<Type>> operator-
1235(
1237 const tmp<fvMatrix<Type>>&
1238);
1239
1240template<class Type>
1241tmp<fvMatrix<Type>> operator-
1242(
1244 const tmp<fvMatrix<Type>>&
1245);
1246
1247
1248template<class Type>
1249tmp<fvMatrix<Type>> operator-
1251 const fvMatrix<Type>&,
1252 const dimensioned<Type>&
1253);
1254
1255template<class Type>
1256tmp<fvMatrix<Type>> operator-
1258 const tmp<fvMatrix<Type>>&,
1259 const dimensioned<Type>&
1260);
1261
1262template<class Type>
1263tmp<fvMatrix<Type>> operator-
1265 const dimensioned<Type>&,
1266 const fvMatrix<Type>&
1267);
1268
1269template<class Type>
1270tmp<fvMatrix<Type>> operator-
1272 const dimensioned<Type>&,
1273 const tmp<fvMatrix<Type>>&
1274);
1275
1276
1277template<class Type>
1279(
1281 const fvMatrix<Type>&
1282);
1283
1284template<class Type>
1286(
1288 const fvMatrix<Type>&
1289);
1290
1291template<class Type>
1293(
1294 const tmp<volScalarField>&,
1295 const fvMatrix<Type>&
1296);
1297
1298template<class Type>
1299tmp<fvMatrix<Type>> operator*
1302 const tmp<fvMatrix<Type>>&
1303);
1304
1305template<class Type>
1306tmp<fvMatrix<Type>> operator*
1309 const tmp<fvMatrix<Type>>&
1310);
1311
1312template<class Type>
1313tmp<fvMatrix<Type>> operator*
1315 const tmp<volScalarField>&,
1316 const tmp<fvMatrix<Type>>&
1317);
1318
1319
1320template<class Type>
1322(
1323 const dimensioned<scalar>&,
1324 const fvMatrix<Type>&
1325);
1326
1327template<class Type>
1328tmp<fvMatrix<Type>> operator*
1330 const dimensioned<scalar>&,
1331 const tmp<fvMatrix<Type>>&
1332);
1333
1334
1335// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1337} // End namespace Foam
1338
1339// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1340
1341#ifdef NoRepository
1342 #include "fvMatrix.C"
1343#endif
1344
1345// Specialisation for scalars
1346#include "fvScalarMatrix.H"
1347
1348// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1349
1350#endif
1351
1352// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Expression wrap of const reference to fvMatrix.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
DimensionedField< scalar, volMesh > Internal
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A dynamically resizable PtrList with allocation management.
Definition PtrDynList.H:58
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
A List with indirect addressing. Like IndirectList but does not store addressing.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Generic dimensioned Type class.
SolverPerformance< Type > solve()
Solve returning the solution statistics.
SolverPerformance< Type > solve(const dictionary &solverControls)
Solve returning the solution statistics.
fvSolver(fvMatrix< Type > &fvMat, autoPtr< lduMatrix::solver > &&sol)
Definition fvMatrix.H:308
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
void operator*=(const tmp< volScalarField > &)
Definition fvMatrix.C:1797
void operator+=(const tmp< DimensionedField< Type, volMesh > > &)
Definition fvMatrix.C:1681
void operator-=(const fvMatrix< Type > &)
Definition fvMatrix.C:1632
autoPtr< fvSolver > solver()
Construct and return the solver.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
void operator+=(const dimensioned< Type > &)
Definition fvMatrix.C:1736
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition fvMatrix.C:1282
tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > flux() const
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
tmp< volScalarField > A() const
Return the central coefficient.
Definition fvMatrix.C:1306
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
Definition fvMatrix.C:351
void manipulateMatrix(direction cmp)
Manipulate matrix.
Definition fvMatrix.C:791
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition fvMatrix.C:41
label globalPatchID(const label fieldi, const label patchi) const
Definition fvMatrix.H:410
void mapContributions(label fieldi, const FieldField< Field, Type > &fluxContrib, FieldField< Field, Type > &contrib, bool internal) const
Add internal and boundary contribution to local patches.
Definition fvMatrix.C:549
void setReferences(const labelUList &cellLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition fvMatrix.C:1020
void operator*=(const dimensioned< scalar > &)
Definition fvMatrix.C:1808
const Field< Type > & source() const noexcept
Definition fvMatrix.H:540
SolverPerformance< Type > solve(const word &name)
Solve returning the solution statistics.
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition fvMatrix.H:549
void transferFvMatrixCoeffs()
Transfer lower, upper, diag and source to this fvMatrix.
Definition fvMatrix.C:814
void relax()
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1242
SolverPerformance< Type > solve()
Solve returning the solution statistics.
void operator+=(const tmp< fvMatrix< Type > > &)
Definition fvMatrix.C:1624
void operator-=(const tmp< DimensionedField< Type, volMesh > > &)
Definition fvMatrix.C:1714
fvMatrix(const fvMatrix< Type > &)
Copy construct.
Definition fvMatrix.C:394
std::unique_ptr< faceFluxFieldType > faceFluxFieldPtrType
Definition fvMatrix.H:584
tmp< volScalarField > H1() const
Return H(1).
Definition fvMatrix.C:1377
void operator+=(const fvMatrix< Type > &)
Definition fvMatrix.C:1595
const fvMatrix< Type > & matrix(const label i) const
Definition fvMatrix.H:399
void createOrUpdateLduPrimitiveAssembly()
Create or update ldu assembly.
Definition fvMatrix.C:907
void operator+=(Foam::zero)
Definition fvMatrix.H:860
autoPtr< fvSolver > solver(const dictionary &)
Construct and return the solver.
faceFluxFieldPtrType & faceFluxCorrectionPtr()
Return pointer to face-flux non-orthogonal correction field.
Definition fvMatrix.H:589
FieldField< Field, Type > & internalCoeffs() noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition fvMatrix.H:558
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition fvMatrix.H:567
GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0)
Definition fvMatrix.H:499
virtual ~fvMatrix()
Destructor.
Definition fvMatrix.C:458
void operator*=(const tmp< volScalarField::Internal > &)
Definition fvMatrix.C:1786
GeometricField< scalar, fvPatchField, volMesh > psiFieldType
Definition fvMatrix.H:128
GeometricField< scalar, fvsPatchField, surfaceMesh > faceFluxFieldType
Definition fvMatrix.H:135
SolverPerformance< Type > solveSegregatedOrCoupled()
Solve segregated or coupled returning the solution statistics.
void setValues(const labelUList &cellLabels, const Type &value)
Set solution in given cells to the specified value and eliminate the corresponding equations from the...
Definition fvMatrix.C:971
label nMatrices() const
Definition fvMatrix.H:394
void operator+=(const DimensionedField< Type, volMesh > &)
Definition fvMatrix.C:1670
const dimensionSet & dimensions() const noexcept
Definition fvMatrix.H:530
const lduPrimitiveMeshAssembly & lduMeshAssembly()
Return optional lduAdressing.
Definition fvMatrix.H:478
void operator*=(const volScalarField::Internal &)
Definition fvMatrix.C:1756
const lduPrimitiveMeshAssembly * lduMeshPtr() const
Const Access to lduPrimitiveMeshAssembly.
Definition fvMatrix.C:894
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition fvMatrix.C:144
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition fvMatrix.C:169
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition fvMatrix.C:1004
void setValues(const labelUList &cellLabels, const UIndirectList< Type > &values)
Set solution in given cells to the specified values and eliminate the corresponding equations from th...
Definition fvMatrix.C:993
void operator=(const tmp< fvMatrix< Type > > &)
Definition fvMatrix.C:1572
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Set reference level for a component of the solution on a given patch face.
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &psi, const Expression::fvMatrixExpression< E, typename E::DiagExpr, typename E::UpperExpr, typename E::LowerExpr, typename E::FaceFluxExpr, typename E::SourceExpr > &expr)
Construct given a field to solve for and expressions for all the components.
Definition fvMatrix.C:2985
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition fvMatrix.C:1260
void setValues(const labelUList &cellLabels, const UList< Type > &values)
Set solution in given cells to the specified values and eliminate the corresponding equations from th...
Definition fvMatrix.C:982
void setInterfaces(lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces)
Set interfaces.
Definition fvMatrix.C:471
void operator-=(const dimensioned< Type > &)
Definition fvMatrix.C:1746
void operator-=(const DimensionedField< Type, volMesh > &)
Definition fvMatrix.C:1703
const dictionary & solverDict(const word &name) const
Return the solver dictionary (from fvSolution) for name.
Definition fvMatrix.C:1513
void operator-=(Foam::zero)
Definition fvMatrix.H:861
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
ClassName("fvMatrix")
void operator-=(const tmp< fvMatrix< Type > > &)
Definition fvMatrix.C:1661
SolverPerformance< Type > solveCoupled(const dictionary &)
Solve coupled returning the solution statistics.
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition fvMatrix.C:1325
tmp< Field< Type > > residual() const
Return the matrix residual.
void faceFluxCorrectionPtr(faceFluxFieldType *flux)
Set pointer to face-flux non-orthogonal correction field.
Definition fvMatrix.H:597
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition fvMatrix.C:116
void negate()
Inplace negate.
Definition fvMatrix.C:1580
tmp< fvMatrix< Type > > clone() const
Construct and return a clone.
Definition fvMatrix.H:376
void setReferences(const labelUList &cellLabels, const UList< Type > &values, const bool forceReference=false)
Set reference level for solution.
Definition fvMatrix.C:1043
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
void clear()
Clear multiple fvMatrices.
Definition fvMatrix.H:523
fvMatrix(const tmp< fvMatrix< Type > > &)
Copy/move construct from tmp<fvMatrix<Type>>.
Definition fvMatrix.C:420
void addFvMatrix(fvMatrix< Type > &matrix)
Add fvMatrix.
Definition fvMatrix.C:1065
tmp< scalarField > D() const
Return the matrix scalar diagonal.
Definition fvMatrix.C:1273
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition fvMatrix.C:79
void addToInternalField(const labelUList &addr, const tmp< Field< Type2 > > &tpf, Field< Type2 > &intf) const
Definition fvMatrix.C:65
const GeometricField< scalar, fvPatchField, volMesh > & psi(const label i=0) const
Definition fvMatrix.H:487
bool checkImplicit(const label fieldi=0)
Name the implicit assembly addressing.
Definition fvMatrix.C:314
Expression::fvMatrixConstRefWrap< fvMatrix< scalar > > expr() const
friend class fvSolver
Declare friendship with the fvSolver class.
Definition fvMatrix.H:205
lduPrimitiveMeshAssembly * lduMeshPtr()
Access to lduPrimitiveMeshAssembly.
Definition fvMatrix.C:881
void operator=(const fvMatrix< Type > &)
Definition fvMatrix.C:1534
void subtractFromInternalField(const labelUList &addr, const tmp< Field< Type2 > > &tpf, Field< Type2 > &intf) const
Definition fvMatrix.C:103
bool hasFaceFluxCorrection() const noexcept
True if face-flux non-orthogonal correction field exists.
Definition fvMatrix.H:605
void setValuesFromList(const labelUList &cellLabels, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition fvMatrix.C:219
FieldField< Field, Type > & boundaryCoeffs() noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition fvMatrix.H:576
Field< scalar > & source() noexcept
Definition fvMatrix.H:535
fvMatrix< Type > & matrix(const label i)
Definition fvMatrix.H:404
const dictionary & solverDict() const
Return the solver dictionary for psi, taking into account finalIteration.
Definition fvMatrix.C:1522
void setBounAndInterCoeffs()
Manipulate boundary/internal coeffs for coupling.
Definition fvMatrix.C:661
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
lduMatrix(const lduMesh &mesh)
Construct (without coefficients) for an LDU addressed mesh.
Definition lduMatrix.C:54
const scalarField & diag() const
Definition lduMatrix.C:195
An assembly of lduMatrix that is specific inter-region coupling through mapped patches.
const labelListList & patchMap() const
Return patchMap.
constexpr refCount() noexcept
Default construct, initializing count to 0.
Definition refCount.H:63
A class for managing temporary objects.
Definition tmp.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
Macro definitions for declaring ClassName(), NamespaceName(), etc.
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
const volScalarField & psi
Expression templates for fvMatrix.
A scalar instance of fvMatrix.
Namespace for OpenFOAM.
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition faMatrix.C:1026
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
uint8_t direction
Definition direction.H:49
const direction noexcept
Definition scalarImpl.H:265
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
volScalarField & alpha
CEqn solve()
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition stdFoam.H:43
Foam::surfaceFields.