fvMatrix.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::fvMatrix
26 
27 Description
28  A special matrix type and solver, designed for finite volume
29  solutions of scalar equations.
30  Face addressing is used to make all matrix assembly
31  and solution loops vectorise.
32 
33 SourceFiles
34  fvMatrix.C
35  fvMatrixSolve.C
36  fvScalarMatrix.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef fvMatrix_H
41 #define fvMatrix_H
42 
43 #include "volFields.H"
44 #include "surfaceFields.H"
45 #include "lduMatrix.H"
46 #include "tmp.H"
47 #include "autoPtr.H"
48 #include "dimensionedTypes.H"
49 #include "zero.H"
50 #include "className.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 // Forward declaration of friend functions and operators
58 
59 template<class Type>
60 class fvMatrix;
61 
62 template<class Type>
63 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
64 (
65  const fvMatrix<Type>&,
66  const DimensionedField<Type, volMesh>&
67 );
68 
69 template<class Type>
70 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
71 (
72  const fvMatrix<Type>&,
73  const tmp<DimensionedField<Type, volMesh>>&
74 );
75 
76 template<class Type>
77 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
78 (
79  const fvMatrix<Type>&,
80  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
81 );
82 
83 template<class Type>
84 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
85 (
86  const tmp<fvMatrix<Type>>&,
87  const DimensionedField<Type, volMesh>&
88 );
89 
90 template<class Type>
91 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
92 (
93  const tmp<fvMatrix<Type>>&,
94  const tmp<DimensionedField<Type, volMesh>>&
95 );
96 
97 template<class Type>
98 tmp<GeometricField<Type, fvPatchField, volMesh>> operator&
99 (
100  const tmp<fvMatrix<Type>>&,
101  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
102 );
103 
104 template<class Type>
105 Ostream& operator<<(Ostream&, const fvMatrix<Type>&);
107 template<class T> class UIndirectList;
108 
109 
110 /*---------------------------------------------------------------------------*\
111  Class fvMatrix Declaration
112 \*---------------------------------------------------------------------------*/
113 
114 template<class Type>
115 class fvMatrix
116 :
117  public tmp<fvMatrix<Type>>::refCount,
118  public lduMatrix
119 {
120  // Private Data
121 
122  //- Const reference to GeometricField<Type, fvPatchField, volMesh>
123  // Converted into a non-const reference at the point of solution.
125 
126  //- Dimension set
127  dimensionSet dimensions_;
128 
129  //- Source term
130  Field<Type> source_;
131 
132  //- Boundary scalar field containing pseudo-matrix coeffs
133  // for internal cells
134  FieldField<Field, Type> internalCoeffs_;
135 
136  //- Boundary scalar field containing pseudo-matrix coeffs
137  // for boundary cells
138  FieldField<Field, Type> boundaryCoeffs_;
139 
140 
141  //- Face flux field for non-orthogonal correction
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  // Matrix manipulation functionality
205 
206  //- Set solution in the given cell to the specified value
207  void setValue
208  (
209  const label celli,
210  const Type& value
211  );
212 
213  //- Set solution in the given cell to the specified value
214  void setValue
215  (
216  const label celli,
217  const Type& value,
218  const scalar fraction,
219  const scalarField& ddtDiag
220  );
221 
222 
223 public:
224 
225  //- Solver class returned by the solver function
226  // used for systems in which it is useful to cache the solver for reuse
227  // e.g. if the solver is potentially expensive to construct (AMG) and can
228  // be used several times (PISO)
229  class fvSolver
230  {
231  fvMatrix<Type>& fvMat_;
232 
234 
235  public:
236 
237  // Constructors
239  fvSolver(fvMatrix<Type>& fvMat, autoPtr<lduMatrix::solver> sol)
240  :
241  fvMat_(fvMat),
242  solver_(sol)
243  {}
244 
245 
246  // Member Functions
247 
248  //- Solve returning the solution statistics.
249  // Use the given solver controls
251 
252  //- Solve returning the solution statistics.
253  // Solver controls read from fvSolution
255  };
256 
257 
258  ClassName("fvMatrix");
259 
260 
261  // Constructors
262 
263  //- Construct given a field to solve for
264  fvMatrix
265  (
267  const dimensionSet&
268  );
269 
270  //- Copy constructor
271  fvMatrix(const fvMatrix<Type>&);
272 
273  //- Copy constructor of tmp<fvMatrix<Type>> deleting argument
274  fvMatrix(const tmp<fvMatrix<Type>>&);
275 
276  //- Construct from Istream given field to solve for
278 
279  //- Clone
280  tmp<fvMatrix<Type>> clone() const;
281 
282 
283  //- Destructor
284  virtual ~fvMatrix();
285 
286 
287  // Member Functions
288 
289  // Access
292  {
293  return psi_;
294  }
296  const dimensionSet& dimensions() const
297  {
298  return dimensions_;
299  }
302  {
303  return source_;
304  }
306  const Field<Type>& source() const
307  {
308  return source_;
309  }
310 
311  //- fvBoundary scalar field containing pseudo-matrix coeffs
312  // for internal cells
314  {
315  return internalCoeffs_;
316  }
317 
318  //- fvBoundary scalar field containing pseudo-matrix coeffs
319  // for boundary cells
321  {
322  return boundaryCoeffs_;
323  }
324 
325 
326  //- Declare return type of the faceFluxCorrectionPtr() function
329 
330  //- Return pointer to face-flux non-orthogonal correction field
331  surfaceTypeFieldPtr& faceFluxCorrectionPtr()
332  {
333  return faceFluxCorrectionPtr_;
334  }
335 
336 
337  // Operations
338 
339  //- Set solution in given cells to the specified values
340  // and eliminate the corresponding equations from the matrix.
341  template<template<class> class ListType>
342  void setValues
343  (
344  const labelUList& cells,
345  const ListType<Type>& values
346  );
347 
348  //- Set solution in given cells to the specified values
349  // and eliminate the corresponding equations from the matrix.
350  template<template<class> class ListType>
351  void setValues
352  (
353  const labelUList& cells,
354  const ListType<Type>& values,
355  const scalarList& fractions,
356  const bool hasDdt = true
357  );
358 
359  //- Set reference level for solution
360  void setReference
361  (
362  const label celli,
363  const Type& value,
364  const bool forceReference = false
365  );
366 
367  //- Set reference level for a component of the solution
368  // on a given patch face
370  (
371  const label patchi,
372  const label facei,
373  const direction cmpt,
374  const scalar value
375  );
376 
377  //- Return the relaxation factor for this equation
378  scalar relaxationFactor() const;
379 
380  //- Relax matrix (for steady-state solution).
381  // alpha = 1 : diagonally equal
382  // alpha < 1 : diagonally dominant
383  // alpha = 0 : do nothing
384  // Note: Requires positive diagonal.
385  void relax(const scalar alpha);
386 
387  //- Relax matrix (for steady-state solution).
388  // alpha is read from controlDict
389  void relax();
390 
391  //- Manipulate based on a boundary field
392  void boundaryManipulate
393  (
395  Boundary& values
396  );
397 
398  //- Construct and return the solver
399  // Use the given solver controls
401 
402  //- Construct and return the solver
403  // Solver controls read from fvSolution
405 
406  //- Solve segregated or coupled returning the solution statistics.
407  // Use the given solver controls
409 
410  //- Solve segregated returning the solution statistics.
411  // Use the given solver controls
413 
414  //- Solve coupled returning the solution statistics.
415  // Use the given solver controls
417 
418  //- Solve segregated or coupled returning the solution statistics.
419  // Solver controls read from fvSolution
421 
422  //- Solve returning the solution statistics.
423  // Solver controls read from fvSolution
425 
426  //- Return the matrix residual
427  tmp<Field<Type>> residual() const;
428 
429  //- Return the matrix scalar diagonal
430  tmp<scalarField> D() const;
431 
432  //- Return the matrix Type diagonal
433  tmp<Field<Type>> DD() const;
434 
435  //- Return the central coefficient
436  tmp<volScalarField> A() const;
437 
438  //- Return the H operation source
440 
441  //- Return H(1)
442  tmp<volScalarField> H1() const;
443 
444  //- Return the face-flux field from the matrix
446  flux() const;
447 
448 
449  // Member Operators
450 
451  void operator=(const fvMatrix<Type>&);
452  void operator=(const tmp<fvMatrix<Type>>&);
453 
454  void negate();
455 
456  void operator+=(const fvMatrix<Type>&);
457  void operator+=(const tmp<fvMatrix<Type>>&);
458 
459  void operator-=(const fvMatrix<Type>&);
460  void operator-=(const tmp<fvMatrix<Type>>&);
461 
462  void operator+=
463  (
465  );
466  void operator+=
467  (
469  );
470  void operator+=
471  (
473  );
474 
475  void operator-=
476  (
478  );
479  void operator-=
480  (
482  );
483  void operator-=
484  (
486  );
487 
488  void operator+=(const dimensioned<Type>&);
489  void operator-=(const dimensioned<Type>&);
490 
491  void operator+=(const zero&);
492  void operator-=(const zero&);
493 
496  void operator*=(const tmp<volScalarField>&);
497  void operator*=(const dimensioned<scalar>&);
498 
501  void operator/=(const tmp<volScalarField>&);
502  void operator/=(const dimensioned<scalar>&);
503 
504 
505  // Friend operators
506 
508  operator& <Type>
509  (
510  const fvMatrix<Type>&,
512  );
513 
515  operator& <Type>
516  (
517  const fvMatrix<Type>&,
519  );
520 
522  operator& <Type>
523  (
524  const tmp<fvMatrix<Type>>&,
526  );
527 
529  operator& <Type>
530  (
531  const tmp<fvMatrix<Type>>&,
533  );
534 
535 
536  // Ostream operator
537 
538  friend Ostream& operator<< <Type>
539  (
540  Ostream&,
541  const fvMatrix<Type>&
542  );
543 };
544 
545 
546 // * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * //
547 
548 template<class Type>
549 void checkMethod
550 (
551  const fvMatrix<Type>&,
552  const fvMatrix<Type>&,
553  const char*
554 );
555 
556 template<class Type>
557 void checkMethod
558 (
559  const fvMatrix<Type>&,
561  const char*
562 );
563 
564 template<class Type>
565 void checkMethod
566 (
567  const fvMatrix<Type>&,
568  const dimensioned<Type>&,
569  const char*
570 );
571 
572 
573 //- Solve returning the solution statistics given convergence tolerance
574 // Use the given solver controls
575 template<class Type>
576 SolverPerformance<Type> solve(fvMatrix<Type>&, const word&);
577 
578 //- Solve returning the solution statistics given convergence tolerance,
579 // deleting temporary matrix after solution.
580 // Use the given solver controls
581 template<class Type>
583 (
584  const tmp<fvMatrix<Type>>&,
585  const word&
586 );
587 
588 //- Solve returning the solution statistics given convergence tolerance
589 // Solver controls read fvSolution
590 template<class Type>
591 SolverPerformance<Type> solve(fvMatrix<Type>&);
592 
593 //- Solve returning the solution statistics given convergence tolerance,
594 // deleting temporary matrix after solution.
595 // Solver controls read fvSolution
596 template<class Type>
597 SolverPerformance<Type> solve(const tmp<fvMatrix<Type>>&);
598 
599 //- Return the correction form of the given matrix
600 // by subtracting the matrix multiplied by the current field
601 template<class Type>
602 tmp<fvMatrix<Type>> correction(const fvMatrix<Type>&);
603 
604 //- Return the correction form of the given temporary matrix
605 // by subtracting the matrix multiplied by the current field
606 template<class Type>
607 tmp<fvMatrix<Type>> correction(const tmp<fvMatrix<Type>>&);
608 
609 
610 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
611 
612 template<class Type>
613 tmp<fvMatrix<Type>> operator==
614 (
615  const fvMatrix<Type>&,
616  const fvMatrix<Type>&
617 );
618 
619 template<class Type>
620 tmp<fvMatrix<Type>> operator==
621 (
622  const tmp<fvMatrix<Type>>&,
623  const fvMatrix<Type>&
624 );
625 
626 template<class Type>
627 tmp<fvMatrix<Type>> operator==
628 (
629  const fvMatrix<Type>&,
630  const tmp<fvMatrix<Type>>&
631 );
632 
633 template<class Type>
634 tmp<fvMatrix<Type>> operator==
635 (
636  const tmp<fvMatrix<Type>>&,
637  const tmp<fvMatrix<Type>>&
638 );
639 
640 
641 template<class Type>
642 tmp<fvMatrix<Type>> operator==
643 (
644  const fvMatrix<Type>&,
646 );
647 
648 template<class Type>
649 tmp<fvMatrix<Type>> operator==
650 (
651  const fvMatrix<Type>&,
653 );
654 
655 template<class Type>
656 tmp<fvMatrix<Type>> operator==
657 (
658  const fvMatrix<Type>&,
660 );
661 
662 template<class Type>
663 tmp<fvMatrix<Type>> operator==
664 (
665  const tmp<fvMatrix<Type>>&,
667 );
668 
669 template<class Type>
670 tmp<fvMatrix<Type>> operator==
671 (
672  const tmp<fvMatrix<Type>>&,
674 );
675 
676 template<class Type>
677 tmp<fvMatrix<Type>> operator==
678 (
679  const tmp<fvMatrix<Type>>&,
681 );
682 
683 template<class Type>
684 tmp<fvMatrix<Type>> operator==
685 (
686  const fvMatrix<Type>&,
687  const dimensioned<Type>&
688 );
689 
690 template<class Type>
691 tmp<fvMatrix<Type>> operator==
692 (
693  const tmp<fvMatrix<Type>>&,
694  const dimensioned<Type>&
695 );
696 
697 
698 template<class Type>
699 tmp<fvMatrix<Type>> operator==
700 (
701  const fvMatrix<Type>&,
702  const zero&
703 );
704 
705 template<class Type>
706 tmp<fvMatrix<Type>> operator==
707 (
708  const tmp<fvMatrix<Type>>&,
709  const zero&
710 );
711 
712 
713 template<class Type>
714 tmp<fvMatrix<Type>> operator-
715 (
716  const fvMatrix<Type>&
717 );
718 
719 template<class Type>
720 tmp<fvMatrix<Type>> operator-
721 (
722  const tmp<fvMatrix<Type>>&
723 );
724 
725 
726 template<class Type>
727 tmp<fvMatrix<Type>> operator+
728 (
729  const fvMatrix<Type>&,
730  const fvMatrix<Type>&
731 );
732 
733 template<class Type>
734 tmp<fvMatrix<Type>> operator+
735 (
736  const tmp<fvMatrix<Type>>&,
737  const fvMatrix<Type>&
738 );
739 
740 template<class Type>
741 tmp<fvMatrix<Type>> operator+
742 (
743  const fvMatrix<Type>&,
744  const tmp<fvMatrix<Type>>&
745 );
746 
747 template<class Type>
748 tmp<fvMatrix<Type>> operator+
749 (
750  const tmp<fvMatrix<Type>>&,
751  const tmp<fvMatrix<Type>>&
752 );
753 
754 
755 template<class Type>
756 tmp<fvMatrix<Type>> operator+
757 (
758  const fvMatrix<Type>&,
760 );
761 
762 template<class Type>
763 tmp<fvMatrix<Type>> operator+
764 (
765  const fvMatrix<Type>&,
767 );
768 
769 template<class Type>
770 tmp<fvMatrix<Type>> operator+
771 (
772  const fvMatrix<Type>&,
774 );
775 
776 template<class Type>
777 tmp<fvMatrix<Type>> operator+
778 (
779  const tmp<fvMatrix<Type>>&,
781 );
782 
783 template<class Type>
784 tmp<fvMatrix<Type>> operator+
785 (
786  const tmp<fvMatrix<Type>>&,
788 );
789 
790 template<class Type>
791 tmp<fvMatrix<Type>> operator+
792 (
793  const tmp<fvMatrix<Type>>&,
795 );
796 
797 template<class Type>
798 tmp<fvMatrix<Type>> operator+
799 (
801  const fvMatrix<Type>&
802 );
803 
804 template<class Type>
805 tmp<fvMatrix<Type>> operator+
806 (
808  const fvMatrix<Type>&
809 );
810 
811 template<class Type>
812 tmp<fvMatrix<Type>> operator+
813 (
815  const fvMatrix<Type>&
816 );
817 
818 template<class Type>
819 tmp<fvMatrix<Type>> operator+
820 (
821  const DimensionedField<Type, volMesh>&,
822  const tmp<fvMatrix<Type>>&
823 );
824 
825 template<class Type>
826 tmp<fvMatrix<Type>> operator+
827 (
828  const tmp<DimensionedField<Type, volMesh>>&,
829  const tmp<fvMatrix<Type>>&
830 );
831 
832 template<class Type>
833 tmp<fvMatrix<Type>> operator+
834 (
835  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
836  const tmp<fvMatrix<Type>>&
837 );
838 
839 
840 template<class Type>
841 tmp<fvMatrix<Type>> operator+
842 (
843  const fvMatrix<Type>&,
844  const dimensioned<Type>&
845 );
846 
847 template<class Type>
848 tmp<fvMatrix<Type>> operator+
849 (
850  const tmp<fvMatrix<Type>>&,
851  const dimensioned<Type>&
852 );
853 
854 template<class Type>
855 tmp<fvMatrix<Type>> operator+
856 (
857  const dimensioned<Type>&,
858  const fvMatrix<Type>&
859 );
860 
861 template<class Type>
862 tmp<fvMatrix<Type>> operator+
863 (
864  const dimensioned<Type>&,
865  const tmp<fvMatrix<Type>>&
866 );
867 
868 
869 template<class Type>
870 tmp<fvMatrix<Type>> operator-
871 (
872  const fvMatrix<Type>&,
873  const fvMatrix<Type>&
874 );
875 
876 template<class Type>
877 tmp<fvMatrix<Type>> operator-
878 (
879  const tmp<fvMatrix<Type>>&,
880  const fvMatrix<Type>&
881 );
882 
883 template<class Type>
884 tmp<fvMatrix<Type>> operator-
885 (
886  const fvMatrix<Type>&,
887  const tmp<fvMatrix<Type>>&
888 );
889 
890 template<class Type>
891 tmp<fvMatrix<Type>> operator-
892 (
893  const tmp<fvMatrix<Type>>&,
894  const tmp<fvMatrix<Type>>&
895 );
896 
897 
898 template<class Type>
899 tmp<fvMatrix<Type>> operator-
900 (
901  const fvMatrix<Type>&,
902  const DimensionedField<Type, volMesh>&
903 );
904 
905 template<class Type>
906 tmp<fvMatrix<Type>> operator-
907 (
908  const fvMatrix<Type>&,
909  const tmp<DimensionedField<Type, volMesh>>&
910 );
911 
912 template<class Type>
913 tmp<fvMatrix<Type>> operator-
914 (
915  const fvMatrix<Type>&,
916  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
917 );
918 
919 template<class Type>
920 tmp<fvMatrix<Type>> operator-
921 (
922  const tmp<fvMatrix<Type>>&,
923  const DimensionedField<Type, volMesh>&
924 );
925 
926 template<class Type>
927 tmp<fvMatrix<Type>> operator-
928 (
929  const tmp<fvMatrix<Type>>&,
930  const tmp<DimensionedField<Type, volMesh>>&
931 );
932 
933 template<class Type>
934 tmp<fvMatrix<Type>> operator-
935 (
936  const tmp<fvMatrix<Type>>&,
937  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
938 );
939 
940 template<class Type>
941 tmp<fvMatrix<Type>> operator-
942 (
943  const DimensionedField<Type, volMesh>&,
944  const fvMatrix<Type>&
945 );
946 
947 template<class Type>
948 tmp<fvMatrix<Type>> operator-
949 (
950  const tmp<DimensionedField<Type, volMesh>>&,
951  const fvMatrix<Type>&
952 );
953 
954 template<class Type>
955 tmp<fvMatrix<Type>> operator-
956 (
957  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
958  const fvMatrix<Type>&
959 );
960 
961 template<class Type>
962 tmp<fvMatrix<Type>> operator-
963 (
964  const DimensionedField<Type, volMesh>&,
965  const tmp<fvMatrix<Type>>&
966 );
967 
968 template<class Type>
969 tmp<fvMatrix<Type>> operator-
970 (
971  const tmp<DimensionedField<Type, volMesh>>&,
972  const tmp<fvMatrix<Type>>&
973 );
974 
975 template<class Type>
976 tmp<fvMatrix<Type>> operator-
977 (
978  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
979  const tmp<fvMatrix<Type>>&
980 );
981 
982 
983 template<class Type>
984 tmp<fvMatrix<Type>> operator-
985 (
986  const fvMatrix<Type>&,
987  const dimensioned<Type>&
988 );
989 
990 template<class Type>
991 tmp<fvMatrix<Type>> operator-
992 (
993  const tmp<fvMatrix<Type>>&,
994  const dimensioned<Type>&
995 );
996 
997 template<class Type>
998 tmp<fvMatrix<Type>> operator-
999 (
1000  const dimensioned<Type>&,
1001  const fvMatrix<Type>&
1002 );
1003 
1004 template<class Type>
1005 tmp<fvMatrix<Type>> operator-
1006 (
1007  const dimensioned<Type>&,
1008  const tmp<fvMatrix<Type>>&
1009 );
1010 
1011 
1012 template<class Type>
1013 tmp<fvMatrix<Type>> operator*
1014 (
1015  const volScalarField::Internal&,
1016  const fvMatrix<Type>&
1017 );
1018 
1019 template<class Type>
1020 tmp<fvMatrix<Type>> operator*
1021 (
1023  const fvMatrix<Type>&
1024 );
1025 
1026 template<class Type>
1027 tmp<fvMatrix<Type>> operator*
1028 (
1029  const tmp<volScalarField>&,
1030  const fvMatrix<Type>&
1031 );
1032 
1033 template<class Type>
1034 tmp<fvMatrix<Type>> operator*
1035 (
1036  const volScalarField::Internal&,
1037  const tmp<fvMatrix<Type>>&
1038 );
1039 
1040 template<class Type>
1041 tmp<fvMatrix<Type>> operator*
1042 (
1043  const tmp<volScalarField::Internal>&,
1044  const tmp<fvMatrix<Type>>&
1045 );
1046 
1047 template<class Type>
1048 tmp<fvMatrix<Type>> operator*
1049 (
1050  const tmp<volScalarField>&,
1051  const tmp<fvMatrix<Type>>&
1052 );
1053 
1054 
1055 template<class Type>
1056 tmp<fvMatrix<Type>> operator*
1057 (
1058  const dimensioned<scalar>&,
1059  const fvMatrix<Type>&
1060 );
1061 
1062 template<class Type>
1063 tmp<fvMatrix<Type>> operator*
1064 (
1065  const dimensioned<scalar>&,
1066  const tmp<fvMatrix<Type>>&
1067 );
1068 
1069 template<class Type>
1070 tmp<fvMatrix<Type>> operator*
1071 (
1072  const volScalarField::Internal&,
1073  const fvMatrix<Type>&
1074 );
1075 
1076 template<class Type>
1077 tmp<fvMatrix<Type>> operator*
1078 (
1079  const tmp<volScalarField::Internal>&,
1080  const fvMatrix<Type>&
1081 );
1082 
1083 template<class Type>
1084 tmp<fvMatrix<Type>> operator*
1085 (
1086  const tmp<volScalarField>&,
1087  const fvMatrix<Type>&
1088 );
1089 
1090 template<class Type>
1091 tmp<fvMatrix<Type>> operator*
1092 (
1093  const volScalarField::Internal&,
1094  const tmp<fvMatrix<Type>>&
1095 );
1096 
1097 template<class Type>
1098 tmp<fvMatrix<Type>> operator*
1099 (
1100  const tmp<volScalarField::Internal>&,
1101  const tmp<fvMatrix<Type>>&
1102 );
1103 
1104 template<class Type>
1105 tmp<fvMatrix<Type>> operator*
1106 (
1107  const tmp<volScalarField>&,
1108  const tmp<fvMatrix<Type>>&
1109 );
1110 
1111 
1112 template<class Type>
1113 tmp<fvMatrix<Type>> operator*
1114 (
1115  const dimensioned<scalar>&,
1116  const fvMatrix<Type>&
1117 );
1118 
1119 template<class Type>
1120 tmp<fvMatrix<Type>> operator*
1121 (
1122  const dimensioned<scalar>&,
1123  const tmp<fvMatrix<Type>>&
1124 );
1125 
1126 
1127 template<class Type>
1128 tmp<fvMatrix<Type>> operator/
1129 (
1130  const fvMatrix<Type>&,
1131  const volScalarField::Internal&
1132 );
1133 
1134 template<class Type>
1135 tmp<fvMatrix<Type>> operator/
1136 (
1137  const fvMatrix<Type>&,
1138  const tmp<volScalarField::Internal>&
1139 );
1140 
1141 template<class Type>
1142 tmp<fvMatrix<Type>> operator/
1143 (
1144  const fvMatrix<Type>&,
1145  const tmp<volScalarField>&
1146 );
1147 
1148 template<class Type>
1149 tmp<fvMatrix<Type>> operator/
1150 (
1151  const tmp<fvMatrix<Type>>&,
1152  const volScalarField::Internal&
1153 );
1154 
1155 template<class Type>
1156 tmp<fvMatrix<Type>> operator/
1157 (
1158  const tmp<fvMatrix<Type>>&,
1159  const tmp<volScalarField::Internal>&
1160 );
1161 
1162 template<class Type>
1163 tmp<fvMatrix<Type>> operator/
1164 (
1165  const tmp<fvMatrix<Type>>&,
1166  const tmp<volScalarField>&
1167 );
1168 
1169 
1170 template<class Type>
1171 tmp<fvMatrix<Type>> operator/
1172 (
1173  const fvMatrix<Type>&,
1174  const dimensioned<scalar>&
1175 );
1176 
1177 template<class Type>
1178 tmp<fvMatrix<Type>> operator/
1179 (
1180  const tmp<fvMatrix<Type>>&,
1181  const dimensioned<scalar>&
1182 );
1183 
1184 
1185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1186 
1187 } // End namespace Foam
1188 
1189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1190 
1191 #ifdef NoRepository
1192  #include "fvMatrix.C"
1193 #endif
1194 
1195 // Specialisation for scalars
1196 #include "fvScalarMatrix.H"
1197 
1198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1199 
1200 #endif
1201 
1202 // ************************************************************************* //
Foam::surfaceFields.
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: fvMatrix.C:40
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
void setValue(const label celli, const Type &value)
Set solution in the given cell to the specified value.
Definition: fvMatrix.C:179
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:113
Reference counter for various OpenFOAM components.
Definition: refCount.H:49
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:755
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void negate()
Definition: fvMatrix.C:1051
FieldField< Field, Type > & internalCoeffs()
fvBoundary scalar field containing pseudo-matrix coeffs
Definition: fvMatrix.H:312
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:290
uint8_t direction
Definition: direction.H:45
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:131
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:563
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:831
void operator-=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1100
scalar relaxationFactor() const
Return the relaxation factor for this equation.
Definition: fvMatrix.C:578
Solver class returned by the solver function.
Definition: fvMatrix.H:228
Generic GeometricField class.
Generic dimensioned Type class.
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:886
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: fvMatrix.C:77
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Generic field type.
Definition: FieldField.H:51
tmp< scalarField > D() const
Return the matrix scalar diagonal.
Definition: fvMatrix.C:776
void operator=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1007
SolverPerformance< Type > solveCoupled(const dictionary &)
Solve coupled returning the solution statistics.
surfaceTypeFieldPtr & faceFluxCorrectionPtr()
Return pointer to face-flux non-orthogonal correction field.
Definition: fvMatrix.H:330
fvSolver(fvMatrix< Type > &fvMat, autoPtr< lduMatrix::solver > sol)
Definition: fvMatrix.H:238
Dimension set for the base types.
Definition: dimensionSet.H:121
SolverPerformance< Type > solve()
Solve returning the solution statistics.
const cellShapeList & cells
void operator/=(const volScalarField::Internal &)
Definition: fvMatrix.C:1306
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
FieldField< Field, Type > & boundaryCoeffs()
fvBoundary scalar field containing pseudo-matrix coeffs
Definition: fvMatrix.H:319
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:482
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
tmp< fvMatrix< Type > > clone() const
Clone.
Definition: fvMatrix.C:470
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Field< Type > & source()
Definition: fvMatrix.H:300
void operator*=(const volScalarField::Internal &)
Definition: fvMatrix.C:1235
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1066
void setValues(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:502
ClassName("fvMatrix")
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:785
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
GeometricField< Type, fvsPatchField, surfaceMesh > * surfaceTypeFieldPtr
Declare return type of the faceFluxCorrectionPtr() function.
Definition: fvMatrix.H:327
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:809
label patchi
A scalar instance of fvMatrix.
Macro definitions for declaring ClassName(), NamespaceName(), etc.
A List with indirect addressing.
Definition: fvMatrix.H:106
autoPtr< fvSolver > solver()
Construct and return the solver.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void checkMethod(const fvMatrix< Type > &, const fvMatrix< Type > &, const char *)
Definition: fvMatrix.C:1379
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:49
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:147
dictionary fractions(initialConditions.subDict("fractions"))
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
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:34
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:927
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:763
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField & diag()
Definition: lduMatrix.C:186
tmp< Field< Type > > residual() const
Return the matrix residual.
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &, const dimensionSet &)
Construct given a field to solve for.
Definition: fvMatrix.C:285
Namespace for OpenFOAM.
const dimensionSet & dimensions() const
Definition: fvMatrix.H:295