3171becb1b958e73911ba29b549421cdb972948f
[OpenFOAM-1.5.x.git] / src / finiteVolume / fvMatrices / fvMatrix / fvMatrix.H
blob3171becb1b958e73911ba29b549421cdb972948f
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 // * * * * * * * * * * * * * * * Global operators  * * * * * * * * * * * * * //
541 template<class Type>
542 tmp<fvMatrix<Type> > operator==
544     const fvMatrix<Type>&,
545     const fvMatrix<Type>&
548 template<class Type>
549 tmp<fvMatrix<Type> > operator==
551     const tmp<fvMatrix<Type> >&,
552     const fvMatrix<Type>&
555 template<class Type>
556 tmp<fvMatrix<Type> > operator==
558     const fvMatrix<Type>&,
559     const tmp<fvMatrix<Type> >&
562 template<class Type>
563 tmp<fvMatrix<Type> > operator==
565     const tmp<fvMatrix<Type> >&,
566     const tmp<fvMatrix<Type> >&
570 template<class Type>
571 tmp<fvMatrix<Type> > operator==
573     const fvMatrix<Type>&,
574     const DimensionedField<Type, volMesh>&
577 template<class Type>
578 tmp<fvMatrix<Type> > operator==
580     const fvMatrix<Type>&,
581     const tmp<DimensionedField<Type, volMesh> >&
584 template<class Type>
585 tmp<fvMatrix<Type> > operator==
587     const fvMatrix<Type>&,
588     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
591 template<class Type>
592 tmp<fvMatrix<Type> > operator==
594     const tmp<fvMatrix<Type> >&,
595     const DimensionedField<Type, volMesh>&
598 template<class Type>
599 tmp<fvMatrix<Type> > operator==
601     const tmp<fvMatrix<Type> >&,
602     const tmp<DimensionedField<Type, volMesh> >&
605 template<class Type>
606 tmp<fvMatrix<Type> > operator==
608     const tmp<fvMatrix<Type> >&,
609     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
612 template<class Type>
613 tmp<fvMatrix<Type> > operator==
615     const fvMatrix<Type>&,
616     const dimensioned<Type>&
619 template<class Type>
620 tmp<fvMatrix<Type> > operator==
622     const tmp<fvMatrix<Type> >&,
623     const dimensioned<Type>&
627 template<class Type>
628 tmp<fvMatrix<Type> > operator==
630     const fvMatrix<Type>&,
631     const zeroField&
634 template<class Type>
635 tmp<fvMatrix<Type> > operator==
637     const tmp<fvMatrix<Type> >&,
638     const zeroField&
642 template<class Type>
643 tmp<fvMatrix<Type> > operator-
645     const fvMatrix<Type>&
648 template<class Type>
649 tmp<fvMatrix<Type> > operator-
651     const tmp<fvMatrix<Type> >&
655 template<class Type>
656 tmp<fvMatrix<Type> > operator+
658     const fvMatrix<Type>&,
659     const fvMatrix<Type>&
662 template<class Type>
663 tmp<fvMatrix<Type> > operator+
665     const tmp<fvMatrix<Type> >&,
666     const fvMatrix<Type>&
669 template<class Type>
670 tmp<fvMatrix<Type> > operator+
672     const fvMatrix<Type>&,
673     const tmp<fvMatrix<Type> >&
676 template<class Type>
677 tmp<fvMatrix<Type> > operator+
679     const tmp<fvMatrix<Type> >&,
680     const tmp<fvMatrix<Type> >&
684 template<class Type>
685 tmp<fvMatrix<Type> > operator+
687     const fvMatrix<Type>&,
688     const DimensionedField<Type, volMesh>&
691 template<class Type>
692 tmp<fvMatrix<Type> > operator+
694     const fvMatrix<Type>&,
695     const tmp<DimensionedField<Type, volMesh> >&
698 template<class Type>
699 tmp<fvMatrix<Type> > operator+
701     const fvMatrix<Type>&,
702     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
705 template<class Type>
706 tmp<fvMatrix<Type> > operator+
708     const tmp<fvMatrix<Type> >&,
709     const DimensionedField<Type, volMesh>&
712 template<class Type>
713 tmp<fvMatrix<Type> > operator+
715     const tmp<fvMatrix<Type> >&,
716     const tmp<DimensionedField<Type, volMesh> >&
719 template<class Type>
720 tmp<fvMatrix<Type> > operator+
722     const tmp<fvMatrix<Type> >&,
723     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
726 template<class Type>
727 tmp<fvMatrix<Type> > operator+
729     const DimensionedField<Type, volMesh>&,
730     const fvMatrix<Type>&
733 template<class Type>
734 tmp<fvMatrix<Type> > operator+
736     const tmp<DimensionedField<Type, volMesh> >&,
737     const fvMatrix<Type>&
740 template<class Type>
741 tmp<fvMatrix<Type> > operator+
743     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
744     const fvMatrix<Type>&
747 template<class Type>
748 tmp<fvMatrix<Type> > operator+
750     const DimensionedField<Type, volMesh>&,
751     const tmp<fvMatrix<Type> >&
754 template<class Type>
755 tmp<fvMatrix<Type> > operator+
757     const tmp<DimensionedField<Type, volMesh> >&,
758     const tmp<fvMatrix<Type> >&
761 template<class Type>
762 tmp<fvMatrix<Type> > operator+
764     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
765     const tmp<fvMatrix<Type> >&
769 template<class Type>
770 tmp<fvMatrix<Type> > operator+
772     const fvMatrix<Type>&,
773     const dimensioned<Type>&
776 template<class Type>
777 tmp<fvMatrix<Type> > operator+
779     const tmp<fvMatrix<Type> >&,
780     const dimensioned<Type>&
783 template<class Type>
784 tmp<fvMatrix<Type> > operator+
786     const dimensioned<Type>&,
787     const fvMatrix<Type>&
790 template<class Type>
791 tmp<fvMatrix<Type> > operator+
793     const dimensioned<Type>&,
794     const tmp<fvMatrix<Type> >&
798 template<class Type>
799 tmp<fvMatrix<Type> > operator-
801     const fvMatrix<Type>&,
802     const fvMatrix<Type>&
805 template<class Type>
806 tmp<fvMatrix<Type> > operator-
808     const tmp<fvMatrix<Type> >&,
809     const fvMatrix<Type>&
812 template<class Type>
813 tmp<fvMatrix<Type> > operator-
815     const fvMatrix<Type>&,
816     const tmp<fvMatrix<Type> >&
819 template<class Type>
820 tmp<fvMatrix<Type> > operator-
822     const tmp<fvMatrix<Type> >&,
823     const tmp<fvMatrix<Type> >&
827 template<class Type>
828 tmp<fvMatrix<Type> > operator-
830     const fvMatrix<Type>&,
831     const DimensionedField<Type, volMesh>&
834 template<class Type>
835 tmp<fvMatrix<Type> > operator-
837     const fvMatrix<Type>&,
838     const tmp<DimensionedField<Type, volMesh> >&
841 template<class Type>
842 tmp<fvMatrix<Type> > operator-
844     const fvMatrix<Type>&,
845     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
848 template<class Type>
849 tmp<fvMatrix<Type> > operator-
851     const tmp<fvMatrix<Type> >&,
852     const DimensionedField<Type, volMesh>&
855 template<class Type>
856 tmp<fvMatrix<Type> > operator-
858     const tmp<fvMatrix<Type> >&,
859     const tmp<DimensionedField<Type, volMesh> >&
862 template<class Type>
863 tmp<fvMatrix<Type> > operator-
865     const tmp<fvMatrix<Type> >&,
866     const tmp<GeometricField<Type, fvPatchField, volMesh> >&
869 template<class Type>
870 tmp<fvMatrix<Type> > operator-
872     const DimensionedField<Type, volMesh>&,
873     const fvMatrix<Type>&
876 template<class Type>
877 tmp<fvMatrix<Type> > operator-
879     const tmp<DimensionedField<Type, volMesh> >&,
880     const fvMatrix<Type>&
883 template<class Type>
884 tmp<fvMatrix<Type> > operator-
886     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
887     const fvMatrix<Type>&
890 template<class Type>
891 tmp<fvMatrix<Type> > operator-
893     const DimensionedField<Type, volMesh>&,
894     const tmp<fvMatrix<Type> >&
897 template<class Type>
898 tmp<fvMatrix<Type> > operator-
900     const tmp<DimensionedField<Type, volMesh> >&,
901     const tmp<fvMatrix<Type> >&
904 template<class Type>
905 tmp<fvMatrix<Type> > operator-
907     const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
908     const tmp<fvMatrix<Type> >&
912 template<class Type>
913 tmp<fvMatrix<Type> > operator-
915     const fvMatrix<Type>&,
916     const dimensioned<Type>&
919 template<class Type>
920 tmp<fvMatrix<Type> > operator-
922     const tmp<fvMatrix<Type> >&,
923     const dimensioned<Type>&
926 template<class Type>
927 tmp<fvMatrix<Type> > operator-
929     const dimensioned<Type>&,
930     const fvMatrix<Type>&
933 template<class Type>
934 tmp<fvMatrix<Type> > operator-
936     const dimensioned<Type>&,
937     const tmp<fvMatrix<Type> >&
941 template<class Type>
942 tmp<fvMatrix<Type> > operator*
944     const DimensionedField<scalar, volMesh>&,
945     const fvMatrix<Type>&
948 template<class Type>
949 tmp<fvMatrix<Type> > operator*
951     const tmp<DimensionedField<scalar, volMesh> >&,
952     const fvMatrix<Type>&
955 template<class Type>
956 tmp<fvMatrix<Type> > operator*
958     const tmp<volScalarField>&,
959     const fvMatrix<Type>&
962 template<class Type>
963 tmp<fvMatrix<Type> > operator*
965     const DimensionedField<scalar, volMesh>&,
966     const tmp<fvMatrix<Type> >&
969 template<class Type>
970 tmp<fvMatrix<Type> > operator*
972     const tmp<DimensionedField<scalar, volMesh> >&,
973     const tmp<fvMatrix<Type> >&
976 template<class Type>
977 tmp<fvMatrix<Type> > operator*
979     const tmp<volScalarField>&,
980     const tmp<fvMatrix<Type> >&
984 template<class Type>
985 tmp<fvMatrix<Type> > operator*
987     const dimensioned<scalar>&,
988     const fvMatrix<Type>&
991 template<class Type>
992 tmp<fvMatrix<Type> > operator*
994     const dimensioned<scalar>&,
995     const tmp<fvMatrix<Type> >&
999 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1001 } // End namespace Foam
1003 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1005 #ifdef NoRepository
1006 #   include "fvMatrix.C"
1007 #endif
1009 // Specialisation for scalars
1010 #include "fvScalarMatrix.H"
1012 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1014 #endif
1016 // ************************************************************************* //