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-2020 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"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 // Forward Declarations
61 template<class Type> class fvMatrix;
62 template<class T> class UIndirectList;
63 
64 template<class Type>
65 Ostream& operator<<(Ostream&, const fvMatrix<Type>&);
66 
67 
68 template<class Type>
69 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
70 (
71  const fvMatrix<Type>&,
72  const DimensionedField<Type, volMesh>&
73 );
74 
75 template<class Type>
76 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
77 (
78  const fvMatrix<Type>&,
79  const tmp<DimensionedField<Type, volMesh>>&
80 );
81 
82 template<class Type>
83 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
84 (
85  const fvMatrix<Type>&,
86  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
87 );
88 
89 template<class Type>
90 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
91 (
92  const tmp<fvMatrix<Type>>&,
93  const DimensionedField<Type, volMesh>&
94 );
95 
96 template<class Type>
97 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
98 (
99  const tmp<fvMatrix<Type>>&,
100  const tmp<DimensionedField<Type, volMesh>>&
101 );
102 
103 template<class Type>
104 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
105 (
106  const tmp<fvMatrix<Type>>&,
107  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
108 );
109 
110 
111 /*---------------------------------------------------------------------------*\
112  Class fvMatrix Declaration
113 \*---------------------------------------------------------------------------*/
114 
115 template<class Type>
116 class fvMatrix
117 :
118  public refCount,
119  public lduMatrix
120 {
121  // Private Data
122 
123  //- Const reference to field
124  // Converted into a non-const reference at the point of solution.
125  const GeometricField<Type, fvPatchField, volMesh>& psi_;
126 
127  //- Dimension set
128  dimensionSet dimensions_;
129 
130  //- Source term
131  Field<Type> source_;
132 
133  //- Boundary scalar field containing pseudo-matrix coeffs
134  //- for internal cells
135  FieldField<Field, Type> internalCoeffs_;
136 
137  //- Boundary scalar field containing pseudo-matrix coeffs
138  //- for boundary cells
139  FieldField<Field, Type> boundaryCoeffs_;
140 
141  //- Face flux field for non-orthogonal correction
142  mutable GeometricField<Type, fvsPatchField, surfaceMesh>
143  *faceFluxCorrectionPtr_;
144 
145 
146 protected:
147 
148  //- Declare friendship with the fvSolver class
149  friend class fvSolver;
150 
151  // Protected Member Functions
152 
153  //- Add patch contribution to internal field
154  template<class Type2>
155  void addToInternalField
156  (
157  const labelUList& addr,
158  const Field<Type2>& pf,
159  Field<Type2>& intf
160  ) const;
161 
162  template<class Type2>
163  void addToInternalField
164  (
165  const labelUList& addr,
166  const tmp<Field<Type2>>& tpf,
167  Field<Type2>& intf
168  ) const;
169 
170  //- Subtract patch contribution from internal field
171  template<class Type2>
173  (
174  const labelUList& addr,
175  const Field<Type2>& pf,
176  Field<Type2>& intf
177  ) const;
178 
179  template<class Type2>
181  (
182  const labelUList& addr,
183  const tmp<Field<Type2>>& tpf,
184  Field<Type2>& intf
185  ) const;
186 
187 
188  // Matrix completion functionality
189 
190  void addBoundaryDiag
191  (
192  scalarField& diag,
193  const direction cmpt
194  ) const;
195 
197 
198  void addBoundarySource
199  (
201  const bool couples=true
202  ) const;
203 
204 
205  // Matrix manipulation functionality
206 
207  //- Set solution in given cells to the specified values
208  template<template<class> class ListType>
209  void setValuesFromList
210  (
211  const labelUList& cellLabels,
212  const ListType<Type>& values
213  );
214 
215 
216 public:
217 
218  //- Solver class returned by the solver function
219  //- used for systems in which it is useful to cache the solver for reuse.
220  // E.g. if the solver is potentially expensive to construct (AMG) and can
221  // be used several times (PISO)
222  class fvSolver
223  {
224  fvMatrix<Type>& fvMat_;
225 
227 
228  public:
229 
230  // Constructors
231 
233  :
234  fvMat_(fvMat),
235  solver_(std::move(sol))
236  {}
237 
238 
239  // Member Functions
240 
241  //- Solve returning the solution statistics.
242  // Solver controls read from dictionary
243  SolverPerformance<Type> solve(const dictionary& solverControls);
244 
245  //- Solve returning the solution statistics.
246  // Solver controls read from fvSolution
248  };
249 
250 
251  ClassName("fvMatrix");
252 
253 
254  // Constructors
255 
256  //- Construct given a field to solve for
257  fvMatrix
258  (
260  const dimensionSet& ds
261  );
262 
263  //- Copy construct
264  fvMatrix(const fvMatrix<Type>&);
265 
266  //- Copy/move construct from tmp<fvMatrix<Type>>
267  fvMatrix(const tmp<fvMatrix<Type>>&);
268 
269  //- Construct from Istream given field to solve for
270  fvMatrix
271  (
273  Istream& is
274  );
275 
276  //- Clone
277  tmp<fvMatrix<Type>> clone() const;
278 
279 
280  //- Destructor
281  virtual ~fvMatrix();
282 
283 
284  // Member Functions
285 
286  // Access
287 
289  {
290  return psi_;
291  }
292 
293  const dimensionSet& dimensions() const
294  {
295  return dimensions_;
296  }
297 
299  {
300  return source_;
301  }
302 
303  const Field<Type>& source() const
304  {
305  return source_;
306  }
307 
308  //- fvBoundary scalar field containing pseudo-matrix coeffs
309  //- for internal cells
311  {
312  return internalCoeffs_;
313  }
314 
315  //- fvBoundary scalar field containing pseudo-matrix coeffs
316  //- for internal cells
318  {
319  return internalCoeffs_;
320  }
321 
322  //- fvBoundary scalar field containing pseudo-matrix coeffs
323  //- for boundary cells
325  {
326  return boundaryCoeffs_;
327  }
328 
329  //- fvBoundary scalar field containing pseudo-matrix coeffs
330  //- for boundary cells
332  {
333  return boundaryCoeffs_;
334  }
335 
336 
337  //- Declare return type of the faceFluxCorrectionPtr() function
340 
341  //- Return pointer to face-flux non-orthogonal correction field
343  {
344  return faceFluxCorrectionPtr_;
345  }
346 
347 
348  // Operations
349 
350  //- Set solution in given cells to the specified value
351  //- and eliminate the corresponding equations from the matrix.
352  void setValues
353  (
354  const labelUList& cellLabels,
355  const Type& value
356  );
357 
358  //- Set solution in given cells to the specified values
359  //- and eliminate the corresponding equations from the matrix.
360  void setValues
361  (
362  const labelUList& cellLabels,
363  const UList<Type>& values
364  );
365 
366  //- Set solution in given cells to the specified values
367  //- and eliminate the corresponding equations from the matrix.
368  void setValues
369  (
370  const labelUList& cellLabels,
372  );
373 
374  //- Set reference level for solution
375  void setReference
376  (
377  const label celli,
378  const Type& value,
379  const bool forceReference = false
380  );
381 
382  //- Set reference level for solution
383  void setReferences
384  (
385  const labelUList& cellLabels,
386  const Type& value,
387  const bool forceReference = false
388  );
389 
390  //- Set reference level for solution
391  void setReferences
392  (
393  const labelUList& cellLabels,
394  const UList<Type>& values,
395  const bool forceReference = false
396  );
397 
398  //- Set reference level for a component of the solution
399  //- on a given patch face
401  (
402  const label patchi,
403  const label facei,
404  const direction cmpt,
405  const scalar value
406  );
407 
408  //- Relax matrix (for steady-state solution).
409  // alpha = 1 : diagonally equal
410  // alpha < 1 : diagonally dominant
411  // alpha = 0 : do nothing
412  // Note: Requires positive diagonal.
413  void relax(const scalar alpha);
414 
415  //- Relax matrix (for steady-state solution).
416  // alpha is read from controlDict
417  void relax();
418 
419  //- Manipulate based on a boundary field
420  void boundaryManipulate
421  (
423  Boundary& values
424  );
425 
426  //- Construct and return the solver
427  // Use the given solver controls
429 
430  //- Construct and return the solver
431  // Solver controls read from fvSolution
433 
434  //- Solve segregated or coupled returning the solution statistics.
435  // Use the given solver controls
437 
438  //- Solve segregated returning the solution statistics.
439  // Use the given solver controls
441 
442  //- Solve coupled returning the solution statistics.
443  // Use the given solver controls
445 
446  //- Solve returning the solution statistics.
447  // Use the given solver controls
449 
450  //- Solve returning the solution statistics.
451  // Solver controls read from fvSolution
453 
454  //- Return the matrix residual
455  tmp<Field<Type>> residual() const;
456 
457  //- Return the matrix scalar diagonal
458  tmp<scalarField> D() const;
459 
460  //- Return the matrix Type diagonal
461  tmp<Field<Type>> DD() const;
462 
463  //- Return the central coefficient
464  tmp<volScalarField> A() const;
465 
466  //- Return the H operation source
468 
469  //- Return H(1)
470  tmp<volScalarField> H1() const;
471 
472  //- Return the face-flux field from the matrix
474  flux() const;
475 
476  //- Return the solver dictionary taking into account finalIteration
477  const dictionary& solverDict() const;
478 
479 
480  // Member Operators
481 
482  void operator=(const fvMatrix<Type>&);
483  void operator=(const tmp<fvMatrix<Type>>&);
484 
485  void negate();
486 
487  void operator+=(const fvMatrix<Type>&);
488  void operator+=(const tmp<fvMatrix<Type>>&);
489 
490  void operator-=(const fvMatrix<Type>&);
491  void operator-=(const tmp<fvMatrix<Type>>&);
492 
493  void operator+=
494  (
496  );
497  void operator+=
498  (
500  );
501  void operator+=
502  (
504  );
505 
506  void operator-=
507  (
509  );
510  void operator-=
511  (
513  );
514  void operator-=
515  (
517  );
518 
519  void operator+=(const dimensioned<Type>&);
520  void operator-=(const dimensioned<Type>&);
521 
522  void operator+=(const zero&);
523  void operator-=(const zero&);
524 
527  void operator*=(const tmp<volScalarField>&);
528 
529  void operator*=(const dimensioned<scalar>&);
530 
531 
532  // Friend Operators
533 
535  operator& <Type>
536  (
537  const fvMatrix<Type>&,
539  );
540 
542  operator& <Type>
543  (
544  const fvMatrix<Type>&,
546  );
547 
549  operator& <Type>
550  (
551  const tmp<fvMatrix<Type>>&,
553  );
554 
556  operator& <Type>
557  (
558  const tmp<fvMatrix<Type>>&,
560  );
561 
562 
563  // Ostream Operator
564 
565  friend Ostream& operator<< <Type>
566  (
567  Ostream&,
568  const fvMatrix<Type>&
569  );
570 };
571 
572 
573 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
574 
575 template<class Type>
576 void checkMethod
577 (
578  const fvMatrix<Type>&,
579  const fvMatrix<Type>&,
580  const char*
581 );
582 
583 template<class Type>
584 void checkMethod
585 (
586  const fvMatrix<Type>&,
587  const DimensionedField<Type, volMesh>&,
588  const char*
589 );
590 
591 template<class Type>
592 void checkMethod
593 (
594  const fvMatrix<Type>&,
595  const dimensioned<Type>&,
596  const char*
597 );
598 
599 
600 //- Solve returning the solution statistics given convergence tolerance
601 // Use the given solver controls
602 template<class Type>
603 SolverPerformance<Type> solve(fvMatrix<Type>&, const dictionary&);
604 
605 
606 //- Solve returning the solution statistics given convergence tolerance,
607 //- deleting temporary matrix after solution.
608 // Use the given solver controls
609 template<class Type>
610 SolverPerformance<Type> solve
611 (
612  const tmp<fvMatrix<Type>>&,
613  const dictionary&
614 );
615 
616 
617 //- Solve returning the solution statistics given convergence tolerance
618 // Solver controls read fvSolution
619 template<class Type>
620 SolverPerformance<Type> solve(fvMatrix<Type>&);
621 
622 
623 //- Solve returning the solution statistics given convergence tolerance,
624 //- deleting temporary matrix after solution.
625 // Solver controls read fvSolution
626 template<class Type>
627 SolverPerformance<Type> solve(const tmp<fvMatrix<Type>>&);
628 
629 
630 //- Return the correction form of the given matrix
631 //- by subtracting the matrix multiplied by the current field
632 template<class Type>
633 tmp<fvMatrix<Type>> correction(const fvMatrix<Type>&);
634 
635 
636 //- Return the correction form of the given temporary matrix
637 //- by subtracting the matrix multiplied by the current field
638 template<class Type>
639 tmp<fvMatrix<Type>> correction(const tmp<fvMatrix<Type>>&);
640 
641 
642 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
643 
644 template<class Type>
645 tmp<fvMatrix<Type>> operator==
646 (
647  const fvMatrix<Type>&,
648  const fvMatrix<Type>&
649 );
650 
651 template<class Type>
652 tmp<fvMatrix<Type>> operator==
653 (
654  const tmp<fvMatrix<Type>>&,
655  const fvMatrix<Type>&
656 );
657 
658 template<class Type>
659 tmp<fvMatrix<Type>> operator==
660 (
661  const fvMatrix<Type>&,
662  const tmp<fvMatrix<Type>>&
663 );
664 
665 template<class Type>
666 tmp<fvMatrix<Type>> operator==
667 (
668  const tmp<fvMatrix<Type>>&,
669  const tmp<fvMatrix<Type>>&
670 );
671 
672 
673 template<class Type>
674 tmp<fvMatrix<Type>> operator==
675 (
676  const fvMatrix<Type>&,
677  const DimensionedField<Type, volMesh>&
678 );
679 
680 template<class Type>
681 tmp<fvMatrix<Type>> operator==
682 (
683  const fvMatrix<Type>&,
684  const tmp<DimensionedField<Type, volMesh>>&
685 );
686 
687 template<class Type>
688 tmp<fvMatrix<Type>> operator==
689 (
690  const fvMatrix<Type>&,
691  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
692 );
693 
694 template<class Type>
695 tmp<fvMatrix<Type>> operator==
696 (
697  const tmp<fvMatrix<Type>>&,
698  const DimensionedField<Type, volMesh>&
699 );
700 
701 template<class Type>
702 tmp<fvMatrix<Type>> operator==
703 (
704  const tmp<fvMatrix<Type>>&,
705  const tmp<DimensionedField<Type, volMesh>>&
706 );
707 
708 template<class Type>
709 tmp<fvMatrix<Type>> operator==
710 (
711  const tmp<fvMatrix<Type>>&,
712  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
713 );
714 
715 template<class Type>
716 tmp<fvMatrix<Type>> operator==
717 (
718  const fvMatrix<Type>&,
719  const dimensioned<Type>&
720 );
721 
722 template<class Type>
723 tmp<fvMatrix<Type>> operator==
724 (
725  const tmp<fvMatrix<Type>>&,
726  const dimensioned<Type>&
727 );
728 
729 
730 template<class Type>
731 tmp<fvMatrix<Type>> operator==
732 (
733  const fvMatrix<Type>&,
734  const zero&
735 );
736 
737 template<class Type>
738 tmp<fvMatrix<Type>> operator==
739 (
740  const tmp<fvMatrix<Type>>&,
741  const zero&
742 );
743 
744 
745 template<class Type>
746 tmp<fvMatrix<Type>> operator-
747 (
748  const fvMatrix<Type>&
749 );
750 
751 template<class Type>
752 tmp<fvMatrix<Type>> operator-
753 (
754  const tmp<fvMatrix<Type>>&
755 );
756 
757 
758 template<class Type>
759 tmp<fvMatrix<Type>> operator+
760 (
761  const fvMatrix<Type>&,
762  const fvMatrix<Type>&
763 );
764 
765 template<class Type>
766 tmp<fvMatrix<Type>> operator+
767 (
768  const tmp<fvMatrix<Type>>&,
769  const fvMatrix<Type>&
770 );
771 
772 template<class Type>
773 tmp<fvMatrix<Type>> operator+
774 (
775  const fvMatrix<Type>&,
776  const tmp<fvMatrix<Type>>&
777 );
778 
779 template<class Type>
780 tmp<fvMatrix<Type>> operator+
781 (
782  const tmp<fvMatrix<Type>>&,
783  const tmp<fvMatrix<Type>>&
784 );
785 
786 
787 template<class Type>
788 tmp<fvMatrix<Type>> operator+
789 (
790  const fvMatrix<Type>&,
791  const DimensionedField<Type, volMesh>&
792 );
793 
794 template<class Type>
795 tmp<fvMatrix<Type>> operator+
796 (
797  const fvMatrix<Type>&,
798  const tmp<DimensionedField<Type, volMesh>>&
799 );
800 
801 template<class Type>
802 tmp<fvMatrix<Type>> operator+
803 (
804  const fvMatrix<Type>&,
805  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
806 );
807 
808 template<class Type>
809 tmp<fvMatrix<Type>> operator+
810 (
811  const tmp<fvMatrix<Type>>&,
812  const DimensionedField<Type, volMesh>&
813 );
814 
815 template<class Type>
816 tmp<fvMatrix<Type>> operator+
817 (
818  const tmp<fvMatrix<Type>>&,
819  const tmp<DimensionedField<Type, volMesh>>&
820 );
821 
822 template<class Type>
823 tmp<fvMatrix<Type>> operator+
824 (
825  const tmp<fvMatrix<Type>>&,
826  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
827 );
828 
829 template<class Type>
830 tmp<fvMatrix<Type>> operator+
831 (
832  const DimensionedField<Type, volMesh>&,
833  const fvMatrix<Type>&
834 );
835 
836 template<class Type>
837 tmp<fvMatrix<Type>> operator+
838 (
839  const tmp<DimensionedField<Type, volMesh>>&,
840  const fvMatrix<Type>&
841 );
842 
843 template<class Type>
844 tmp<fvMatrix<Type>> operator+
845 (
846  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
847  const fvMatrix<Type>&
848 );
849 
850 template<class Type>
851 tmp<fvMatrix<Type>> operator+
852 (
853  const DimensionedField<Type, volMesh>&,
854  const tmp<fvMatrix<Type>>&
855 );
856 
857 template<class Type>
858 tmp<fvMatrix<Type>> operator+
859 (
860  const tmp<DimensionedField<Type, volMesh>>&,
861  const tmp<fvMatrix<Type>>&
862 );
863 
864 template<class Type>
865 tmp<fvMatrix<Type>> operator+
866 (
867  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
868  const tmp<fvMatrix<Type>>&
869 );
870 
871 
872 template<class Type>
873 tmp<fvMatrix<Type>> operator+
874 (
875  const fvMatrix<Type>&,
876  const dimensioned<Type>&
877 );
878 
879 template<class Type>
880 tmp<fvMatrix<Type>> operator+
881 (
882  const tmp<fvMatrix<Type>>&,
883  const dimensioned<Type>&
884 );
885 
886 template<class Type>
887 tmp<fvMatrix<Type>> operator+
888 (
889  const dimensioned<Type>&,
890  const fvMatrix<Type>&
891 );
892 
893 template<class Type>
894 tmp<fvMatrix<Type>> operator+
895 (
896  const dimensioned<Type>&,
897  const tmp<fvMatrix<Type>>&
898 );
899 
900 
901 template<class Type>
902 tmp<fvMatrix<Type>> operator-
903 (
904  const fvMatrix<Type>&,
905  const fvMatrix<Type>&
906 );
907 
908 template<class Type>
909 tmp<fvMatrix<Type>> operator-
910 (
911  const tmp<fvMatrix<Type>>&,
912  const fvMatrix<Type>&
913 );
914 
915 template<class Type>
916 tmp<fvMatrix<Type>> operator-
917 (
918  const fvMatrix<Type>&,
919  const tmp<fvMatrix<Type>>&
920 );
921 
922 template<class Type>
923 tmp<fvMatrix<Type>> operator-
924 (
925  const tmp<fvMatrix<Type>>&,
926  const tmp<fvMatrix<Type>>&
927 );
928 
929 
930 template<class Type>
931 tmp<fvMatrix<Type>> operator-
932 (
933  const fvMatrix<Type>&,
934  const DimensionedField<Type, volMesh>&
935 );
936 
937 template<class Type>
938 tmp<fvMatrix<Type>> operator-
939 (
940  const fvMatrix<Type>&,
941  const tmp<DimensionedField<Type, volMesh>>&
942 );
943 
944 template<class Type>
945 tmp<fvMatrix<Type>> operator-
946 (
947  const fvMatrix<Type>&,
948  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
949 );
950 
951 template<class Type>
952 tmp<fvMatrix<Type>> operator-
953 (
954  const tmp<fvMatrix<Type>>&,
955  const DimensionedField<Type, volMesh>&
956 );
957 
958 template<class Type>
959 tmp<fvMatrix<Type>> operator-
960 (
961  const tmp<fvMatrix<Type>>&,
962  const tmp<DimensionedField<Type, volMesh>>&
963 );
964 
965 template<class Type>
966 tmp<fvMatrix<Type>> operator-
967 (
968  const tmp<fvMatrix<Type>>&,
969  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
970 );
971 
972 template<class Type>
973 tmp<fvMatrix<Type>> operator-
974 (
975  const DimensionedField<Type, volMesh>&,
976  const fvMatrix<Type>&
977 );
978 
979 template<class Type>
980 tmp<fvMatrix<Type>> operator-
981 (
982  const tmp<DimensionedField<Type, volMesh>>&,
983  const fvMatrix<Type>&
984 );
985 
986 template<class Type>
987 tmp<fvMatrix<Type>> operator-
988 (
989  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
990  const fvMatrix<Type>&
991 );
992 
993 template<class Type>
994 tmp<fvMatrix<Type>> operator-
995 (
996  const DimensionedField<Type, volMesh>&,
997  const tmp<fvMatrix<Type>>&
998 );
999 
1000 template<class Type>
1001 tmp<fvMatrix<Type>> operator-
1002 (
1003  const tmp<DimensionedField<Type, volMesh>>&,
1004  const tmp<fvMatrix<Type>>&
1005 );
1006 
1007 template<class Type>
1008 tmp<fvMatrix<Type>> operator-
1009 (
1010  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
1011  const tmp<fvMatrix<Type>>&
1012 );
1013 
1014 
1015 template<class Type>
1016 tmp<fvMatrix<Type>> operator-
1017 (
1018  const fvMatrix<Type>&,
1019  const dimensioned<Type>&
1020 );
1021 
1022 template<class Type>
1023 tmp<fvMatrix<Type>> operator-
1024 (
1025  const tmp<fvMatrix<Type>>&,
1026  const dimensioned<Type>&
1027 );
1028 
1029 template<class Type>
1030 tmp<fvMatrix<Type>> operator-
1031 (
1032  const dimensioned<Type>&,
1033  const fvMatrix<Type>&
1034 );
1035 
1036 template<class Type>
1037 tmp<fvMatrix<Type>> operator-
1038 (
1039  const dimensioned<Type>&,
1040  const tmp<fvMatrix<Type>>&
1041 );
1042 
1043 
1044 template<class Type>
1045 tmp<fvMatrix<Type>> operator*
1046 (
1047  const volScalarField::Internal&,
1048  const fvMatrix<Type>&
1049 );
1050 
1051 template<class Type>
1052 tmp<fvMatrix<Type>> operator*
1053 (
1054  const tmp<volScalarField::Internal>&,
1055  const fvMatrix<Type>&
1056 );
1057 
1058 template<class Type>
1059 tmp<fvMatrix<Type>> operator*
1060 (
1061  const tmp<volScalarField>&,
1062  const fvMatrix<Type>&
1063 );
1064 
1065 template<class Type>
1066 tmp<fvMatrix<Type>> operator*
1067 (
1068  const volScalarField::Internal&,
1069  const tmp<fvMatrix<Type>>&
1070 );
1071 
1072 template<class Type>
1073 tmp<fvMatrix<Type>> operator*
1074 (
1075  const tmp<volScalarField::Internal>&,
1076  const tmp<fvMatrix<Type>>&
1077 );
1078 
1079 template<class Type>
1080 tmp<fvMatrix<Type>> operator*
1081 (
1082  const tmp<volScalarField>&,
1083  const tmp<fvMatrix<Type>>&
1084 );
1085 
1086 
1087 template<class Type>
1088 tmp<fvMatrix<Type>> operator*
1089 (
1090  const dimensioned<scalar>&,
1091  const fvMatrix<Type>&
1092 );
1093 
1094 template<class Type>
1095 tmp<fvMatrix<Type>> operator*
1096 (
1097  const dimensioned<scalar>&,
1098  const tmp<fvMatrix<Type>>&
1099 );
1100 
1101 
1102 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1103 
1104 } // End namespace Foam
1105 
1106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1107 
1108 #ifdef NoRepository
1109  #include "fvMatrix.C"
1110 #endif
1111 
1112 // Specialisation for scalars
1113 #include "fvScalarMatrix.H"
1114 
1115 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1116 
1117 #endif
1118 
1119 // ************************************************************************* //
Foam::checkMethod
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1035
Foam::fvMatrix::residual
tmp< Field< Type > > residual() const
Return the matrix residual.
Definition: fvMatrixSolve.C:332
volFields.H
Foam::fvMatrix::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:325
Foam::fvMatrix::H
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:799
Foam::fvMatrix::addBoundarySource
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:152
Foam::fvMatrix::boundaryManipulate
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:724
Foam::fvMatrix::addCmptAvBoundaryDiag
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:136
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::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:737
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:81
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:1012
Foam::fvMatrix::fvSolver::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:318
Foam::fvMatrix::dimensions
const dimensionSet & dimensions() const
Definition: fvMatrix.H:292
Foam::fvMatrix::relax
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:707
Foam::fvMatrix::negate
void negate()
Definition: fvMatrix.C:1054
Foam::fvMatrix::internalCoeffs
FieldField< Field, Type > & internalCoeffs()
Definition: fvMatrix.H:316
Foam::fvMatrix::boundaryCoeffs
const FieldField< Field, Type > & boundaryCoeffs() const
Definition: fvMatrix.H:323
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:185
Foam::fvMatrix::addBoundaryDiag
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:118
lduMatrix.H
Foam::fvMatrix::fvSolver::fvSolver
fvSolver(fvMatrix< Type > &fvMat, autoPtr< lduMatrix::solver > &&sol)
Definition: fvMatrix.H:231
Foam::fvMatrix::setReference
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:498
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:287
Foam::fvMatrix::faceFluxCorrectionPtr
surfaceTypeFieldPtr & faceFluxCorrectionPtr()
Return pointer to face-flux non-orthogonal correction field.
Definition: fvMatrix.H:341
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::fvMatrix::solveCoupled
SolverPerformance< Type > solveCoupled(const dictionary &)
Solve coupled returning the solution statistics.
Definition: fvMatrixSolve.C:241
Foam::fvMatrix::operator-=
void operator-=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1103
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::operator*=
void operator*=(const volScalarField::Internal &)
Definition: fvMatrix.C:1238
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:861
className.H
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Foam::fvMatrix::~fvMatrix
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:452
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:514
zero.H
Foam::fvMatrix::ClassName
ClassName("fvMatrix")
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
fvScalarMatrix.H
A scalar instance of fvMatrix.
Foam::fvMatrix::DD
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:746
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:770
fvMatrix.C
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvMatrix::surfaceTypeFieldPtr
GeometricField< Type, fvsPatchField, surfaceMesh > * surfaceTypeFieldPtr
Declare return type of the faceFluxCorrectionPtr() function.
Definition: fvMatrix.H:338
Foam::fvMatrix::fvMatrix
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
Definition: fvMatrix.C:268
Foam::fvMatrix::solver
autoPtr< fvSolver > solver()
Construct and return the solver.
Definition: fvMatrixSolve.C:311
Foam::fvMatrix::setValues
void setValues(const labelUList &cellLabels, const Type &value)
Definition: fvMatrix.C:465
Foam::fvMatrix::boundaryCoeffs
FieldField< Field, Type > & boundaryCoeffs()
Definition: fvMatrix.H:330
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
tmp.H
Foam::fvMatrix::source
const Field< Type > & source() const
Definition: fvMatrix.H:302
Foam::fvMatrix::operator+=
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1069
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:297
Foam::fvMatrix::solverDict
const dictionary & solverDict() const
Return the solver dictionary taking into account finalIteration.
Definition: fvMatrix.C:996
Foam::UList< label >
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::fvMatrix::clone
tmp< fvMatrix< Type > > clone() const
Clone.
Definition: fvMatrix.C:440
dimensionedTypes.H
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fvMatrix::fvSolver
Definition: fvMatrix.H:221
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::fvMatrix::solveSegregated
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
Definition: fvMatrixSolve.C:115
Foam::fvMatrix::internalCoeffs
const FieldField< Field, Type > & internalCoeffs() const
Definition: fvMatrix.H:309
Foam::GeometricField< Type, fvPatchField, volMesh >
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:43
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
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:909
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