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 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::fvMatrix
29 
30 Description
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 
36 SourceFiles
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 // Forward Declarations
62 template<class Type> class fvMatrix;
63 template<class T> class UIndirectList;
64 
65 template<class Type>
66 Ostream& operator<<(Ostream&, const fvMatrix<Type>&);
67 
68 
69 template<class Type>
70 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
71 (
72  const fvMatrix<Type>&,
73  const DimensionedField<Type, volMesh>&
74 );
75 
76 template<class Type>
77 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
78 (
79  const fvMatrix<Type>&,
80  const tmp<DimensionedField<Type, volMesh>>&
81 );
82 
83 template<class Type>
84 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
85 (
86  const fvMatrix<Type>&,
87  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
88 );
89 
90 template<class Type>
91 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
92 (
93  const tmp<fvMatrix<Type>>&,
94  const DimensionedField<Type, volMesh>&
95 );
96 
97 template<class Type>
98 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
99 (
100  const tmp<fvMatrix<Type>>&,
101  const tmp<DimensionedField<Type, volMesh>>&
102 );
103 
104 template<class Type>
105 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
106 (
107  const tmp<fvMatrix<Type>>&,
108  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
109 );
110 
111 
112 /*---------------------------------------------------------------------------*\
113  Class fvMatrix Declaration
114 \*---------------------------------------------------------------------------*/
115 
116 template<class Type>
117 class fvMatrix
118 :
119  public refCount,
120  public lduMatrix
121 {
122 public:
123 
124  // Public Types
125 
126  //- The geometric field type for psi
127  typedef
129  psiFieldType;
130 
131  //- Field type for face flux (for non-orthogonal correction)
132  typedef
135 
138 
139 
140 private:
141 
142  // Private Data
143 
144  //- Const reference to field
145  // Converted into a non-const reference at the point of solution.
146  const psiFieldType& psi_;
147 
148  //- Originating fvMatrices when assembling matrices. Empty if not used.
149  PtrList<fvMatrix<Type>> subMatrices_;
150 
151  //- Is the fvMatrix using implicit formulation
152  bool useImplicit_;
153 
154  //- Name of the lduAssembly
155  word lduAssemblyName_;
156 
157  //- Number of fvMatrices added to this
158  label nMatrix_;
159 
160  //- Dimension set
161  dimensionSet dimensions_;
162 
163  //- Source term
164  Field<Type> source_;
165 
166  //- Boundary scalar field containing pseudo-matrix coeffs
167  //- for internal cells
168  FieldField<Field, Type> internalCoeffs_;
169 
170  //- Boundary scalar field containing pseudo-matrix coeffs
171  //- for boundary cells
172  FieldField<Field, Type> boundaryCoeffs_;
173 
174  //- Face flux field for non-orthogonal correction
175  mutable faceFluxFieldType* faceFluxCorrectionPtr_;
176 
177 
178 protected:
179 
180  //- Declare friendship with the fvSolver class
181  friend class fvSolver;
182 
183 
184  // Protected Member Functions
185 
186  //- Add patch contribution to internal field
187  template<class Type2>
188  void addToInternalField
189  (
190  const labelUList& addr,
191  const Field<Type2>& pf,
192  Field<Type2>& intf
193  ) const;
194 
195  template<class Type2>
196  void addToInternalField
197  (
198  const labelUList& addr,
199  const tmp<Field<Type2>>& tpf,
200  Field<Type2>& intf
201  ) const;
202 
203  //- Subtract patch contribution from internal field
204  template<class Type2>
206  (
207  const labelUList& addr,
208  const Field<Type2>& pf,
209  Field<Type2>& intf
210  ) const;
211 
212  template<class Type2>
214  (
215  const labelUList& addr,
216  const tmp<Field<Type2>>& tpf,
217  Field<Type2>& intf
218  ) const;
219 
220 
221  // Implicit helper functions
222 
223  //- Name the implicit assembly addressing
224  // \return true if assembly is used
225  bool checkImplicit(const label fieldi = 0);
226 
227 
228  // Matrix completion functionality
229 
230  void addBoundaryDiag
231  (
232  scalarField& diag,
233  const direction cmpt
234  ) const;
235 
237 
238  void addBoundarySource
239  (
241  const bool couples=true
242  ) const;
243 
244 
245  // Matrix manipulation functionality
246 
247  //- Set solution in given cells to the specified values
248  template<template<class> class ListType>
249  void setValuesFromList
250  (
251  const labelUList& cellLabels,
252  const ListType<Type>& values
253  );
254 
255 
256 public:
257 
258  //- Solver class returned by the solver function
259  //- used for systems in which it is useful to cache the solver for reuse.
260  // E.g. if the solver is potentially expensive to construct (AMG) and can
261  // be used several times (PISO)
262  class fvSolver
263  {
264  fvMatrix<Type>& fvMat_;
265 
267 
268  public:
269 
270  // Constructors
271 
273  :
274  fvMat_(fvMat),
275  solver_(std::move(sol))
276  {}
277 
278 
279  // Member Functions
280 
281  //- Solve returning the solution statistics.
282  // Solver controls read from dictionary
283  SolverPerformance<Type> solve(const dictionary& solverControls);
284 
285  //- Solve returning the solution statistics.
286  // Solver controls read from fvSolution
288  };
289 
290 
291  // Runtime information
292  ClassName("fvMatrix");
293 
294 
295  // Constructors
297  //- Construct given a field to solve for
298  fvMatrix
299  (
301  const dimensionSet& ds
302  );
303 
304  //- Copy construct
305  fvMatrix(const fvMatrix<Type>&);
307  //- Copy/move construct from tmp<fvMatrix<Type>>
308  fvMatrix(const tmp<fvMatrix<Type>>&);
309 
310  //- Deprecated(2022-05) - construct with dimensionSet instead
311  // \deprecated(2022-05) - construct with dimensionSet instead
312  FOAM_DEPRECATED_FOR(2022-05, "Construct with dimensionSet")
313  fvMatrix
314  (
315  const GeometricField<Type, fvPatchField, volMesh>& psi,
316  Istream& is
317  )
318  :
319  fvMatrix<Type>(psi, dimensionSet(is))
320  {}
321 
322  //- Construct and return a clone
323  tmp<fvMatrix<Type>> clone() const
324  {
325  return tmp<fvMatrix<Type>>::New(*this);
326  }
327 
328 
329  //- Destructor
330  virtual ~fvMatrix();
331 
332 
333  // Member Functions
334 
335  // Access
336 
337  // Coupling
338 
339  label nMatrices() const
340  {
341  return (nMatrix_ == 0 ? 1 : nMatrix_);
342  }
343 
344  const fvMatrix<Type>& matrix(const label i) const
345  {
346  return (nMatrix_ == 0 ? *this : subMatrices_[i]);
347  }
348 
349  fvMatrix<Type>& matrix(const label i)
350  {
351  return (nMatrix_ == 0 ? *this : subMatrices_[i]);
352  }
353 
354  label globalPatchID
355  (
356  const label fieldi,
357  const label patchi
358  ) const
359  {
360  if (!lduMeshPtr())
361  {
362  return patchi;
363  }
364  else
365  {
366  return lduMeshPtr()->patchMap()[fieldi][patchi];
367  }
368  }
369 
370  //- Transfer lower, upper, diag and source to this fvMatrix
371  void transferFvMatrixCoeffs();
372 
373  //- Create or update ldu assembly
375 
376  //- Access to lduPrimitiveMeshAssembly
378 
379  //- Const Access to lduPrimitiveMeshAssembly
380  const lduPrimitiveMeshAssembly* lduMeshPtr() const;
381 
382  //- Manipulate matrix
383  void manipulateMatrix(direction cmp);
384 
385  //- Manipulate boundary/internal coeffs for coupling
386  void setBounAndInterCoeffs();
387 
388  //- Set interfaces
389  void setInterfaces
390  (
393  );
394 
395  //- Add internal and boundary contribution to local patches
396  void mapContributions
397  (
398  label fieldi,
399  const FieldField<Field, Type>& fluxContrib,
400  FieldField<Field, Type>& contrib,
401  bool internal
402  ) const;
403 
404  //- Return optional lduAdressing
406  {
407  return *lduMeshPtr();
408  }
409 
410  //- Return psi
412  (
413  const label i = 0
414  ) const
415  {
416  return
417  (
418  (i == 0 && nMatrix_ == 0) ? psi_ : matrix(i).psi()
419  );
420  }
421 
422 
423  GeometricField<Type, fvPatchField, volMesh>& psi
424  (
425  const label i = 0
426  )
427  {
428  return
429  (
430  (i == 0 && nMatrix_ == 0)
431  ? const_cast
432  <
433  GeometricField<Type, fvPatchField, volMesh>&
434  >(psi_)
435  : const_cast
436  <
437  GeometricField<Type, fvPatchField, volMesh>&
438  >
439  (
440  matrix(i).psi()
441  )
442  );
443  }
444 
445  //- Clear multiple fvMatrices
446  void clear()
447  {
448  subMatrices_.clear();
449  nMatrix_ = 0;
450  }
451 
452 
453  const dimensionSet& dimensions() const noexcept
454  {
455  return dimensions_;
456  }
457 
458  Field<Type>& source() noexcept
459  {
460  return source_;
461  }
462 
463  const Field<Type>& source() const noexcept
464  {
465  return source_;
466  }
467 
468  //- fvBoundary scalar field containing pseudo-matrix coeffs
469  //- for internal cells
470  const FieldField<Field, Type>& internalCoeffs() const noexcept
471  {
472  return internalCoeffs_;
473  }
474 
475  //- fvBoundary scalar field containing pseudo-matrix coeffs
476  //- for internal cells
478  {
479  return internalCoeffs_;
480  }
481 
482  //- fvBoundary scalar field containing pseudo-matrix coeffs
483  //- for boundary cells
485  {
486  return boundaryCoeffs_;
487  }
488 
489  //- fvBoundary scalar field containing pseudo-matrix coeffs
490  //- for boundary cells
492  {
493  return boundaryCoeffs_;
494  }
495 
496  //- Declare return type of the faceFluxCorrectionPtr() function
499 
500  //- Return pointer to face-flux non-orthogonal correction field
502  {
503  return faceFluxCorrectionPtr_;
504  }
505 
506  //- True if face-flux non-orthogonal correction field exists
507  bool hasFaceFluxCorrection() const noexcept
508  {
509  return bool(faceFluxCorrectionPtr_);
510  }
511 
512 
513  // Operations
514 
515  //- Set solution in given cells to the specified value
516  //- and eliminate the corresponding equations from the matrix.
517  void setValues
518  (
519  const labelUList& cellLabels,
520  const Type& value
521  );
522 
523  //- Set solution in given cells to the specified values
524  //- and eliminate the corresponding equations from the matrix.
525  void setValues
526  (
527  const labelUList& cellLabels,
529  );
530 
531  //- Set solution in given cells to the specified values
532  //- and eliminate the corresponding equations from the matrix.
534  (
535  const labelUList& cellLabels,
537  );
539  //- Set reference level for solution
540  void setReference
541  (
542  const label celli,
543  const Type& value,
544  const bool forceReference = false
545  );
546 
547  //- Set reference level for solution
548  void setReferences
549  (
550  const labelUList& cellLabels,
551  const Type& value,
552  const bool forceReference = false
553  );
554 
555  //- Set reference level for solution
557  (
558  const labelUList& cellLabels,
559  const UList<Type>& values,
560  const bool forceReference = false
561  );
562 
563  //- Set reference level for a component of the solution
564  //- on a given patch face
566  (
567  const label patchi,
568  const label facei,
569  const direction cmpt,
570  const scalar value
571  );
572 
573  //- Add fvMatrix
575 
576  //- Relax matrix (for steady-state solution).
577  // alpha = 1 : diagonally equal
578  // alpha < 1 : diagonally dominant
579  // alpha = 0 : do nothing
580  // Note: Requires positive diagonal.
581  void relax(const scalar alpha);
582 
583  //- Relax matrix (for steady-state solution).
584  // alpha is read from controlDict
585  void relax();
586 
587  //- Manipulate based on a boundary field
589  (
591  Boundary& values
592  );
593 
594  //- Construct and return the solver
595  // Use the given solver controls
597 
598  //- Construct and return the solver
599  // Solver controls read from fvSolution
601 
602  //- Solve segregated or coupled returning the solution statistics.
603  // Use the given solver controls
605 
606  //- Solve segregated or coupled returning the solution statistics.
607  // Solver controls read from fvSolution
609 
610  //- Solve segregated returning the solution statistics.
611  // Use the given solver controls
613 
614  //- Solve coupled returning the solution statistics.
615  // Use the given solver controls
617 
618  //- Solve returning the solution statistics.
619  // Use the given solver controls
621 
622  //- Solve returning the solution statistics.
623  // Uses \p name solver controls from fvSolution
625 
626  //- Solve returning the solution statistics.
627  // Solver controls read from fvSolution
629 
630  //- Return the matrix residual
631  tmp<Field<Type>> residual() const;
632 
633  //- Return the matrix scalar diagonal
634  tmp<scalarField> D() const;
635 
636  //- Return the matrix Type diagonal
637  tmp<Field<Type>> DD() const;
638 
639  //- Return the central coefficient
640  tmp<volScalarField> A() const;
641 
642  //- Return the H operation source
644 
645  //- Return H(1)
646  tmp<volScalarField> H1() const;
647 
648  //- Return the face-flux field from the matrix
650  flux() const;
651 
652  //- Return the solver dictionary (from fvSolution) for \p name
653  const dictionary& solverDict(const word& name) const;
654 
655  //- Return the solver dictionary for psi,
656  //- taking into account finalIteration
657  const dictionary& solverDict() const;
658 
659 
660  // Member Operators
661 
662  void operator=(const fvMatrix<Type>&);
663  void operator=(const tmp<fvMatrix<Type>>&);
664 
665  //- Inplace negate
666  void negate();
667 
668  void operator+=(const fvMatrix<Type>&);
669  void operator+=(const tmp<fvMatrix<Type>>&);
670 
671  void operator-=(const fvMatrix<Type>&);
672  void operator-=(const tmp<fvMatrix<Type>>&);
673 
676  void operator+=
677  (
679  );
680 
683  void operator-=
684  (
686  );
687 
688  void operator+=(const dimensioned<Type>&);
689  void operator-=(const dimensioned<Type>&);
690 
691  void operator+=(const Foam::zero) {}
692  void operator-=(const Foam::zero) {}
693 
695  void operator*=(const tmp<volScalarField::Internal>&);
696  void operator*=(const tmp<volScalarField>&);
697 
698  void operator*=(const dimensioned<scalar>&);
699 
700 
701  // Friend Operators
702 
703  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
704  operator& <Type>
705  (
706  const fvMatrix<Type>&,
707  const DimensionedField<Type, volMesh>&
708  );
709 
710  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
711  operator& <Type>
712  (
713  const fvMatrix<Type>&,
714  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
715  );
716 
717  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
718  operator& <Type>
719  (
720  const tmp<fvMatrix<Type>>&,
721  const DimensionedField<Type, volMesh>&
722  );
723 
724  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
725  operator& <Type>
726  (
727  const tmp<fvMatrix<Type>>&,
728  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
729  );
730 
731 
732  // Ostream Operator
733 
734  friend Ostream& operator<< <Type>
735  (
736  Ostream&,
737  const fvMatrix<Type>&
738  );
739 };
740 
741 
742 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
743 
744 template<class Type>
745 void checkMethod
746 (
747  const fvMatrix<Type>&,
748  const fvMatrix<Type>&,
749  const char*
750 );
751 
752 template<class Type>
753 void checkMethod
754 (
755  const fvMatrix<Type>&,
756  const DimensionedField<Type, volMesh>&,
757  const char*
758 );
759 
760 template<class Type>
761 void checkMethod
762 (
763  const fvMatrix<Type>&,
764  const dimensioned<Type>&,
765  const char*
766 );
767 
768 
769 //- Solve returning the solution statistics given convergence tolerance
770 // Use the given solver controls
771 template<class Type>
772 SolverPerformance<Type> solve
773 (
774  fvMatrix<Type>&,
775  const dictionary& solverControls
776 );
777 
778 //- Solve returning the solution statistics given convergence tolerance,
779 //- deleting temporary matrix after solution.
780 // Use the given solver controls
781 template<class Type>
782 SolverPerformance<Type> solve
783 (
784  const tmp<fvMatrix<Type>>&,
785  const dictionary& solverControls
786 );
787 
788 
789 //- Solve returning the solution statistics given convergence tolerance
790 // Uses \p name solver controls from fvSolution
791 template<class Type>
792 SolverPerformance<Type> solve(fvMatrix<Type>&, const word& name);
793 
794 //- Solve returning the solution statistics given convergence tolerance,
795 //- deleting temporary matrix after solution.
796 // Uses \p name solver controls from fvSolution
797 template<class Type>
798 SolverPerformance<Type> solve(const tmp<fvMatrix<Type>>&, const word& name);
799 
800 
801 //- Solve returning the solution statistics given convergence tolerance
802 // Uses solver controls from fvSolution
803 template<class Type>
804 SolverPerformance<Type> solve(fvMatrix<Type>&);
805 
806 //- Solve returning the solution statistics given convergence tolerance,
807 //- deleting temporary matrix after solution.
808 // Uses solver controls from fvSolution
809 template<class Type>
810 SolverPerformance<Type> solve(const tmp<fvMatrix<Type>>&);
811 
812 
813 //- Return the correction form of the given matrix
814 //- by subtracting the matrix multiplied by the current field
815 template<class Type>
816 tmp<fvMatrix<Type>> correction(const fvMatrix<Type>&);
817 
818 
819 //- Return the correction form of the given temporary matrix
820 //- by subtracting the matrix multiplied by the current field
821 template<class Type>
822 tmp<fvMatrix<Type>> correction(const tmp<fvMatrix<Type>>&);
823 
824 
825 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
826 
827 template<class Type>
828 tmp<fvMatrix<Type>> operator==
829 (
830  const fvMatrix<Type>&,
831  const fvMatrix<Type>&
832 );
833 
834 template<class Type>
835 tmp<fvMatrix<Type>> operator==
836 (
837  const tmp<fvMatrix<Type>>&,
838  const fvMatrix<Type>&
839 );
840 
841 template<class Type>
842 tmp<fvMatrix<Type>> operator==
843 (
844  const fvMatrix<Type>&,
845  const tmp<fvMatrix<Type>>&
846 );
847 
848 template<class Type>
849 tmp<fvMatrix<Type>> operator==
850 (
853 );
854 
855 
856 template<class Type>
857 tmp<fvMatrix<Type>> operator==
858 (
859  const fvMatrix<Type>&,
861 );
862 
863 template<class Type>
864 tmp<fvMatrix<Type>> operator==
865 (
866  const fvMatrix<Type>&,
868 );
869 
870 template<class Type>
871 tmp<fvMatrix<Type>> operator==
872 (
873  const fvMatrix<Type>&,
875 );
876 
877 template<class Type>
878 tmp<fvMatrix<Type>> operator==
879 (
880  const tmp<fvMatrix<Type>>&,
882 );
883 
884 template<class Type>
885 tmp<fvMatrix<Type>> operator==
886 (
887  const tmp<fvMatrix<Type>>&,
889 );
890 
891 template<class Type>
892 tmp<fvMatrix<Type>> operator==
893 (
894  const tmp<fvMatrix<Type>>&,
896 );
897 
898 template<class Type>
899 tmp<fvMatrix<Type>> operator==
900 (
901  const fvMatrix<Type>&,
902  const dimensioned<Type>&
903 );
904 
905 template<class Type>
906 tmp<fvMatrix<Type>> operator==
907 (
908  const tmp<fvMatrix<Type>>&,
909  const dimensioned<Type>&
910 );
911 
912 
913 template<class Type>
914 tmp<fvMatrix<Type>> operator==
915 (
916  const fvMatrix<Type>&,
917  const Foam::zero
918 );
919 
920 template<class Type>
921 tmp<fvMatrix<Type>> operator==
922 (
923  const tmp<fvMatrix<Type>>&,
924  const Foam::zero
925 );
926 
927 
928 //- Unary negation
929 template<class Type>
930 tmp<fvMatrix<Type>> operator-
931 (
932  const fvMatrix<Type>&
933 );
934 
935 //- Unary negation
936 template<class Type>
937 tmp<fvMatrix<Type>> operator-
938 (
939  const tmp<fvMatrix<Type>>&
940 );
941 
942 
943 template<class Type>
944 tmp<fvMatrix<Type>> operator+
945 (
946  const fvMatrix<Type>&,
947  const fvMatrix<Type>&
948 );
949 
950 template<class Type>
951 tmp<fvMatrix<Type>> operator+
952 (
953  const tmp<fvMatrix<Type>>&,
954  const fvMatrix<Type>&
955 );
956 
957 template<class Type>
958 tmp<fvMatrix<Type>> operator+
959 (
960  const fvMatrix<Type>&,
961  const tmp<fvMatrix<Type>>&
962 );
963 
964 template<class Type>
965 tmp<fvMatrix<Type>> operator+
966 (
967  const tmp<fvMatrix<Type>>&,
968  const tmp<fvMatrix<Type>>&
969 );
970 
971 
972 template<class Type>
973 tmp<fvMatrix<Type>> operator+
974 (
975  const fvMatrix<Type>&,
977 );
978 
979 template<class Type>
980 tmp<fvMatrix<Type>> operator+
981 (
982  const fvMatrix<Type>&,
984 );
985 
986 template<class Type>
987 tmp<fvMatrix<Type>> operator+
988 (
989  const fvMatrix<Type>&,
991 );
992 
993 template<class Type>
994 tmp<fvMatrix<Type>> operator+
995 (
996  const tmp<fvMatrix<Type>>&,
998 );
999 
1000 template<class Type>
1001 tmp<fvMatrix<Type>> operator+
1002 (
1003  const tmp<fvMatrix<Type>>&,
1005 );
1006 
1007 template<class Type>
1008 tmp<fvMatrix<Type>> operator+
1009 (
1010  const tmp<fvMatrix<Type>>&,
1012 );
1013 
1014 template<class Type>
1015 tmp<fvMatrix<Type>> operator+
1016 (
1018  const fvMatrix<Type>&
1019 );
1020 
1021 template<class Type>
1022 tmp<fvMatrix<Type>> operator+
1023 (
1025  const fvMatrix<Type>&
1026 );
1027 
1028 template<class Type>
1029 tmp<fvMatrix<Type>> operator+
1030 (
1032  const fvMatrix<Type>&
1033 );
1034 
1035 template<class Type>
1036 tmp<fvMatrix<Type>> operator+
1037 (
1039  const tmp<fvMatrix<Type>>&
1040 );
1041 
1042 template<class Type>
1043 tmp<fvMatrix<Type>> operator+
1044 (
1046  const tmp<fvMatrix<Type>>&
1047 );
1048 
1049 template<class Type>
1050 tmp<fvMatrix<Type>> operator+
1051 (
1053  const tmp<fvMatrix<Type>>&
1054 );
1055 
1056 
1057 template<class Type>
1058 tmp<fvMatrix<Type>> operator+
1059 (
1060  const fvMatrix<Type>&,
1061  const dimensioned<Type>&
1062 );
1063 
1064 template<class Type>
1065 tmp<fvMatrix<Type>> operator+
1066 (
1067  const tmp<fvMatrix<Type>>&,
1068  const dimensioned<Type>&
1069 );
1070 
1071 template<class Type>
1072 tmp<fvMatrix<Type>> operator+
1073 (
1074  const dimensioned<Type>&,
1075  const fvMatrix<Type>&
1076 );
1077 
1078 template<class Type>
1079 tmp<fvMatrix<Type>> operator+
1080 (
1081  const dimensioned<Type>&,
1082  const tmp<fvMatrix<Type>>&
1083 );
1084 
1085 
1086 template<class Type>
1087 tmp<fvMatrix<Type>> operator-
1088 (
1089  const fvMatrix<Type>&,
1090  const fvMatrix<Type>&
1091 );
1092 
1093 template<class Type>
1094 tmp<fvMatrix<Type>> operator-
1095 (
1096  const tmp<fvMatrix<Type>>&,
1097  const fvMatrix<Type>&
1098 );
1099 
1100 template<class Type>
1101 tmp<fvMatrix<Type>> operator-
1102 (
1103  const fvMatrix<Type>&,
1104  const tmp<fvMatrix<Type>>&
1105 );
1106 
1107 template<class Type>
1108 tmp<fvMatrix<Type>> operator-
1109 (
1110  const tmp<fvMatrix<Type>>&,
1111  const tmp<fvMatrix<Type>>&
1112 );
1113 
1114 
1115 template<class Type>
1116 tmp<fvMatrix<Type>> operator-
1117 (
1118  const fvMatrix<Type>&,
1120 );
1121 
1122 template<class Type>
1123 tmp<fvMatrix<Type>> operator-
1124 (
1125  const fvMatrix<Type>&,
1127 );
1128 
1129 template<class Type>
1130 tmp<fvMatrix<Type>> operator-
1131 (
1132  const fvMatrix<Type>&,
1134 );
1135 
1136 template<class Type>
1137 tmp<fvMatrix<Type>> operator-
1138 (
1139  const tmp<fvMatrix<Type>>&,
1141 );
1142 
1143 template<class Type>
1144 tmp<fvMatrix<Type>> operator-
1145 (
1146  const tmp<fvMatrix<Type>>&,
1148 );
1149 
1150 template<class Type>
1151 tmp<fvMatrix<Type>> operator-
1152 (
1153  const tmp<fvMatrix<Type>>&,
1155 );
1156 
1157 template<class Type>
1158 tmp<fvMatrix<Type>> operator-
1159 (
1161  const fvMatrix<Type>&
1162 );
1163 
1164 template<class Type>
1165 tmp<fvMatrix<Type>> operator-
1166 (
1168  const fvMatrix<Type>&
1169 );
1170 
1171 template<class Type>
1172 tmp<fvMatrix<Type>> operator-
1173 (
1175  const fvMatrix<Type>&
1176 );
1177 
1178 template<class Type>
1179 tmp<fvMatrix<Type>> operator-
1180 (
1182  const tmp<fvMatrix<Type>>&
1183 );
1184 
1185 template<class Type>
1186 tmp<fvMatrix<Type>> operator-
1187 (
1189  const tmp<fvMatrix<Type>>&
1190 );
1191 
1192 template<class Type>
1193 tmp<fvMatrix<Type>> operator-
1194 (
1196  const tmp<fvMatrix<Type>>&
1197 );
1198 
1199 
1200 template<class Type>
1201 tmp<fvMatrix<Type>> operator-
1202 (
1203  const fvMatrix<Type>&,
1204  const dimensioned<Type>&
1205 );
1206 
1207 template<class Type>
1208 tmp<fvMatrix<Type>> operator-
1209 (
1210  const tmp<fvMatrix<Type>>&,
1211  const dimensioned<Type>&
1212 );
1213 
1214 template<class Type>
1215 tmp<fvMatrix<Type>> operator-
1216 (
1217  const dimensioned<Type>&,
1218  const fvMatrix<Type>&
1219 );
1220 
1221 template<class Type>
1222 tmp<fvMatrix<Type>> operator-
1223 (
1224  const dimensioned<Type>&,
1225  const tmp<fvMatrix<Type>>&
1226 );
1227 
1228 
1229 template<class Type>
1230 tmp<fvMatrix<Type>> operator*
1231 (
1232  const volScalarField::Internal&,
1233  const fvMatrix<Type>&
1234 );
1235 
1236 template<class Type>
1237 tmp<fvMatrix<Type>> operator*
1238 (
1240  const fvMatrix<Type>&
1241 );
1242 
1243 template<class Type>
1244 tmp<fvMatrix<Type>> operator*
1245 (
1246  const tmp<volScalarField>&,
1247  const fvMatrix<Type>&
1248 );
1249 
1250 template<class Type>
1251 tmp<fvMatrix<Type>> operator*
1252 (
1253  const volScalarField::Internal&,
1254  const tmp<fvMatrix<Type>>&
1255 );
1256 
1257 template<class Type>
1258 tmp<fvMatrix<Type>> operator*
1259 (
1261  const tmp<fvMatrix<Type>>&
1262 );
1263 
1264 template<class Type>
1265 tmp<fvMatrix<Type>> operator*
1266 (
1267  const tmp<volScalarField>&,
1268  const tmp<fvMatrix<Type>>&
1269 );
1270 
1271 
1272 template<class Type>
1273 tmp<fvMatrix<Type>> operator*
1274 (
1275  const dimensioned<scalar>&,
1276  const fvMatrix<Type>&
1277 );
1278 
1279 template<class Type>
1280 tmp<fvMatrix<Type>> operator*
1281 (
1282  const dimensioned<scalar>&,
1283  const tmp<fvMatrix<Type>>&
1284 );
1285 
1286 
1287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1288 
1289 } // End namespace Foam
1290 
1291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1292 
1293 #ifdef NoRepository
1294  #include "fvMatrix.C"
1295 #endif
1296 
1297 // Specialisation for scalars
1298 #include "fvScalarMatrix.H"
1299 
1300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1301 
1302 #endif
1303 
1304 // ************************************************************************* //
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:556
Foam::surfaceFields.
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: fvMatrix.C:42
GeometricField< Type, fvsPatchField, surfaceMesh > faceFluxFieldType
Field type for face flux (for non-orthogonal correction)
Definition: fvMatrix.H:133
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
tmp< fvMatrix< Type > > clone() const
Construct and return a clone.
Definition: fvMatrix.H:374
void setValuesFromList(const labelUList &cellLabels, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:220
uint8_t direction
Definition: direction.H:46
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:117
void createOrUpdateLduPrimitiveAssembly()
Create or update ldu assembly.
Definition: fvMatrix.C:914
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1249
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void negate()
Inplace negate.
Definition: fvMatrix.C:1592
const lduPrimitiveMeshAssembly & lduMeshAssembly()
Return optional lduAdressing.
Definition: fvMatrix.H:476
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:145
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:1011
void manipulateMatrix(direction cmp)
Manipulate matrix.
Definition: fvMatrix.C:798
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1332
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.
void operator-=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1645
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Solver class returned by the solver function used for systems in which it is useful to cache the solv...
Definition: fvMatrix.H:296
bool checkImplicit(const label fieldi=0)
Name the implicit assembly addressing.
Definition: fvMatrix.C:315
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Generic dimensioned Type class.
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:1384
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: fvMatrix.C:80
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
const labelListList & patchMap() const
Return patchMap.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
void clear()
Clear multiple fvMatrices.
Definition: fvMatrix.H:521
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
tmp< scalarField > D() const
Return the matrix scalar diagonal.
Definition: fvMatrix.C:1280
void operator=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1547
class FOAM_DEPRECATED_FOR(2017-05, "Foam::Enum") NamedEnum
Definition: NamedEnum.H:65
SolverPerformance< Type > solveCoupled(const dictionary &)
Solve coupled returning the solution statistics.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
GeometricField< Type, fvPatchField, volMesh > psiFieldType
The geometric field type for psi.
Definition: fvMatrix.H:126
label nMatrices() const
Definition: fvMatrix.H:392
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:45
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
SolverPerformance< Type > solve()
Solve returning the solution statistics.
void addFvMatrix(fvMatrix< Type > &matrix)
Add fvMatrix.
Definition: fvMatrix.C:1072
Generic templated field type.
Definition: Field.H:62
void setInterfaces(lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces)
Set interfaces.
Definition: fvMatrix.C:478
A class for handling words, derived from Foam::string.
Definition: word.H:63
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
void setReferences(const labelUList &cellLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:1027
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1043
bool hasFaceFluxCorrection() const noexcept
True if face-flux non-orthogonal correction field exists.
Definition: fvMatrix.H:596
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:464
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
void transferFvMatrixCoeffs()
Transfer lower, upper, diag and source to this fvMatrix.
Definition: fvMatrix.C:821
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
const direction noexcept
Definition: Scalar.H:258
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
Definition: fvMatrix.C:352
fvSolver(fvMatrix< Type > &fvMat, autoPtr< lduMatrix::solver > &&sol)
Definition: fvMatrix.H:306
SolverPerformance< Type > solveSegregatedOrCoupled()
Solve segregated or coupled returning the solution statistics.
void operator*=(const volScalarField::Internal &)
Definition: fvMatrix.C:1768
An assembly of lduMatrix that is specific inter-region coupling through mapped patches.
A dynamically resizable PtrList with allocation management.
Definition: PtrDynList.H:48
faceFluxFieldPtrType & faceFluxCorrectionPtr()
Return pointer to face-flux non-orthogonal correction field.
Definition: fvMatrix.H:588
lduPrimitiveMeshAssembly * lduMeshPtr()
Access to lduPrimitiveMeshAssembly.
Definition: fvMatrix.C:888
GeometricField< Type, fvsPatchField, surfaceMesh > * faceFluxFieldPtrType
Declare return type of the faceFluxCorrectionPtr() function.
Definition: fvMatrix.H:583
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1607
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
ClassName("fvMatrix")
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:1289
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1313
const dimensionSet & dimensions() const noexcept
Definition: fvMatrix.H:528
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
A scalar instance of fvMatrix.
Macro definitions for declaring ClassName(), NamespaceName(), etc.
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
autoPtr< fvSolver > solver()
Construct and return the solver.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:81
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:170
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition: fvMatrix.H:547
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
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.
Definition: fvMatrixSolve.C:31
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1422
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition: fvMatrix.H:565
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated 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:978
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1267
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const fvMatrix< Type > & matrix(const label i) const
Definition: fvMatrix.H:397
scalarField & diag()
Definition: lduMatrix.C:197
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< Field< Type > > residual() const
Return the matrix residual.
SolverPerformance< Type > solve()
Solve returning the solution statistics.
void setBounAndInterCoeffs()
Manipulate boundary/internal coeffs for coupling.
Definition: fvMatrix.C:668
label globalPatchID(const label fieldi, const label patchi) const
Definition: fvMatrix.H:408
Namespace for OpenFOAM.
const dictionary & solverDict() const
Return the solver dictionary for psi, taking into account finalIteration.
Definition: fvMatrix.C:1535