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-2020 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 given cells to the specified values
207  template<template<class> class ListType>
208  void setValuesFromList
209  (
210  const labelUList& cells,
211  const ListType<Type>& values
212  );
213 
214 
215 public:
216 
217  //- Solver class returned by the solver function
218  // used for systems in which it is useful to cache the solver for reuse
219  // e.g. if the solver is potentially expensive to construct (AMG) and can
220  // be used several times (PISO)
221  class fvSolver
222  {
223  fvMatrix<Type>& fvMat_;
224 
226 
227  public:
228 
229  // Constructors
231  fvSolver(fvMatrix<Type>& fvMat, autoPtr<lduMatrix::solver> sol)
232  :
233  fvMat_(fvMat),
234  solver_(sol)
235  {}
236 
237 
238  // Member Functions
239 
240  //- Solve returning the solution statistics.
241  // Use the given solver controls
243 
244  //- Solve returning the solution statistics.
245  // Solver controls read from fvSolution
247  };
248 
249 
250  ClassName("fvMatrix");
251 
252 
253  // Constructors
254 
255  //- Construct given a field to solve for
256  fvMatrix
257  (
259  const dimensionSet&
260  );
261 
262  //- Copy constructor
263  fvMatrix(const fvMatrix<Type>&);
264 
265  //- Copy constructor of tmp<fvMatrix<Type>> deleting argument
266  fvMatrix(const tmp<fvMatrix<Type>>&);
267 
268  //- Construct from Istream given field to solve for
270 
271  //- Clone
272  tmp<fvMatrix<Type>> clone() const;
273 
274 
275  //- Destructor
276  virtual ~fvMatrix();
277 
278 
279  // Member Functions
280 
281  // Access
284  {
285  return psi_;
286  }
288  const dimensionSet& dimensions() const
289  {
290  return dimensions_;
291  }
294  {
295  return source_;
296  }
298  const Field<Type>& source() const
299  {
300  return source_;
301  }
302 
303  //- fvBoundary scalar field containing pseudo-matrix coeffs
304  // for internal cells
306  {
307  return internalCoeffs_;
308  }
309 
310  //- fvBoundary scalar field containing pseudo-matrix coeffs
311  // for boundary cells
313  {
314  return boundaryCoeffs_;
315  }
316 
317 
318  //- Declare return type of the faceFluxCorrectionPtr() function
321 
322  //- Return pointer to face-flux non-orthogonal correction field
323  surfaceTypeFieldPtr& faceFluxCorrectionPtr()
324  {
325  return faceFluxCorrectionPtr_;
326  }
327 
328 
329  // Operations
330 
331  //- Set solution in given cells to the specified values
332  // and eliminate the corresponding equations from the matrix.
333  void setValues
334  (
335  const labelUList& cells,
336  const UList<Type>& values
337  );
338 
339  //- Set solution in given cells to the specified values
340  // and eliminate the corresponding equations from the matrix.
341  void setValues
342  (
343  const labelUList& cells,
344  const UIndirectList<Type>& values
345  );
346 
347  //- Set reference level for solution
348  void setReference
349  (
350  const label celli,
351  const Type& value,
352  const bool forceReference = false
353  );
354 
355  //- Set reference level for a component of the solution
356  // on a given patch face
358  (
359  const label patchi,
360  const label facei,
361  const direction cmpt,
362  const scalar value
363  );
364 
365  //- Relax matrix (for steady-state solution).
366  // alpha = 1 : diagonally equal
367  // alpha < 1 : diagonally dominant
368  // alpha = 0 : do nothing
369  // Note: Requires positive diagonal.
370  void relax(const scalar alpha);
371 
372  //- Relax matrix (for steady-state solution).
373  // alpha is read from controlDict
374  void relax();
375 
376  //- Manipulate based on a boundary field
377  void boundaryManipulate
378  (
380  Boundary& values
381  );
382 
383  //- Construct and return the solver
384  // Use the given solver controls
386 
387  //- Construct and return the solver
388  // Solver controls read from fvSolution
390 
391  //- Solve segregated or coupled returning the solution statistics.
392  // Use the given solver controls
394 
395  //- Solve segregated returning the solution statistics.
396  // Use the given solver controls
398 
399  //- Solve coupled returning the solution statistics.
400  // Use the given solver controls
402 
403  //- Solve segregated or coupled returning the solution statistics.
404  // Solver controls read from fvSolution
406 
407  //- Solve returning the solution statistics.
408  // Solver controls read from fvSolution
410 
411  //- Return the matrix residual
412  tmp<Field<Type>> residual() const;
413 
414  //- Return the matrix scalar diagonal
415  tmp<scalarField> D() const;
416 
417  //- Return the matrix Type diagonal
418  tmp<Field<Type>> DD() const;
419 
420  //- Return the central coefficient
421  tmp<volScalarField> A() const;
422 
423  //- Return the H operation source
425 
426  //- Return H(1)
427  tmp<volScalarField> H1() const;
428 
429  //- Return the face-flux field from the matrix
431  flux() const;
432 
433 
434  // Member Operators
435 
436  void operator=(const fvMatrix<Type>&);
437  void operator=(const tmp<fvMatrix<Type>>&);
438 
439  void negate();
440 
441  void operator+=(const fvMatrix<Type>&);
442  void operator+=(const tmp<fvMatrix<Type>>&);
443 
444  void operator-=(const fvMatrix<Type>&);
445  void operator-=(const tmp<fvMatrix<Type>>&);
446 
447  void operator+=
448  (
450  );
451  void operator+=
452  (
454  );
455  void operator+=
456  (
458  );
459 
460  void operator-=
461  (
463  );
464  void operator-=
465  (
467  );
468  void operator-=
469  (
471  );
472 
473  void operator+=(const dimensioned<Type>&);
474  void operator-=(const dimensioned<Type>&);
475 
476  void operator+=(const zero&);
477  void operator-=(const zero&);
478 
481  void operator*=(const tmp<volScalarField>&);
482  void operator*=(const dimensioned<scalar>&);
483 
486  void operator/=(const tmp<volScalarField>&);
487  void operator/=(const dimensioned<scalar>&);
488 
489 
490  // Friend operators
491 
493  operator& <Type>
494  (
495  const fvMatrix<Type>&,
497  );
498 
500  operator& <Type>
501  (
502  const fvMatrix<Type>&,
504  );
505 
507  operator& <Type>
508  (
509  const tmp<fvMatrix<Type>>&,
511  );
512 
514  operator& <Type>
515  (
516  const tmp<fvMatrix<Type>>&,
518  );
519 
520 
521  // Ostream operator
522 
523  friend Ostream& operator<< <Type>
524  (
525  Ostream&,
526  const fvMatrix<Type>&
527  );
528 };
529 
530 
531 // * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * //
532 
533 template<class Type>
534 void checkMethod
535 (
536  const fvMatrix<Type>&,
537  const fvMatrix<Type>&,
538  const char*
539 );
540 
541 template<class Type>
542 void checkMethod
543 (
544  const fvMatrix<Type>&,
546  const char*
547 );
548 
549 template<class Type>
550 void checkMethod
551 (
552  const fvMatrix<Type>&,
553  const dimensioned<Type>&,
554  const char*
555 );
556 
557 
558 //- Solve returning the solution statistics given convergence tolerance
559 // Use the given solver controls
560 template<class Type>
561 SolverPerformance<Type> solve(fvMatrix<Type>&, const word&);
562 
563 //- Solve returning the solution statistics given convergence tolerance,
564 // deleting temporary matrix after solution.
565 // Use the given solver controls
566 template<class Type>
568 (
569  const tmp<fvMatrix<Type>>&,
570  const word&
571 );
572 
573 //- Solve returning the solution statistics given convergence tolerance
574 // Solver controls read fvSolution
575 template<class Type>
576 SolverPerformance<Type> solve(fvMatrix<Type>&);
577 
578 //- Solve returning the solution statistics given convergence tolerance,
579 // deleting temporary matrix after solution.
580 // Solver controls read fvSolution
581 template<class Type>
582 SolverPerformance<Type> solve(const tmp<fvMatrix<Type>>&);
583 
584 //- Return the correction form of the given matrix
585 // by subtracting the matrix multiplied by the current field
586 template<class Type>
587 tmp<fvMatrix<Type>> correction(const fvMatrix<Type>&);
588 
589 //- Return the correction form of the given temporary matrix
590 // by subtracting the matrix multiplied by the current field
591 template<class Type>
592 tmp<fvMatrix<Type>> correction(const tmp<fvMatrix<Type>>&);
593 
594 
595 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
596 
597 template<class Type>
598 tmp<fvMatrix<Type>> operator==
599 (
600  const fvMatrix<Type>&,
601  const fvMatrix<Type>&
602 );
603 
604 template<class Type>
605 tmp<fvMatrix<Type>> operator==
606 (
607  const tmp<fvMatrix<Type>>&,
608  const fvMatrix<Type>&
609 );
610 
611 template<class Type>
612 tmp<fvMatrix<Type>> operator==
613 (
614  const fvMatrix<Type>&,
615  const tmp<fvMatrix<Type>>&
616 );
617 
618 template<class Type>
619 tmp<fvMatrix<Type>> operator==
620 (
621  const tmp<fvMatrix<Type>>&,
622  const tmp<fvMatrix<Type>>&
623 );
624 
625 
626 template<class Type>
627 tmp<fvMatrix<Type>> operator==
628 (
629  const fvMatrix<Type>&,
631 );
632 
633 template<class Type>
634 tmp<fvMatrix<Type>> operator==
635 (
636  const fvMatrix<Type>&,
638 );
639 
640 template<class Type>
641 tmp<fvMatrix<Type>> operator==
642 (
643  const fvMatrix<Type>&,
645 );
646 
647 template<class Type>
648 tmp<fvMatrix<Type>> operator==
649 (
650  const tmp<fvMatrix<Type>>&,
652 );
653 
654 template<class Type>
655 tmp<fvMatrix<Type>> operator==
656 (
657  const tmp<fvMatrix<Type>>&,
659 );
660 
661 template<class Type>
662 tmp<fvMatrix<Type>> operator==
663 (
664  const tmp<fvMatrix<Type>>&,
666 );
667 
668 template<class Type>
669 tmp<fvMatrix<Type>> operator==
670 (
671  const fvMatrix<Type>&,
672  const dimensioned<Type>&
673 );
674 
675 template<class Type>
676 tmp<fvMatrix<Type>> operator==
677 (
678  const tmp<fvMatrix<Type>>&,
679  const dimensioned<Type>&
680 );
681 
682 
683 template<class Type>
684 tmp<fvMatrix<Type>> operator==
685 (
686  const fvMatrix<Type>&,
687  const zero&
688 );
689 
690 template<class Type>
691 tmp<fvMatrix<Type>> operator==
692 (
693  const tmp<fvMatrix<Type>>&,
694  const zero&
695 );
696 
697 
698 template<class Type>
699 tmp<fvMatrix<Type>> operator-
700 (
701  const fvMatrix<Type>&
702 );
703 
704 template<class Type>
705 tmp<fvMatrix<Type>> operator-
706 (
707  const tmp<fvMatrix<Type>>&
708 );
709 
710 
711 template<class Type>
712 tmp<fvMatrix<Type>> operator+
713 (
714  const fvMatrix<Type>&,
715  const fvMatrix<Type>&
716 );
717 
718 template<class Type>
719 tmp<fvMatrix<Type>> operator+
720 (
721  const tmp<fvMatrix<Type>>&,
722  const fvMatrix<Type>&
723 );
724 
725 template<class Type>
726 tmp<fvMatrix<Type>> operator+
727 (
728  const fvMatrix<Type>&,
729  const tmp<fvMatrix<Type>>&
730 );
731 
732 template<class Type>
733 tmp<fvMatrix<Type>> operator+
734 (
735  const tmp<fvMatrix<Type>>&,
736  const tmp<fvMatrix<Type>>&
737 );
738 
739 
740 template<class Type>
741 tmp<fvMatrix<Type>> operator+
742 (
743  const fvMatrix<Type>&,
745 );
746 
747 template<class Type>
748 tmp<fvMatrix<Type>> operator+
749 (
750  const fvMatrix<Type>&,
752 );
753 
754 template<class Type>
755 tmp<fvMatrix<Type>> operator+
756 (
757  const fvMatrix<Type>&,
759 );
760 
761 template<class Type>
762 tmp<fvMatrix<Type>> operator+
763 (
764  const tmp<fvMatrix<Type>>&,
766 );
767 
768 template<class Type>
769 tmp<fvMatrix<Type>> operator+
770 (
771  const tmp<fvMatrix<Type>>&,
773 );
774 
775 template<class Type>
776 tmp<fvMatrix<Type>> operator+
777 (
778  const tmp<fvMatrix<Type>>&,
780 );
781 
782 template<class Type>
783 tmp<fvMatrix<Type>> operator+
784 (
786  const fvMatrix<Type>&
787 );
788 
789 template<class Type>
790 tmp<fvMatrix<Type>> operator+
791 (
793  const fvMatrix<Type>&
794 );
795 
796 template<class Type>
797 tmp<fvMatrix<Type>> operator+
798 (
800  const fvMatrix<Type>&
801 );
802 
803 template<class Type>
804 tmp<fvMatrix<Type>> operator+
805 (
806  const DimensionedField<Type, volMesh>&,
807  const tmp<fvMatrix<Type>>&
808 );
809 
810 template<class Type>
811 tmp<fvMatrix<Type>> operator+
812 (
813  const tmp<DimensionedField<Type, volMesh>>&,
814  const tmp<fvMatrix<Type>>&
815 );
816 
817 template<class Type>
818 tmp<fvMatrix<Type>> operator+
819 (
820  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
821  const tmp<fvMatrix<Type>>&
822 );
823 
824 
825 template<class Type>
826 tmp<fvMatrix<Type>> operator+
827 (
828  const fvMatrix<Type>&,
829  const dimensioned<Type>&
830 );
831 
832 template<class Type>
833 tmp<fvMatrix<Type>> operator+
834 (
835  const tmp<fvMatrix<Type>>&,
836  const dimensioned<Type>&
837 );
838 
839 template<class Type>
840 tmp<fvMatrix<Type>> operator+
841 (
842  const dimensioned<Type>&,
843  const fvMatrix<Type>&
844 );
845 
846 template<class Type>
847 tmp<fvMatrix<Type>> operator+
848 (
849  const dimensioned<Type>&,
850  const tmp<fvMatrix<Type>>&
851 );
852 
853 
854 template<class Type>
855 tmp<fvMatrix<Type>> operator-
856 (
857  const fvMatrix<Type>&,
858  const fvMatrix<Type>&
859 );
860 
861 template<class Type>
862 tmp<fvMatrix<Type>> operator-
863 (
864  const tmp<fvMatrix<Type>>&,
865  const fvMatrix<Type>&
866 );
867 
868 template<class Type>
869 tmp<fvMatrix<Type>> operator-
870 (
871  const fvMatrix<Type>&,
872  const tmp<fvMatrix<Type>>&
873 );
874 
875 template<class Type>
876 tmp<fvMatrix<Type>> operator-
877 (
878  const tmp<fvMatrix<Type>>&,
879  const tmp<fvMatrix<Type>>&
880 );
881 
882 
883 template<class Type>
884 tmp<fvMatrix<Type>> operator-
885 (
886  const fvMatrix<Type>&,
887  const DimensionedField<Type, volMesh>&
888 );
889 
890 template<class Type>
891 tmp<fvMatrix<Type>> operator-
892 (
893  const fvMatrix<Type>&,
894  const tmp<DimensionedField<Type, volMesh>>&
895 );
896 
897 template<class Type>
898 tmp<fvMatrix<Type>> operator-
899 (
900  const fvMatrix<Type>&,
901  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
902 );
903 
904 template<class Type>
905 tmp<fvMatrix<Type>> operator-
906 (
907  const tmp<fvMatrix<Type>>&,
908  const DimensionedField<Type, volMesh>&
909 );
910 
911 template<class Type>
912 tmp<fvMatrix<Type>> operator-
913 (
914  const tmp<fvMatrix<Type>>&,
915  const tmp<DimensionedField<Type, volMesh>>&
916 );
917 
918 template<class Type>
919 tmp<fvMatrix<Type>> operator-
920 (
921  const tmp<fvMatrix<Type>>&,
922  const tmp<GeometricField<Type, fvPatchField, volMesh>>&
923 );
924 
925 template<class Type>
926 tmp<fvMatrix<Type>> operator-
927 (
928  const DimensionedField<Type, volMesh>&,
929  const fvMatrix<Type>&
930 );
931 
932 template<class Type>
933 tmp<fvMatrix<Type>> operator-
934 (
935  const tmp<DimensionedField<Type, volMesh>>&,
936  const fvMatrix<Type>&
937 );
938 
939 template<class Type>
940 tmp<fvMatrix<Type>> operator-
941 (
942  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
943  const fvMatrix<Type>&
944 );
945 
946 template<class Type>
947 tmp<fvMatrix<Type>> operator-
948 (
949  const DimensionedField<Type, volMesh>&,
950  const tmp<fvMatrix<Type>>&
951 );
952 
953 template<class Type>
954 tmp<fvMatrix<Type>> operator-
955 (
956  const tmp<DimensionedField<Type, volMesh>>&,
957  const tmp<fvMatrix<Type>>&
958 );
959 
960 template<class Type>
961 tmp<fvMatrix<Type>> operator-
962 (
963  const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
964  const tmp<fvMatrix<Type>>&
965 );
966 
967 
968 template<class Type>
969 tmp<fvMatrix<Type>> operator-
970 (
971  const fvMatrix<Type>&,
972  const dimensioned<Type>&
973 );
974 
975 template<class Type>
976 tmp<fvMatrix<Type>> operator-
977 (
978  const tmp<fvMatrix<Type>>&,
979  const dimensioned<Type>&
980 );
981 
982 template<class Type>
983 tmp<fvMatrix<Type>> operator-
984 (
985  const dimensioned<Type>&,
986  const fvMatrix<Type>&
987 );
988 
989 template<class Type>
990 tmp<fvMatrix<Type>> operator-
991 (
992  const dimensioned<Type>&,
993  const tmp<fvMatrix<Type>>&
994 );
995 
996 
997 template<class Type>
998 tmp<fvMatrix<Type>> operator*
999 (
1000  const volScalarField::Internal&,
1001  const fvMatrix<Type>&
1002 );
1003 
1004 template<class Type>
1005 tmp<fvMatrix<Type>> operator*
1006 (
1008  const fvMatrix<Type>&
1009 );
1010 
1011 template<class Type>
1012 tmp<fvMatrix<Type>> operator*
1013 (
1014  const tmp<volScalarField>&,
1015  const fvMatrix<Type>&
1016 );
1017 
1018 template<class Type>
1019 tmp<fvMatrix<Type>> operator*
1020 (
1021  const volScalarField::Internal&,
1022  const tmp<fvMatrix<Type>>&
1023 );
1024 
1025 template<class Type>
1026 tmp<fvMatrix<Type>> operator*
1027 (
1028  const tmp<volScalarField::Internal>&,
1029  const tmp<fvMatrix<Type>>&
1030 );
1031 
1032 template<class Type>
1033 tmp<fvMatrix<Type>> operator*
1034 (
1035  const tmp<volScalarField>&,
1036  const tmp<fvMatrix<Type>>&
1037 );
1038 
1039 
1040 template<class Type>
1041 tmp<fvMatrix<Type>> operator*
1042 (
1043  const dimensioned<scalar>&,
1044  const fvMatrix<Type>&
1045 );
1046 
1047 template<class Type>
1048 tmp<fvMatrix<Type>> operator*
1049 (
1050  const dimensioned<scalar>&,
1051  const tmp<fvMatrix<Type>>&
1052 );
1053 
1054 template<class Type>
1055 tmp<fvMatrix<Type>> operator*
1056 (
1057  const volScalarField::Internal&,
1058  const fvMatrix<Type>&
1059 );
1060 
1061 template<class Type>
1062 tmp<fvMatrix<Type>> operator*
1063 (
1064  const tmp<volScalarField::Internal>&,
1065  const fvMatrix<Type>&
1066 );
1067 
1068 template<class Type>
1069 tmp<fvMatrix<Type>> operator*
1070 (
1071  const tmp<volScalarField>&,
1072  const fvMatrix<Type>&
1073 );
1074 
1075 template<class Type>
1076 tmp<fvMatrix<Type>> operator*
1077 (
1078  const volScalarField::Internal&,
1079  const tmp<fvMatrix<Type>>&
1080 );
1081 
1082 template<class Type>
1083 tmp<fvMatrix<Type>> operator*
1084 (
1085  const tmp<volScalarField::Internal>&,
1086  const tmp<fvMatrix<Type>>&
1087 );
1088 
1089 template<class Type>
1090 tmp<fvMatrix<Type>> operator*
1091 (
1092  const tmp<volScalarField>&,
1093  const tmp<fvMatrix<Type>>&
1094 );
1095 
1096 
1097 template<class Type>
1098 tmp<fvMatrix<Type>> operator*
1099 (
1100  const dimensioned<scalar>&,
1101  const fvMatrix<Type>&
1102 );
1103 
1104 template<class Type>
1105 tmp<fvMatrix<Type>> operator*
1106 (
1107  const dimensioned<scalar>&,
1108  const tmp<fvMatrix<Type>>&
1109 );
1110 
1111 
1112 template<class Type>
1113 tmp<fvMatrix<Type>> operator/
1114 (
1115  const fvMatrix<Type>&,
1116  const volScalarField::Internal&
1117 );
1118 
1119 template<class Type>
1120 tmp<fvMatrix<Type>> operator/
1121 (
1122  const fvMatrix<Type>&,
1123  const tmp<volScalarField::Internal>&
1124 );
1125 
1126 template<class Type>
1127 tmp<fvMatrix<Type>> operator/
1128 (
1129  const fvMatrix<Type>&,
1130  const tmp<volScalarField>&
1131 );
1132 
1133 template<class Type>
1134 tmp<fvMatrix<Type>> operator/
1135 (
1136  const tmp<fvMatrix<Type>>&,
1137  const volScalarField::Internal&
1138 );
1139 
1140 template<class Type>
1141 tmp<fvMatrix<Type>> operator/
1142 (
1143  const tmp<fvMatrix<Type>>&,
1144  const tmp<volScalarField::Internal>&
1145 );
1146 
1147 template<class Type>
1148 tmp<fvMatrix<Type>> operator/
1149 (
1150  const tmp<fvMatrix<Type>>&,
1151  const tmp<volScalarField>&
1152 );
1153 
1154 
1155 template<class Type>
1156 tmp<fvMatrix<Type>> operator/
1157 (
1158  const fvMatrix<Type>&,
1159  const dimensioned<scalar>&
1160 );
1161 
1162 template<class Type>
1163 tmp<fvMatrix<Type>> operator/
1164 (
1165  const tmp<fvMatrix<Type>>&,
1166  const dimensioned<scalar>&
1167 );
1168 
1169 
1170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1171 
1172 } // End namespace Foam
1173 
1174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1175 
1176 #ifdef NoRepository
1177  #include "fvMatrix.C"
1178 #endif
1179 
1180 // Specialisation for scalars
1181 #include "fvScalarMatrix.H"
1182 
1183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1184 
1185 #endif
1186 
1187 // ************************************************************************* //
Foam::surfaceFields.
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: fvMatrix.C:38
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
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:59
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:111
Reference counter for various OpenFOAM components.
Definition: refCount.H:49
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:672
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:980
FieldField< Field, Type > & internalCoeffs()
fvBoundary scalar field containing pseudo-matrix coeffs
Definition: fvMatrix.H:304
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
uint8_t direction
Definition: direction.H:45
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:129
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:506
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:760
void operator-=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1029
Solver class returned by the solver function.
Definition: fvMatrix.H:220
Generic GeometricField class.
Generic dimensioned Type class.
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:815
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: fvMatrix.C:75
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:705
void operator=(const fvMatrix< Type > &)
Definition: fvMatrix.C:936
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:322
fvSolver(fvMatrix< Type > &fvMat, autoPtr< lduMatrix::solver > sol)
Definition: fvMatrix.H:230
Dimension set for the base types.
Definition: dimensionSet.H:120
SolverPerformance< Type > solve()
Solve returning the solution statistics.
const cellShapeList & cells
void operator/=(const volScalarField::Internal &)
Definition: fvMatrix.C:1235
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...
void setValues(const labelUList &cells, const UList< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:484
FieldField< Field, Type > & boundaryCoeffs()
fvBoundary scalar field containing pseudo-matrix coeffs
Definition: fvMatrix.H:311
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:465
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:453
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:292
void operator*=(const volScalarField::Internal &)
Definition: fvMatrix.C:1164
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:995
ClassName("fvMatrix")
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:714
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:319
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:738
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 setValuesFromList(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:178
void checkMethod(const fvMatrix< Type > &, const fvMatrix< Type > &, const char *)
Definition: fvMatrix.C:1308
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:145
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:856
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:692
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:268
Namespace for OpenFOAM.
const dimensionSet & dimensions() const
Definition: fvMatrix.H:287