initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / fvMatrices / fvMatrix / fvMatrix.H
blobf9c82a2d58cf1228795bc7112a124fd6ffc8eaa1
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 cell,
330                 const Type& value
331             );
333             //- Set reference level for a component of the solution
334             //  on a given patch face
335             void setComponentReference
336             (
337                 const label patchi,
338                 const label facei,
339                 const direction cmpt,
340                 const scalar value
341             );
343             //- Relax matrix (for steady-state solution).
344             //  alpha = 1 : diagonally equal
345             //  alpha < 1 : diagonally dominant
346             //  alpha = 0 : do nothing
347             void relax(const scalar alpha);
349             //- Relax matrix (for steady-state solution).
350             //  alpha is read from controlDict
351             void relax();
353             //- Construct and return the solver
354             //  Solver controls read from Istream
355             autoPtr<fvSolver> solver(Istream&);
357             //- Construct and return the solver
358             //  Solver controls read from fvSolution
359             autoPtr<fvSolver> solver();
361             //- Solve returning the solution statistics.
362             //  Solver controls read from Istream
363             lduMatrix::solverPerformance solve(Istream&);
365             //- Solve returning the solution statistics.
366             //  Solver controls read from fvSolution
367             lduMatrix::solverPerformance solve();
369             //- Return the matrix residual
370             tmp<Field<Type> > residual() const;
372             //- Return the matrix scalar diagonal
373             tmp<scalarField> D() const;
375             //- Return the matrix Type diagonal
376             tmp<Field<Type> > DD() const;
378             //- Return the central coefficient
379             tmp<volScalarField> A() const;
381             //- Return the H operation source
382             tmp<GeometricField<Type, fvPatchField, volMesh> > H() const;
384             //- Return H(1)
385             tmp<volScalarField> H1() const;
387             //- Return the face-flux field from the matrix
388             tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
389             flux() const;
392     // Member operators
394         void operator=(const fvMatrix<Type>&);
395         void operator=(const tmp<fvMatrix<Type> >&);
397         void negate();
399         void operator+=(const fvMatrix<Type>&);
400         void operator+=(const tmp<fvMatrix<Type> >&);
402         void operator-=(const fvMatrix<Type>&);
403         void operator-=(const tmp<fvMatrix<Type> >&);
405         void operator+=
406         (
407             const DimensionedField<Type, volMesh>&
408         );
409         void operator+=
410         (
411             const tmp<DimensionedField<Type, volMesh> >&
412         );
413         void operator+=
414         (
415             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
416         );
418         void operator-=
419         (
420             const DimensionedField<Type, volMesh>&
421         );
422         void operator-=
423         (
424             const tmp<DimensionedField<Type, volMesh> >&
425         );
426         void operator-=
427         (
428             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
429         );
431         void operator+=(const dimensioned<Type>&);
432         void operator-=(const dimensioned<Type>&);
434         void operator+=(const zeroField&);
435         void operator-=(const zeroField&);
437         void operator*=(const DimensionedField<scalar, volMesh>&);
438         void operator*=(const tmp<DimensionedField<scalar, volMesh> >&);
439         void operator*=(const tmp<volScalarField>&);
441         void operator*=(const dimensioned<scalar>&);
444     // Friend operators
446         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
447         operator& <Type>
448         (
449             const fvMatrix<Type>&,
450             const DimensionedField<Type, volMesh>&
451         );
453         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
454         operator& <Type>
455         (
456             const fvMatrix<Type>&,
457             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
458         );
460         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
461         operator& <Type>
462         (
463             const tmp<fvMatrix<Type> >&,
464             const DimensionedField<Type, volMesh>&
465         );
467         friend tmp<GeometricField<Type, fvPatchField, volMesh> >
468         operator& <Type>
469         (
470             const tmp<fvMatrix<Type> >&,
471             const tmp<GeometricField<Type, fvPatchField, volMesh> >&
472         );
475     // Ostream operator
477         friend Ostream& operator<< <Type>
478         (
479             Ostream&,
480             const fvMatrix<Type>&
481         );
485 // * * * * * * * * * * * * * * * Global functions  * * * * * * * * * * * * * //
487 template<class Type>
488 void checkMethod
490     const fvMatrix<Type>&,
491     const fvMatrix<Type>&,
492     const char*
495 template<class Type>
496 void checkMethod
498     const fvMatrix<Type>&,
499     const DimensionedField<Type, volMesh>&,
500     const char*
503 template<class Type>
504 void checkMethod
506     const fvMatrix<Type>&,
507     const dimensioned<Type>&,
508     const char*
512 //- Solve returning the solution statistics given convergence tolerance
513 //  Solver controls read Istream
514 template<class Type>
515 lduMatrix::solverPerformance solve(fvMatrix<Type>&, Istream&);
518 //- Solve returning the solution statistics given convergence tolerance,
519 //  deleting temporary matrix after solution.
520 //  Solver controls read Istream
521 template<class Type>
522 lduMatrix::solverPerformance solve(const tmp<fvMatrix<Type> >&, Istream&);
525 //- Solve returning the solution statistics given convergence tolerance
526 //  Solver controls read fvSolution
527 template<class Type>
528 lduMatrix::solverPerformance solve(fvMatrix<Type>&);
531 //- Solve returning the solution statistics given convergence tolerance,
532 //  deleting temporary matrix after solution.
533 //  Solver controls read fvSolution
534 template<class Type>
535 lduMatrix::solverPerformance solve(const tmp<fvMatrix<Type> >&);
538 // * * * * * * * * * * * * * * * Global operators  * * * * * * * * * * * * * //
540 template<class Type>
541 tmp<fvMatrix<Type> > operator==
543     const fvMatrix<Type>&,
544     const fvMatrix<Type>&
547 template<class Type>
548 tmp<fvMatrix<Type> > operator==
550     const tmp<fvMatrix<Type> >&,
551     const fvMatrix<Type>&
554 template<class Type>
555 tmp<fvMatrix<Type> > operator==
557     const fvMatrix<Type>&,
558     const tmp<fvMatrix<Type> >&
561 template<class Type>
562 tmp<fvMatrix<Type> > operator==
564     const tmp<fvMatrix<Type> >&,
565     const tmp<fvMatrix<Type> >&
569 template<class Type>
570 tmp<fvMatrix<Type> > operator==
572     const fvMatrix<Type>&,
573     const DimensionedField<Type, volMesh>&
576 template<class Type>
577 tmp<fvMatrix<Type> > operator==
579     const fvMatrix<Type>&,
580     const tmp<DimensionedField<Type, volMesh> >&
583 template<class Type>
584 tmp<fvMatrix<Type> > operator==
586     const fvMatrix<Type>&,
587     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
590 template<class Type>
591 tmp<fvMatrix<Type> > operator==
593     const tmp<fvMatrix<Type> >&,
594     const DimensionedField<Type, volMesh>&
597 template<class Type>
598 tmp<fvMatrix<Type> > operator==
600     const tmp<fvMatrix<Type> >&,
601     const tmp<DimensionedField<Type, volMesh> >&
604 template<class Type>
605 tmp<fvMatrix<Type> > operator==
607     const tmp<fvMatrix<Type> >&,
608     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
611 template<class Type>
612 tmp<fvMatrix<Type> > operator==
614     const fvMatrix<Type>&,
615     const dimensioned<Type>&
618 template<class Type>
619 tmp<fvMatrix<Type> > operator==
621     const tmp<fvMatrix<Type> >&,
622     const dimensioned<Type>&
626 template<class Type>
627 tmp<fvMatrix<Type> > operator==
629     const fvMatrix<Type>&,
630     const zeroField&
633 template<class Type>
634 tmp<fvMatrix<Type> > operator==
636     const tmp<fvMatrix<Type> >&,
637     const zeroField&
641 template<class Type>
642 tmp<fvMatrix<Type> > operator-
644     const fvMatrix<Type>&
647 template<class Type>
648 tmp<fvMatrix<Type> > operator-
650     const tmp<fvMatrix<Type> >&
654 template<class Type>
655 tmp<fvMatrix<Type> > operator+
657     const fvMatrix<Type>&,
658     const fvMatrix<Type>&
661 template<class Type>
662 tmp<fvMatrix<Type> > operator+
664     const tmp<fvMatrix<Type> >&,
665     const fvMatrix<Type>&
668 template<class Type>
669 tmp<fvMatrix<Type> > operator+
671     const fvMatrix<Type>&,
672     const tmp<fvMatrix<Type> >&
675 template<class Type>
676 tmp<fvMatrix<Type> > operator+
678     const tmp<fvMatrix<Type> >&,
679     const tmp<fvMatrix<Type> >&
683 template<class Type>
684 tmp<fvMatrix<Type> > operator+
686     const fvMatrix<Type>&,
687     const DimensionedField<Type, volMesh>&
690 template<class Type>
691 tmp<fvMatrix<Type> > operator+
693     const fvMatrix<Type>&,
694     const tmp<DimensionedField<Type, volMesh> >&
697 template<class Type>
698 tmp<fvMatrix<Type> > operator+
700     const fvMatrix<Type>&,
701     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
704 template<class Type>
705 tmp<fvMatrix<Type> > operator+
707     const tmp<fvMatrix<Type> >&,
708     const DimensionedField<Type, volMesh>&
711 template<class Type>
712 tmp<fvMatrix<Type> > operator+
714     const tmp<fvMatrix<Type> >&,
715     const tmp<DimensionedField<Type, volMesh> >&
718 template<class Type>
719 tmp<fvMatrix<Type> > operator+
721     const tmp<fvMatrix<Type> >&,
722     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
725 template<class Type>
726 tmp<fvMatrix<Type> > operator+
728     const DimensionedField<Type, volMesh>&,
729     const fvMatrix<Type>&
732 template<class Type>
733 tmp<fvMatrix<Type> > operator+
735     const tmp<DimensionedField<Type, volMesh> >&,
736     const fvMatrix<Type>&
739 template<class Type>
740 tmp<fvMatrix<Type> > operator+
742     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
743     const fvMatrix<Type>&
746 template<class Type>
747 tmp<fvMatrix<Type> > operator+
749     const DimensionedField<Type, volMesh>&,
750     const tmp<fvMatrix<Type> >&
753 template<class Type>
754 tmp<fvMatrix<Type> > operator+
756     const tmp<DimensionedField<Type, volMesh> >&,
757     const tmp<fvMatrix<Type> >&
760 template<class Type>
761 tmp<fvMatrix<Type> > operator+
763     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
764     const tmp<fvMatrix<Type> >&
768 template<class Type>
769 tmp<fvMatrix<Type> > operator+
771     const fvMatrix<Type>&,
772     const dimensioned<Type>&
775 template<class Type>
776 tmp<fvMatrix<Type> > operator+
778     const tmp<fvMatrix<Type> >&,
779     const dimensioned<Type>&
782 template<class Type>
783 tmp<fvMatrix<Type> > operator+
785     const dimensioned<Type>&,
786     const fvMatrix<Type>&
789 template<class Type>
790 tmp<fvMatrix<Type> > operator+
792     const dimensioned<Type>&,
793     const tmp<fvMatrix<Type> >&
797 template<class Type>
798 tmp<fvMatrix<Type> > operator-
800     const fvMatrix<Type>&,
801     const fvMatrix<Type>&
804 template<class Type>
805 tmp<fvMatrix<Type> > operator-
807     const tmp<fvMatrix<Type> >&,
808     const fvMatrix<Type>&
811 template<class Type>
812 tmp<fvMatrix<Type> > operator-
814     const fvMatrix<Type>&,
815     const tmp<fvMatrix<Type> >&
818 template<class Type>
819 tmp<fvMatrix<Type> > operator-
821     const tmp<fvMatrix<Type> >&,
822     const tmp<fvMatrix<Type> >&
826 template<class Type>
827 tmp<fvMatrix<Type> > operator-
829     const fvMatrix<Type>&,
830     const DimensionedField<Type, volMesh>&
833 template<class Type>
834 tmp<fvMatrix<Type> > operator-
836     const fvMatrix<Type>&,
837     const tmp<DimensionedField<Type, volMesh> >&
840 template<class Type>
841 tmp<fvMatrix<Type> > operator-
843     const fvMatrix<Type>&,
844     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
847 template<class Type>
848 tmp<fvMatrix<Type> > operator-
850     const tmp<fvMatrix<Type> >&,
851     const DimensionedField<Type, volMesh>&
854 template<class Type>
855 tmp<fvMatrix<Type> > operator-
857     const tmp<fvMatrix<Type> >&,
858     const tmp<DimensionedField<Type, volMesh> >&
861 template<class Type>
862 tmp<fvMatrix<Type> > operator-
864     const tmp<fvMatrix<Type> >&,
865     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
868 template<class Type>
869 tmp<fvMatrix<Type> > operator-
871     const DimensionedField<Type, volMesh>&,
872     const fvMatrix<Type>&
875 template<class Type>
876 tmp<fvMatrix<Type> > operator-
878     const tmp<DimensionedField<Type, volMesh> >&,
879     const fvMatrix<Type>&
882 template<class Type>
883 tmp<fvMatrix<Type> > operator-
885     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
886     const fvMatrix<Type>&
889 template<class Type>
890 tmp<fvMatrix<Type> > operator-
892     const DimensionedField<Type, volMesh>&,
893     const tmp<fvMatrix<Type> >&
896 template<class Type>
897 tmp<fvMatrix<Type> > operator-
899     const tmp<DimensionedField<Type, volMesh> >&,
900     const tmp<fvMatrix<Type> >&
903 template<class Type>
904 tmp<fvMatrix<Type> > operator-
906     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
907     const tmp<fvMatrix<Type> >&
911 template<class Type>
912 tmp<fvMatrix<Type> > operator-
914     const fvMatrix<Type>&,
915     const dimensioned<Type>&
918 template<class Type>
919 tmp<fvMatrix<Type> > operator-
921     const tmp<fvMatrix<Type> >&,
922     const dimensioned<Type>&
925 template<class Type>
926 tmp<fvMatrix<Type> > operator-
928     const dimensioned<Type>&,
929     const fvMatrix<Type>&
932 template<class Type>
933 tmp<fvMatrix<Type> > operator-
935     const dimensioned<Type>&,
936     const tmp<fvMatrix<Type> >&
940 template<class Type>
941 tmp<fvMatrix<Type> > operator*
943     const DimensionedField<scalar, volMesh>&,
944     const fvMatrix<Type>&
947 template<class Type>
948 tmp<fvMatrix<Type> > operator*
950     const tmp<DimensionedField<scalar, volMesh> >&,
951     const fvMatrix<Type>&
954 template<class Type>
955 tmp<fvMatrix<Type> > operator*
957     const tmp<volScalarField>&,
958     const fvMatrix<Type>&
961 template<class Type>
962 tmp<fvMatrix<Type> > operator*
964     const DimensionedField<scalar, volMesh>&,
965     const tmp<fvMatrix<Type> >&
968 template<class Type>
969 tmp<fvMatrix<Type> > operator*
971     const tmp<DimensionedField<scalar, volMesh> >&,
972     const tmp<fvMatrix<Type> >&
975 template<class Type>
976 tmp<fvMatrix<Type> > operator*
978     const tmp<volScalarField>&,
979     const tmp<fvMatrix<Type> >&
983 template<class Type>
984 tmp<fvMatrix<Type> > operator*
986     const dimensioned<scalar>&,
987     const fvMatrix<Type>&
990 template<class Type>
991 tmp<fvMatrix<Type> > operator*
993     const dimensioned<scalar>&,
994     const tmp<fvMatrix<Type> >&
998 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1000 } // End namespace Foam
1002 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1004 #ifdef NoRepository
1005 #   include "fvMatrix.C"
1006 #endif
1008 // Specialisation for scalars
1009 #include "fvScalarMatrix.H"
1011 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1013 #endif
1015 // ************************************************************************* //