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 std::unique_ptr<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
497  typedef std::unique_ptr<faceFluxFieldType> faceFluxFieldPtrType;
498 
499  //- Return pointer to face-flux non-orthogonal correction field
501  {
502  return faceFluxCorrectionPtr_;
503  }
504 
505  //- Set pointer to face-flux non-orthogonal correction field
507  {
508  faceFluxCorrectionPtr_.reset(flux);
509  }
510 
511  //- True if face-flux non-orthogonal correction field exists
512  bool hasFaceFluxCorrection() const noexcept
513  {
514  return bool(faceFluxCorrectionPtr_);
515  }
516 
517 
518  // Operations
519 
520  //- Set solution in given cells to the specified value
521  //- and eliminate the corresponding equations from the matrix.
522  void setValues
523  (
524  const labelUList& cellLabels,
525  const Type& value
526  );
527 
528  //- Set solution in given cells to the specified values
529  //- and eliminate the corresponding equations from the matrix.
530  void setValues
531  (
532  const labelUList& cellLabels,
534  );
535 
536  //- Set solution in given cells to the specified values
537  //- and eliminate the corresponding equations from the matrix.
539  (
540  const labelUList& cellLabels,
542  );
543 
544  //- Set reference level for solution
545  void setReference
546  (
547  const label celli,
548  const Type& value,
549  const bool forceReference = false
550  );
551 
552  //- Set reference level for solution
553  void setReferences
554  (
555  const labelUList& cellLabels,
556  const Type& value,
557  const bool forceReference = false
558  );
559 
560  //- Set reference level for solution
561  void setReferences
562  (
563  const labelUList& cellLabels,
564  const UList<Type>& values,
565  const bool forceReference = false
566  );
567 
568  //- Set reference level for a component of the solution
569  //- on a given patch face
571  (
572  const label patchi,
573  const label facei,
574  const direction cmpt,
575  const scalar value
576  );
577 
578  //- Add fvMatrix
580 
581  //- Relax matrix (for steady-state solution).
582  // alpha = 1 : diagonally equal
583  // alpha < 1 : diagonally dominant
584  // alpha = 0 : do nothing
585  // Note: Requires positive diagonal.
586  void relax(const scalar alpha);
588  //- Relax matrix (for steady-state solution).
589  // alpha is read from controlDict
590  void relax();
591 
592  //- Manipulate based on a boundary field
593  void boundaryManipulate
594  (
596  Boundary& values
597  );
598 
599  //- Construct and return the solver
600  // Use the given solver controls
602 
603  //- Construct and return the solver
604  // Solver controls read from fvSolution
606 
607  //- Solve segregated or coupled returning the solution statistics.
608  // Use the given solver controls
610 
611  //- Solve segregated or coupled returning the solution statistics.
612  // Solver controls read from fvSolution
614 
615  //- Solve segregated returning the solution statistics.
616  // Use the given solver controls
618 
619  //- Solve coupled returning the solution statistics.
620  // Use the given solver controls
622 
623  //- Solve returning the solution statistics.
624  // Use the given solver controls
626 
627  //- Solve returning the solution statistics.
628  // Uses \p name solver controls from fvSolution
630 
631  //- Solve returning the solution statistics.
632  // Solver controls read from fvSolution
634 
635  //- Return the matrix residual
636  tmp<Field<Type>> residual() const;
637 
638  //- Return the matrix scalar diagonal
639  tmp<scalarField> D() const;
640 
641  //- Return the matrix Type diagonal
642  tmp<Field<Type>> DD() const;
643 
644  //- Return the central coefficient
645  tmp<volScalarField> A() const;
646 
647  //- Return the H operation source
649 
650  //- Return H(1)
651  tmp<volScalarField> H1() const;
652 
653  //- Return the face-flux field from the matrix
655  flux() const;
656 
657  //- Return the solver dictionary (from fvSolution) for \p name
658  const dictionary& solverDict(const word& name) const;
659 
660  //- Return the solver dictionary for psi,
661  //- taking into account finalIteration
662  const dictionary& solverDict() const;
663 
664 
665  // Member Operators
666 
667  void operator=(const fvMatrix<Type>&);
668  void operator=(const tmp<fvMatrix<Type>>&);
669 
670  //- Inplace negate
671  void negate();
672 
673  void operator+=(const fvMatrix<Type>&);
674  void operator+=(const tmp<fvMatrix<Type>>&);
675 
676  void operator-=(const fvMatrix<Type>&);
677  void operator-=(const tmp<fvMatrix<Type>>&);
678 
681  void operator+=
682  (
684  );
685 
688  void operator-=
689  (
691  );
692 
693  void operator+=(const dimensioned<Type>&);
694  void operator-=(const dimensioned<Type>&);
695 
696  void operator+=(const Foam::zero) {}
697  void operator-=(const Foam::zero) {}
698 
700  void operator*=(const tmp<volScalarField::Internal>&);
701  void operator*=(const tmp<volScalarField>&);
702 
703  void operator*=(const dimensioned<scalar>&);
704 
705 
706  // Friend Operators
707 
708  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
709  operator& <Type>
710  (
711  const fvMatrix<Type>&,
712  const DimensionedField<Type, volMesh>&
713  );
714 
715  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
716  operator& <Type>
717  (
718  const fvMatrix<Type>&,
719  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
720  );
721 
722  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
723  operator& <Type>
724  (
725  const tmp<fvMatrix<Type>>&,
726  const DimensionedField<Type, volMesh>&
727  );
728 
729  friend tmp<GeometricField<Type, fvPatchField, volMesh>>
730  operator& <Type>
731  (
732  const tmp<fvMatrix<Type>>&,
733  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
734  );
735 
736 
737  // Ostream Operator
738 
739  friend Ostream& operator<< <Type>
740  (
741  Ostream&,
742  const fvMatrix<Type>&
743  );
744 };
745 
746 
747 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
748 
749 template<class Type>
750 void checkMethod
751 (
752  const fvMatrix<Type>&,
753  const fvMatrix<Type>&,
754  const char*
755 );
756 
757 template<class Type>
758 void checkMethod
759 (
760  const fvMatrix<Type>&,
761  const DimensionedField<Type, volMesh>&,
762  const char*
763 );
764 
765 template<class Type>
766 void checkMethod
767 (
768  const fvMatrix<Type>&,
769  const dimensioned<Type>&,
770  const char*
771 );
772 
773 
774 //- Solve returning the solution statistics given convergence tolerance
775 // Use the given solver controls
776 template<class Type>
777 SolverPerformance<Type> solve
778 (
779  fvMatrix<Type>&,
780  const dictionary& solverControls
781 );
782 
783 //- Solve returning the solution statistics given convergence tolerance,
784 //- deleting temporary matrix after solution.
785 // Use the given solver controls
786 template<class Type>
787 SolverPerformance<Type> solve
788 (
789  const tmp<fvMatrix<Type>>&,
790  const dictionary& solverControls
791 );
792 
793 
794 //- Solve returning the solution statistics given convergence tolerance
795 // Uses \p name solver controls from fvSolution
796 template<class Type>
797 SolverPerformance<Type> solve(fvMatrix<Type>&, const word& name);
798 
799 //- Solve returning the solution statistics given convergence tolerance,
800 //- deleting temporary matrix after solution.
801 // Uses \p name solver controls from fvSolution
802 template<class Type>
803 SolverPerformance<Type> solve(const tmp<fvMatrix<Type>>&, const word& name);
804 
805 
806 //- Solve returning the solution statistics given convergence tolerance
807 // Uses solver controls from fvSolution
808 template<class Type>
809 SolverPerformance<Type> solve(fvMatrix<Type>&);
810 
811 //- Solve returning the solution statistics given convergence tolerance,
812 //- deleting temporary matrix after solution.
813 // Uses solver controls from fvSolution
814 template<class Type>
815 SolverPerformance<Type> solve(const tmp<fvMatrix<Type>>&);
816 
817 
818 //- Return the correction form of the given matrix
819 //- by subtracting the matrix multiplied by the current field
820 template<class Type>
821 tmp<fvMatrix<Type>> correction(const fvMatrix<Type>&);
822 
823 
824 //- Return the correction form of the given temporary matrix
825 //- by subtracting the matrix multiplied by the current field
826 template<class Type>
827 tmp<fvMatrix<Type>> correction(const tmp<fvMatrix<Type>>&);
828 
829 
830 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
831 
832 template<class Type>
833 tmp<fvMatrix<Type>> operator==
834 (
835  const fvMatrix<Type>&,
836  const fvMatrix<Type>&
837 );
838 
839 template<class Type>
840 tmp<fvMatrix<Type>> operator==
841 (
842  const tmp<fvMatrix<Type>>&,
843  const fvMatrix<Type>&
844 );
845 
846 template<class Type>
847 tmp<fvMatrix<Type>> operator==
848 (
849  const fvMatrix<Type>&,
850  const tmp<fvMatrix<Type>>&
851 );
852 
853 template<class Type>
854 tmp<fvMatrix<Type>> operator==
855 (
856  const tmp<fvMatrix<Type>>&,
857  const tmp<fvMatrix<Type>>&
858 );
860 
861 template<class Type>
862 tmp<fvMatrix<Type>> operator==
863 (
864  const fvMatrix<Type>&,
866 );
867 
868 template<class Type>
869 tmp<fvMatrix<Type>> operator==
870 (
871  const fvMatrix<Type>&,
873 );
874 
875 template<class Type>
876 tmp<fvMatrix<Type>> operator==
877 (
878  const fvMatrix<Type>&,
880 );
881 
882 template<class Type>
883 tmp<fvMatrix<Type>> operator==
884 (
885  const tmp<fvMatrix<Type>>&,
887 );
888 
889 template<class Type>
890 tmp<fvMatrix<Type>> operator==
891 (
892  const tmp<fvMatrix<Type>>&,
894 );
895 
896 template<class Type>
897 tmp<fvMatrix<Type>> operator==
898 (
899  const tmp<fvMatrix<Type>>&,
901 );
902 
903 template<class Type>
904 tmp<fvMatrix<Type>> operator==
905 (
906  const fvMatrix<Type>&,
907  const dimensioned<Type>&
908 );
909 
910 template<class Type>
911 tmp<fvMatrix<Type>> operator==
912 (
913  const tmp<fvMatrix<Type>>&,
914  const dimensioned<Type>&
915 );
916 
917 
918 template<class Type>
919 tmp<fvMatrix<Type>> operator==
920 (
921  const fvMatrix<Type>&,
922  const Foam::zero
923 );
924 
925 template<class Type>
926 tmp<fvMatrix<Type>> operator==
927 (
928  const tmp<fvMatrix<Type>>&,
929  const Foam::zero
930 );
931 
932 
933 //- Unary negation
934 template<class Type>
935 tmp<fvMatrix<Type>> operator-
936 (
937  const fvMatrix<Type>&
938 );
939 
940 //- Unary negation
941 template<class Type>
942 tmp<fvMatrix<Type>> operator-
943 (
944  const tmp<fvMatrix<Type>>&
945 );
946 
947 
948 template<class Type>
949 tmp<fvMatrix<Type>> operator+
950 (
951  const fvMatrix<Type>&,
952  const fvMatrix<Type>&
953 );
954 
955 template<class Type>
956 tmp<fvMatrix<Type>> operator+
957 (
958  const tmp<fvMatrix<Type>>&,
959  const fvMatrix<Type>&
960 );
961 
962 template<class Type>
963 tmp<fvMatrix<Type>> operator+
964 (
965  const fvMatrix<Type>&,
966  const tmp<fvMatrix<Type>>&
967 );
968 
969 template<class Type>
970 tmp<fvMatrix<Type>> operator+
971 (
972  const tmp<fvMatrix<Type>>&,
973  const tmp<fvMatrix<Type>>&
974 );
975 
976 
977 template<class Type>
978 tmp<fvMatrix<Type>> operator+
979 (
980  const fvMatrix<Type>&,
982 );
983 
984 template<class Type>
985 tmp<fvMatrix<Type>> operator+
986 (
987  const fvMatrix<Type>&,
989 );
990 
991 template<class Type>
992 tmp<fvMatrix<Type>> operator+
993 (
994  const fvMatrix<Type>&,
996 );
997 
998 template<class Type>
999 tmp<fvMatrix<Type>> operator+
1000 (
1001  const tmp<fvMatrix<Type>>&,
1003 );
1004 
1005 template<class Type>
1006 tmp<fvMatrix<Type>> operator+
1007 (
1008  const tmp<fvMatrix<Type>>&,
1010 );
1011 
1012 template<class Type>
1013 tmp<fvMatrix<Type>> operator+
1014 (
1015  const tmp<fvMatrix<Type>>&,
1017 );
1018 
1019 template<class Type>
1020 tmp<fvMatrix<Type>> operator+
1021 (
1023  const fvMatrix<Type>&
1024 );
1025 
1026 template<class Type>
1027 tmp<fvMatrix<Type>> operator+
1028 (
1030  const fvMatrix<Type>&
1031 );
1032 
1033 template<class Type>
1034 tmp<fvMatrix<Type>> operator+
1035 (
1037  const fvMatrix<Type>&
1038 );
1039 
1040 template<class Type>
1041 tmp<fvMatrix<Type>> operator+
1042 (
1044  const tmp<fvMatrix<Type>>&
1045 );
1046 
1047 template<class Type>
1048 tmp<fvMatrix<Type>> operator+
1049 (
1051  const tmp<fvMatrix<Type>>&
1052 );
1053 
1054 template<class Type>
1055 tmp<fvMatrix<Type>> operator+
1056 (
1058  const tmp<fvMatrix<Type>>&
1059 );
1060 
1061 
1062 template<class Type>
1063 tmp<fvMatrix<Type>> operator+
1064 (
1065  const fvMatrix<Type>&,
1066  const dimensioned<Type>&
1067 );
1068 
1069 template<class Type>
1070 tmp<fvMatrix<Type>> operator+
1071 (
1072  const tmp<fvMatrix<Type>>&,
1073  const dimensioned<Type>&
1074 );
1075 
1076 template<class Type>
1077 tmp<fvMatrix<Type>> operator+
1078 (
1079  const dimensioned<Type>&,
1080  const fvMatrix<Type>&
1081 );
1082 
1083 template<class Type>
1084 tmp<fvMatrix<Type>> operator+
1085 (
1086  const dimensioned<Type>&,
1087  const tmp<fvMatrix<Type>>&
1088 );
1089 
1090 
1091 template<class Type>
1092 tmp<fvMatrix<Type>> operator-
1093 (
1094  const fvMatrix<Type>&,
1095  const fvMatrix<Type>&
1096 );
1097 
1098 template<class Type>
1099 tmp<fvMatrix<Type>> operator-
1100 (
1101  const tmp<fvMatrix<Type>>&,
1102  const fvMatrix<Type>&
1103 );
1104 
1105 template<class Type>
1106 tmp<fvMatrix<Type>> operator-
1107 (
1108  const fvMatrix<Type>&,
1109  const tmp<fvMatrix<Type>>&
1110 );
1111 
1112 template<class Type>
1113 tmp<fvMatrix<Type>> operator-
1114 (
1115  const tmp<fvMatrix<Type>>&,
1116  const tmp<fvMatrix<Type>>&
1117 );
1118 
1119 
1120 template<class Type>
1121 tmp<fvMatrix<Type>> operator-
1122 (
1123  const fvMatrix<Type>&,
1125 );
1126 
1127 template<class Type>
1128 tmp<fvMatrix<Type>> operator-
1129 (
1130  const fvMatrix<Type>&,
1132 );
1133 
1134 template<class Type>
1135 tmp<fvMatrix<Type>> operator-
1136 (
1137  const fvMatrix<Type>&,
1139 );
1140 
1141 template<class Type>
1142 tmp<fvMatrix<Type>> operator-
1143 (
1144  const tmp<fvMatrix<Type>>&,
1146 );
1147 
1148 template<class Type>
1149 tmp<fvMatrix<Type>> operator-
1150 (
1151  const tmp<fvMatrix<Type>>&,
1153 );
1154 
1155 template<class Type>
1156 tmp<fvMatrix<Type>> operator-
1157 (
1158  const tmp<fvMatrix<Type>>&,
1160 );
1161 
1162 template<class Type>
1163 tmp<fvMatrix<Type>> operator-
1164 (
1166  const fvMatrix<Type>&
1167 );
1168 
1169 template<class Type>
1170 tmp<fvMatrix<Type>> operator-
1171 (
1173  const fvMatrix<Type>&
1174 );
1175 
1176 template<class Type>
1177 tmp<fvMatrix<Type>> operator-
1178 (
1180  const fvMatrix<Type>&
1181 );
1182 
1183 template<class Type>
1184 tmp<fvMatrix<Type>> operator-
1185 (
1187  const tmp<fvMatrix<Type>>&
1188 );
1189 
1190 template<class Type>
1191 tmp<fvMatrix<Type>> operator-
1192 (
1194  const tmp<fvMatrix<Type>>&
1195 );
1196 
1197 template<class Type>
1198 tmp<fvMatrix<Type>> operator-
1199 (
1201  const tmp<fvMatrix<Type>>&
1202 );
1203 
1204 
1205 template<class Type>
1206 tmp<fvMatrix<Type>> operator-
1207 (
1208  const fvMatrix<Type>&,
1209  const dimensioned<Type>&
1210 );
1211 
1212 template<class Type>
1213 tmp<fvMatrix<Type>> operator-
1214 (
1215  const tmp<fvMatrix<Type>>&,
1216  const dimensioned<Type>&
1217 );
1218 
1219 template<class Type>
1220 tmp<fvMatrix<Type>> operator-
1221 (
1222  const dimensioned<Type>&,
1223  const fvMatrix<Type>&
1224 );
1225 
1226 template<class Type>
1227 tmp<fvMatrix<Type>> operator-
1228 (
1229  const dimensioned<Type>&,
1230  const tmp<fvMatrix<Type>>&
1231 );
1232 
1233 
1234 template<class Type>
1235 tmp<fvMatrix<Type>> operator*
1236 (
1237  const volScalarField::Internal&,
1238  const fvMatrix<Type>&
1239 );
1240 
1241 template<class Type>
1242 tmp<fvMatrix<Type>> operator*
1243 (
1245  const fvMatrix<Type>&
1246 );
1247 
1248 template<class Type>
1249 tmp<fvMatrix<Type>> operator*
1250 (
1251  const tmp<volScalarField>&,
1252  const fvMatrix<Type>&
1253 );
1254 
1255 template<class Type>
1256 tmp<fvMatrix<Type>> operator*
1257 (
1258  const volScalarField::Internal&,
1259  const tmp<fvMatrix<Type>>&
1260 );
1261 
1262 template<class Type>
1263 tmp<fvMatrix<Type>> operator*
1264 (
1266  const tmp<fvMatrix<Type>>&
1267 );
1268 
1269 template<class Type>
1270 tmp<fvMatrix<Type>> operator*
1271 (
1272  const tmp<volScalarField>&,
1273  const tmp<fvMatrix<Type>>&
1274 );
1275 
1276 
1277 template<class Type>
1278 tmp<fvMatrix<Type>> operator*
1279 (
1280  const dimensioned<scalar>&,
1281  const fvMatrix<Type>&
1282 );
1283 
1284 template<class Type>
1285 tmp<fvMatrix<Type>> operator*
1286 (
1287  const dimensioned<scalar>&,
1288  const tmp<fvMatrix<Type>>&
1289 );
1290 
1291 
1292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1293 
1294 } // End namespace Foam
1295 
1296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1297 
1298 #ifdef NoRepository
1299  #include "fvMatrix.C"
1300 #endif
1301 
1302 // Specialisation for scalars
1303 #include "fvScalarMatrix.H"
1304 
1305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1306 
1307 #endif
1308 
1309 // ************************************************************************* //
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:550
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...
const scalarField & diag() const
Definition: lduMatrix.C:163
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:908
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1243
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:1587
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:1005
void manipulateMatrix(direction cmp)
Manipulate matrix.
Definition: fvMatrix.C:792
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1326
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:1639
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.
Generic dimensioned Type class.
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:1378
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.
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:1274
void operator=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1541
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:1066
Generic templated field type.
Definition: Field.H:62
void setInterfaces(lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces)
Set interfaces.
Definition: fvMatrix.C:472
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:1021
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1034
bool hasFaceFluxCorrection() const noexcept
True if face-flux non-orthogonal correction field exists.
Definition: fvMatrix.H:603
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:459
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:815
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:1763
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:587
lduPrimitiveMeshAssembly * lduMeshPtr()
Access to lduPrimitiveMeshAssembly.
Definition: fvMatrix.C:882
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1602
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
ClassName("fvMatrix")
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:1283
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1307
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:52
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...
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:1416
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:972
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1261
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const fvMatrix< Type > & matrix(const label i) const
Definition: fvMatrix.H:397
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:662
label globalPatchID(const label fieldi, const label patchi) const
Definition: fvMatrix.H:408
std::unique_ptr< faceFluxFieldType > faceFluxFieldPtrType
Declare return type of the faceFluxCorrectionPtr() function.
Definition: fvMatrix.H:582
Namespace for OpenFOAM.
const dictionary & solverDict() const
Return the solver dictionary for psi, taking into account finalIteration.
Definition: fvMatrix.C:1529