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