added matrix correction operation
[OpenFOAM-1.5.x.git] / src / finiteVolume / fvMatrices / fvMatrix / fvMatrix.H
blob056ce560963542489acba70578645007a5fe2141
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::fvMatrix
28 Description
29     Finite-Volume matrix.
31 SourceFiles
32     fvMatrix.C
33     fvMatrixSolve.C
34     fvScalarMatrix.C
36 \*---------------------------------------------------------------------------*/
38 #ifndef fvMatrix_H
39 #define fvMatrix_H
41 #include "volFields.H"
42 #include "surfaceFields.H"
43 #include "lduMatrix.H"
44 #include "tmp.H"
45 #include "autoPtr.H"
46 #include "dimensionedTypes.H"
47 #include "zeroField.H"
48 #include "className.H"
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 namespace Foam
55 // Forward declaration of friend functions and operators
57 template<class Type>
58 class fvMatrix;
60 template<class Type>
61 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
63     const fvMatrix<Type>&,
64     const DimensionedField<Type, volMesh>&
67 template<class Type>
68 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
70     const fvMatrix<Type>&,
71     const tmp<DimensionedField<Type, volMesh> >&
74 template<class Type>
75 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
77     const fvMatrix<Type>&,
78     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
81 template<class Type>
82 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
84     const tmp<fvMatrix<Type> >&,
85     const DimensionedField<Type, volMesh>&
88 template<class Type>
89 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
91     const tmp<fvMatrix<Type> >&,
92     const tmp<DimensionedField<Type, volMesh> >&
95 template<class Type>
96 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
98     const tmp<fvMatrix<Type> >&,
99     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
102 template<class Type>
103 Ostream& operator<<(Ostream&, const fvMatrix<Type>&);
106 /*---------------------------------------------------------------------------*\
107                            Class fvMatrix Declaration
108 \*---------------------------------------------------------------------------*/
110 template<class Type>
111 class fvMatrix
113     public refCount,
114     public lduMatrix
116 public:
118     // Private data
120         // Reference to GeometricField<Type, fvPatchField, volMesh>
121         GeometricField<Type, fvPatchField, volMesh>& psi_;
123         //- Dimension set
124         dimensionSet dimensions_;
126         //- Source term
127         Field<Type> source_;
129         //- Boundary scalar field containing pseudo-matrix coeffs
130         //  for internal cells
131         FieldField<Field, Type> internalCoeffs_;
133         //- Boundary scalar field containing pseudo-matrix coeffs
134         //  for boundary cells
135         FieldField<Field, Type> boundaryCoeffs_;
138         //- Face flux field for non-orthogonal correction
139         mutable GeometricField<Type, fvsPatchField, surfaceMesh>
140             *faceFluxCorrectionPtr_;
143     // Private member functions
145         //- Add patch contribution to internal field
146         template<class Type2>
147         void addToInternalField
148         (
149             const unallocLabelList& addr,
150             const Field<Type2>& pf,
151             Field<Type2>& intf
152         ) const;
154         template<class Type2>
155         void addToInternalField
156         (
157             const unallocLabelList& addr,
158             const tmp<Field<Type2> >& tpf,
159             Field<Type2>& intf
160         ) const;
162         //- Subtract patch contribution from internal field
163         template<class Type2>
164         void subtractFromInternalField
165         (
166             const unallocLabelList& addr,
167             const Field<Type2>& pf,
168             Field<Type2>& intf
169         ) const;
171         template<class Type2>
172         void subtractFromInternalField
173         (
174             const unallocLabelList& addr,
175             const tmp<Field<Type2> >& tpf,
176             Field<Type2>& intf
177         ) const;
180         // Matrix completion functionality
182             void addBoundaryDiag
183             (
184                 scalarField& diag,
185                 const direction cmpt
186             ) const;
188             void addCmptAvBoundaryDiag(scalarField& diag) const;
190             void addBoundarySource
191             (
192                 Field<Type>& source,
193                 const bool couples=true
194             ) const;
197 public:
199     //- Solver class returned by the solver function
200     //  used for systems in which it is useful to cache the solver for reuse
201     //  e.g. is the solver is potentialy expensive to construct (AMG) and can
202     //  be used several times (PISO)
203     class fvSolver
204     {
205         fvMatrix<Type>& fvMat_;
207         autoPtr<lduMatrix::solver> solver_;
209     public:
211         // Constructors
213             fvSolver(fvMatrix<Type>& fvMat, autoPtr<lduMatrix::solver> sol)
214             :
215                 fvMat_(fvMat),
216                 solver_(sol)
217             {}
220         // Member functions
222             //- Solve returning the solution statistics.
223             //  Solver controls read from Istream
224             lduMatrix::solverPerformance solve(Istream&);
226             //- Solve returning the solution statistics.
227             //  Solver controls read from fvSolution
228             lduMatrix::solverPerformance solve();
229     };
232     ClassName("fvMatrix");
235     // Constructors
237         //- Construct given a field to solve for
238         fvMatrix
239         (
240             GeometricField<Type, fvPatchField, volMesh>&,
241             const dimensionSet&
242         );
244         //- Construct as copy
245         fvMatrix(const fvMatrix<Type>&);
247         //- Construct as copy of tmp<fvMatrix<Type> > deleting argument
248 #       ifdef ConstructFromTmp
249         fvMatrix(const tmp<fvMatrix<Type> >&);
250 #       endif
252         //- Construct from Istream given field to solve for
253         fvMatrix(GeometricField<Type, fvPatchField, volMesh>&, Istream&);
256     // Destructor
258         virtual ~fvMatrix();
261     // Member functions
263         // Access
265             const GeometricField<Type, fvPatchField, volMesh>& psi() const
266             {
267                 return psi_;
268             }
270             GeometricField<Type, fvPatchField, volMesh>& psi()
271             {
272                 return psi_;
273             }
275             const dimensionSet& dimensions() const
276             {
277                 return dimensions_;
278             }
280             Field<Type>& source()
281             {
282                 return source_;
283             }
285             const Field<Type>& source() const
286             {
287                 return source_;
288             }
290             //- fvBoundary scalar field containing pseudo-matrix coeffs
291             //  for internal cells
292             FieldField<Field, Type>& internalCoeffs()
293             {
294                 return internalCoeffs_;
295             }
297             //- fvBoundary scalar field containing pseudo-matrix coeffs
298             //  for boundary cells
299             FieldField<Field, Type>& boundaryCoeffs()
300             {
301                 return boundaryCoeffs_;
302             }
305             //- Declare return type of the faceFluxCorrectionPtr() function
306             typedef GeometricField<Type, fvsPatchField, surfaceMesh>
307                 *surfaceTypeFieldPtr;
309             //- Return pointer to face-flux non-orthogonal correction field
310             surfaceTypeFieldPtr& faceFluxCorrectionPtr()
311             {
312                 return faceFluxCorrectionPtr_;
313             }
316         // Operations
318             //- Set solution in given cells and eliminate corresponding
319             //  equations from the matrix
320             void setValues
321             (
322                 const labelList& cells,
323                 const Field<Type>& values
324             );
326             //- Set reference level for solution
327             void setReference
328             (
329                 const label celli,
330                 const Type& value,
331                 const bool forceReference = false
332             );
334             //- Set reference level for a component of the solution
335             //  on a given patch face
336             void setComponentReference
337             (
338                 const label patchi,
339                 const label facei,
340                 const direction cmpt,
341                 const scalar value
342             );
344             //- Relax matrix (for steady-state solution).
345             //  alpha = 1 : diagonally equal
346             //  alpha < 1 : diagonally dominant
347             //  alpha = 0 : do nothing
348             void relax(const scalar alpha);
350             //- Relax matrix (for steady-state solution).
351             //  alpha is read from controlDict
352             void relax();
354             //- Construct and return the solver
355             //  Solver controls read from Istream
356             autoPtr<fvSolver> solver(Istream&);
358             //- Construct and return the solver
359             //  Solver controls read from fvSolution
360             autoPtr<fvSolver> solver();
362             //- Solve returning the solution statistics.
363             //  Solver controls read from Istream
364             lduMatrix::solverPerformance solve(Istream&);
366             //- Solve returning the solution statistics.
367             //  Solver controls read from fvSolution
368             lduMatrix::solverPerformance solve();
370             //- Return the matrix residual
371             tmp<Field<Type> > residual() const;
373             //- Return the matrix scalar diagonal
374             tmp<scalarField> D() const;
376             //- Return the matrix Type diagonal
377             tmp<Field<Type> > DD() const;
379             //- Return the central coefficient
380             tmp<volScalarField> A() const;
382             //- Return the H operation source
383             tmp<GeometricField<Type, fvPatchField, volMesh> > H() const;
385             //- Return H(1)
386             tmp<volScalarField> H1() const;
388             //- Return the face-flux field from the matrix
389             tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
390             flux() const;
393     // Member operators
395         void operator=(const fvMatrix<Type>&);
396         void operator=(const tmp<fvMatrix<Type> >&);
398         void negate();
400         void operator+=(const fvMatrix<Type>&);
401         void operator+=(const tmp<fvMatrix<Type> >&);
403         void operator-=(const fvMatrix<Type>&);
404         void operator-=(const tmp<fvMatrix<Type> >&);
406         void operator+=
407         (
408             const DimensionedField<Type, volMesh>&
409         );
410         void operator+=
411         (
412             const tmp<DimensionedField<Type, volMesh> >&
413         );
414         void operator+=
415         (
416             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
417         );
419         void operator-=
420         (
421             const DimensionedField<Type, volMesh>&
422         );
423         void operator-=
424         (
425             const tmp<DimensionedField<Type, volMesh> >&
426         );
427         void operator-=
428         (
429             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
430         );
432         void operator+=(const dimensioned<Type>&);
433         void operator-=(const dimensioned<Type>&);
435         void operator+=(const zeroField&);
436         void operator-=(const zeroField&);
438         void operator*=(const DimensionedField<scalar, volMesh>&);
439         void operator*=(const tmp<DimensionedField<scalar, volMesh> >&);
440         void operator*=(const tmp<volScalarField>&);
442         void operator*=(const dimensioned<scalar>&);
445     // Friend operators
447         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
448         operator& <Type>
449         (
450             const fvMatrix<Type>&,
451             const DimensionedField<Type, volMesh>&
452         );
454         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
455         operator& <Type>
456         (
457             const fvMatrix<Type>&,
458             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
459         );
461         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
462         operator& <Type>
463         (
464             const tmp<fvMatrix<Type> >&,
465             const DimensionedField<Type, volMesh>&
466         );
468         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
469         operator& <Type>
470         (
471             const tmp<fvMatrix<Type> >&,
472             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
473         );
476     // Ostream operator
478         friend Ostream& operator<< <Type>
479         (
480             Ostream&,
481             const fvMatrix<Type>&
482         );
486 // * * * * * * * * * * * * * * * Global functions  * * * * * * * * * * * * * //
488 template<class Type>
489 void checkMethod
491     const fvMatrix<Type>&,
492     const fvMatrix<Type>&,
493     const char*
496 template<class Type>
497 void checkMethod
499     const fvMatrix<Type>&,
500     const DimensionedField<Type, volMesh>&,
501     const char*
504 template<class Type>
505 void checkMethod
507     const fvMatrix<Type>&,
508     const dimensioned<Type>&,
509     const char*
513 //- Solve returning the solution statistics given convergence tolerance
514 //  Solver controls read Istream
515 template<class Type>
516 lduMatrix::solverPerformance solve(fvMatrix<Type>&, Istream&);
519 //- Solve returning the solution statistics given convergence tolerance,
520 //  deleting temporary matrix after solution.
521 //  Solver controls read Istream
522 template<class Type>
523 lduMatrix::solverPerformance solve(const tmp<fvMatrix<Type> >&, Istream&);
526 //- Solve returning the solution statistics given convergence tolerance
527 //  Solver controls read fvSolution
528 template<class Type>
529 lduMatrix::solverPerformance solve(fvMatrix<Type>&);
532 //- Solve returning the solution statistics given convergence tolerance,
533 //  deleting temporary matrix after solution.
534 //  Solver controls read fvSolution
535 template<class Type>
536 lduMatrix::solverPerformance solve(const tmp<fvMatrix<Type> >&);
539 //- Return the correction form of the given matrix
540 //  by subtracting the matrix multiplied by the current field
541 template<class Type>
542 tmp<fvMatrix<Type> > correction(const fvMatrix<Type>&);
545 //- Return the correction form of the given temporary matrix
546 //  by subtracting the matrix multiplied by the current field
547 template<class Type>
548 tmp<fvMatrix<Type> > correction(const tmp<fvMatrix<Type> >&);
551 // * * * * * * * * * * * * * * * Global operators  * * * * * * * * * * * * * //
553 template<class Type>
554 tmp<fvMatrix<Type> > operator==
556     const fvMatrix<Type>&,
557     const fvMatrix<Type>&
560 template<class Type>
561 tmp<fvMatrix<Type> > operator==
563     const tmp<fvMatrix<Type> >&,
564     const fvMatrix<Type>&
567 template<class Type>
568 tmp<fvMatrix<Type> > operator==
570     const fvMatrix<Type>&,
571     const tmp<fvMatrix<Type> >&
574 template<class Type>
575 tmp<fvMatrix<Type> > operator==
577     const tmp<fvMatrix<Type> >&,
578     const tmp<fvMatrix<Type> >&
582 template<class Type>
583 tmp<fvMatrix<Type> > operator==
585     const fvMatrix<Type>&,
586     const DimensionedField<Type, volMesh>&
589 template<class Type>
590 tmp<fvMatrix<Type> > operator==
592     const fvMatrix<Type>&,
593     const tmp<DimensionedField<Type, volMesh> >&
596 template<class Type>
597 tmp<fvMatrix<Type> > operator==
599     const fvMatrix<Type>&,
600     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
603 template<class Type>
604 tmp<fvMatrix<Type> > operator==
606     const tmp<fvMatrix<Type> >&,
607     const DimensionedField<Type, volMesh>&
610 template<class Type>
611 tmp<fvMatrix<Type> > operator==
613     const tmp<fvMatrix<Type> >&,
614     const tmp<DimensionedField<Type, volMesh> >&
617 template<class Type>
618 tmp<fvMatrix<Type> > operator==
620     const tmp<fvMatrix<Type> >&,
621     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
624 template<class Type>
625 tmp<fvMatrix<Type> > operator==
627     const fvMatrix<Type>&,
628     const dimensioned<Type>&
631 template<class Type>
632 tmp<fvMatrix<Type> > operator==
634     const tmp<fvMatrix<Type> >&,
635     const dimensioned<Type>&
639 template<class Type>
640 tmp<fvMatrix<Type> > operator==
642     const fvMatrix<Type>&,
643     const zeroField&
646 template<class Type>
647 tmp<fvMatrix<Type> > operator==
649     const tmp<fvMatrix<Type> >&,
650     const zeroField&
654 template<class Type>
655 tmp<fvMatrix<Type> > operator-
657     const fvMatrix<Type>&
660 template<class Type>
661 tmp<fvMatrix<Type> > operator-
663     const tmp<fvMatrix<Type> >&
667 template<class Type>
668 tmp<fvMatrix<Type> > operator+
670     const fvMatrix<Type>&,
671     const fvMatrix<Type>&
674 template<class Type>
675 tmp<fvMatrix<Type> > operator+
677     const tmp<fvMatrix<Type> >&,
678     const fvMatrix<Type>&
681 template<class Type>
682 tmp<fvMatrix<Type> > operator+
684     const fvMatrix<Type>&,
685     const tmp<fvMatrix<Type> >&
688 template<class Type>
689 tmp<fvMatrix<Type> > operator+
691     const tmp<fvMatrix<Type> >&,
692     const tmp<fvMatrix<Type> >&
696 template<class Type>
697 tmp<fvMatrix<Type> > operator+
699     const fvMatrix<Type>&,
700     const DimensionedField<Type, volMesh>&
703 template<class Type>
704 tmp<fvMatrix<Type> > operator+
706     const fvMatrix<Type>&,
707     const tmp<DimensionedField<Type, volMesh> >&
710 template<class Type>
711 tmp<fvMatrix<Type> > operator+
713     const fvMatrix<Type>&,
714     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
717 template<class Type>
718 tmp<fvMatrix<Type> > operator+
720     const tmp<fvMatrix<Type> >&,
721     const DimensionedField<Type, volMesh>&
724 template<class Type>
725 tmp<fvMatrix<Type> > operator+
727     const tmp<fvMatrix<Type> >&,
728     const tmp<DimensionedField<Type, volMesh> >&
731 template<class Type>
732 tmp<fvMatrix<Type> > operator+
734     const tmp<fvMatrix<Type> >&,
735     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
738 template<class Type>
739 tmp<fvMatrix<Type> > operator+
741     const DimensionedField<Type, volMesh>&,
742     const fvMatrix<Type>&
745 template<class Type>
746 tmp<fvMatrix<Type> > operator+
748     const tmp<DimensionedField<Type, volMesh> >&,
749     const fvMatrix<Type>&
752 template<class Type>
753 tmp<fvMatrix<Type> > operator+
755     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
756     const fvMatrix<Type>&
759 template<class Type>
760 tmp<fvMatrix<Type> > operator+
762     const DimensionedField<Type, volMesh>&,
763     const tmp<fvMatrix<Type> >&
766 template<class Type>
767 tmp<fvMatrix<Type> > operator+
769     const tmp<DimensionedField<Type, volMesh> >&,
770     const tmp<fvMatrix<Type> >&
773 template<class Type>
774 tmp<fvMatrix<Type> > operator+
776     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
777     const tmp<fvMatrix<Type> >&
781 template<class Type>
782 tmp<fvMatrix<Type> > operator+
784     const fvMatrix<Type>&,
785     const dimensioned<Type>&
788 template<class Type>
789 tmp<fvMatrix<Type> > operator+
791     const tmp<fvMatrix<Type> >&,
792     const dimensioned<Type>&
795 template<class Type>
796 tmp<fvMatrix<Type> > operator+
798     const dimensioned<Type>&,
799     const fvMatrix<Type>&
802 template<class Type>
803 tmp<fvMatrix<Type> > operator+
805     const dimensioned<Type>&,
806     const tmp<fvMatrix<Type> >&
810 template<class Type>
811 tmp<fvMatrix<Type> > operator-
813     const fvMatrix<Type>&,
814     const fvMatrix<Type>&
817 template<class Type>
818 tmp<fvMatrix<Type> > operator-
820     const tmp<fvMatrix<Type> >&,
821     const fvMatrix<Type>&
824 template<class Type>
825 tmp<fvMatrix<Type> > operator-
827     const fvMatrix<Type>&,
828     const tmp<fvMatrix<Type> >&
831 template<class Type>
832 tmp<fvMatrix<Type> > operator-
834     const tmp<fvMatrix<Type> >&,
835     const tmp<fvMatrix<Type> >&
839 template<class Type>
840 tmp<fvMatrix<Type> > operator-
842     const fvMatrix<Type>&,
843     const DimensionedField<Type, volMesh>&
846 template<class Type>
847 tmp<fvMatrix<Type> > operator-
849     const fvMatrix<Type>&,
850     const tmp<DimensionedField<Type, volMesh> >&
853 template<class Type>
854 tmp<fvMatrix<Type> > operator-
856     const fvMatrix<Type>&,
857     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
860 template<class Type>
861 tmp<fvMatrix<Type> > operator-
863     const tmp<fvMatrix<Type> >&,
864     const DimensionedField<Type, volMesh>&
867 template<class Type>
868 tmp<fvMatrix<Type> > operator-
870     const tmp<fvMatrix<Type> >&,
871     const tmp<DimensionedField<Type, volMesh> >&
874 template<class Type>
875 tmp<fvMatrix<Type> > operator-
877     const tmp<fvMatrix<Type> >&,
878     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
881 template<class Type>
882 tmp<fvMatrix<Type> > operator-
884     const DimensionedField<Type, volMesh>&,
885     const fvMatrix<Type>&
888 template<class Type>
889 tmp<fvMatrix<Type> > operator-
891     const tmp<DimensionedField<Type, volMesh> >&,
892     const fvMatrix<Type>&
895 template<class Type>
896 tmp<fvMatrix<Type> > operator-
898     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
899     const fvMatrix<Type>&
902 template<class Type>
903 tmp<fvMatrix<Type> > operator-
905     const DimensionedField<Type, volMesh>&,
906     const tmp<fvMatrix<Type> >&
909 template<class Type>
910 tmp<fvMatrix<Type> > operator-
912     const tmp<DimensionedField<Type, volMesh> >&,
913     const tmp<fvMatrix<Type> >&
916 template<class Type>
917 tmp<fvMatrix<Type> > operator-
919     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
920     const tmp<fvMatrix<Type> >&
924 template<class Type>
925 tmp<fvMatrix<Type> > operator-
927     const fvMatrix<Type>&,
928     const dimensioned<Type>&
931 template<class Type>
932 tmp<fvMatrix<Type> > operator-
934     const tmp<fvMatrix<Type> >&,
935     const dimensioned<Type>&
938 template<class Type>
939 tmp<fvMatrix<Type> > operator-
941     const dimensioned<Type>&,
942     const fvMatrix<Type>&
945 template<class Type>
946 tmp<fvMatrix<Type> > operator-
948     const dimensioned<Type>&,
949     const tmp<fvMatrix<Type> >&
953 template<class Type>
954 tmp<fvMatrix<Type> > operator*
956     const DimensionedField<scalar, volMesh>&,
957     const fvMatrix<Type>&
960 template<class Type>
961 tmp<fvMatrix<Type> > operator*
963     const tmp<DimensionedField<scalar, volMesh> >&,
964     const fvMatrix<Type>&
967 template<class Type>
968 tmp<fvMatrix<Type> > operator*
970     const tmp<volScalarField>&,
971     const fvMatrix<Type>&
974 template<class Type>
975 tmp<fvMatrix<Type> > operator*
977     const DimensionedField<scalar, volMesh>&,
978     const tmp<fvMatrix<Type> >&
981 template<class Type>
982 tmp<fvMatrix<Type> > operator*
984     const tmp<DimensionedField<scalar, volMesh> >&,
985     const tmp<fvMatrix<Type> >&
988 template<class Type>
989 tmp<fvMatrix<Type> > operator*
991     const tmp<volScalarField>&,
992     const tmp<fvMatrix<Type> >&
996 template<class Type>
997 tmp<fvMatrix<Type> > operator*
999     const dimensioned<scalar>&,
1000     const fvMatrix<Type>&
1003 template<class Type>
1004 tmp<fvMatrix<Type> > operator*
1006     const dimensioned<scalar>&,
1007     const tmp<fvMatrix<Type> >&
1011 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1013 } // End namespace Foam
1015 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1017 #ifdef NoRepository
1018 #   include "fvMatrix.C"
1019 #endif
1021 // Specialisation for scalars
1022 #include "fvScalarMatrix.H"
1024 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1026 #endif
1028 // ************************************************************************* //