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