1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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
36 \*---------------------------------------------------------------------------*/
41 #include "volFields.H"
42 #include "surfaceFields.H"
43 #include "lduMatrix.H"
46 #include "dimensionedTypes.H"
47 #include "zeroField.H"
48 #include "className.H"
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 // Forward declaration of friend functions and operators
61 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
63 const fvMatrix<Type>&,
64 const DimensionedField<Type, volMesh>&
68 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
70 const fvMatrix<Type>&,
71 const tmp<DimensionedField<Type, volMesh> >&
75 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
77 const fvMatrix<Type>&,
78 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
82 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
84 const tmp<fvMatrix<Type> >&,
85 const DimensionedField<Type, volMesh>&
89 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
91 const tmp<fvMatrix<Type> >&,
92 const tmp<DimensionedField<Type, volMesh> >&
96 tmp<GeometricField<Type, fvPatchField, volMesh> > operator&
98 const tmp<fvMatrix<Type> >&,
99 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
103 Ostream& operator<<(Ostream&, const fvMatrix<Type>&);
106 /*---------------------------------------------------------------------------*\
107 Class fvMatrix Declaration
108 \*---------------------------------------------------------------------------*/
120 // Reference to GeometricField<Type, fvPatchField, volMesh>
121 GeometricField<Type, fvPatchField, volMesh>& psi_;
124 dimensionSet dimensions_;
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
149 const unallocLabelList& addr,
150 const Field<Type2>& pf,
154 template<class Type2>
155 void addToInternalField
157 const unallocLabelList& addr,
158 const tmp<Field<Type2> >& tpf,
162 //- Subtract patch contribution from internal field
163 template<class Type2>
164 void subtractFromInternalField
166 const unallocLabelList& addr,
167 const Field<Type2>& pf,
171 template<class Type2>
172 void subtractFromInternalField
174 const unallocLabelList& addr,
175 const tmp<Field<Type2> >& tpf,
180 // Matrix completion functionality
188 void addCmptAvBoundaryDiag(scalarField& diag) const;
190 void addBoundarySource
193 const bool couples=true
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)
205 fvMatrix<Type>& fvMat_;
207 autoPtr<lduMatrix::solver> solver_;
213 fvSolver(fvMatrix<Type>& fvMat, autoPtr<lduMatrix::solver> sol)
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();
232 ClassName("fvMatrix");
237 //- Construct given a field to solve for
240 GeometricField<Type, fvPatchField, volMesh>&,
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> >&);
252 //- Construct from Istream given field to solve for
253 fvMatrix(GeometricField<Type, fvPatchField, volMesh>&, Istream&);
265 const GeometricField<Type, fvPatchField, volMesh>& psi() const
270 GeometricField<Type, fvPatchField, volMesh>& psi()
275 const dimensionSet& dimensions() const
280 Field<Type>& source()
285 const Field<Type>& source() const
290 //- fvBoundary scalar field containing pseudo-matrix coeffs
291 // for internal cells
292 FieldField<Field, Type>& internalCoeffs()
294 return internalCoeffs_;
297 //- fvBoundary scalar field containing pseudo-matrix coeffs
298 // for boundary cells
299 FieldField<Field, Type>& boundaryCoeffs()
301 return boundaryCoeffs_;
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()
312 return faceFluxCorrectionPtr_;
318 //- Set solution in given cells and eliminate corresponding
319 // equations from the matrix
322 const labelList& cells,
323 const Field<Type>& values
326 //- Set reference level for solution
333 //- Set reference level for a component of the solution
334 // on a given patch face
335 void setComponentReference
339 const direction cmpt,
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
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;
385 tmp<volScalarField> H1() const;
387 //- Return the face-flux field from the matrix
388 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
394 void operator=(const fvMatrix<Type>&);
395 void operator=(const tmp<fvMatrix<Type> >&);
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> >&);
407 const DimensionedField<Type, volMesh>&
411 const tmp<DimensionedField<Type, volMesh> >&
415 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
420 const DimensionedField<Type, volMesh>&
424 const tmp<DimensionedField<Type, volMesh> >&
428 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
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>&);
446 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
449 const fvMatrix<Type>&,
450 const DimensionedField<Type, volMesh>&
453 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
456 const fvMatrix<Type>&,
457 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
460 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
463 const tmp<fvMatrix<Type> >&,
464 const DimensionedField<Type, volMesh>&
467 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
470 const tmp<fvMatrix<Type> >&,
471 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
477 friend Ostream& operator<< <Type>
480 const fvMatrix<Type>&
485 // * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * //
490 const fvMatrix<Type>&,
491 const fvMatrix<Type>&,
498 const fvMatrix<Type>&,
499 const DimensionedField<Type, volMesh>&,
506 const fvMatrix<Type>&,
507 const dimensioned<Type>&,
512 //- Solve returning the solution statistics given convergence tolerance
513 // Solver controls read Istream
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
522 lduMatrix::solverPerformance solve(const tmp<fvMatrix<Type> >&, Istream&);
525 //- Solve returning the solution statistics given convergence tolerance
526 // Solver controls read fvSolution
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
535 lduMatrix::solverPerformance solve(const tmp<fvMatrix<Type> >&);
538 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
541 tmp<fvMatrix<Type> > operator==
543 const fvMatrix<Type>&,
544 const fvMatrix<Type>&
548 tmp<fvMatrix<Type> > operator==
550 const tmp<fvMatrix<Type> >&,
551 const fvMatrix<Type>&
555 tmp<fvMatrix<Type> > operator==
557 const fvMatrix<Type>&,
558 const tmp<fvMatrix<Type> >&
562 tmp<fvMatrix<Type> > operator==
564 const tmp<fvMatrix<Type> >&,
565 const tmp<fvMatrix<Type> >&
570 tmp<fvMatrix<Type> > operator==
572 const fvMatrix<Type>&,
573 const DimensionedField<Type, volMesh>&
577 tmp<fvMatrix<Type> > operator==
579 const fvMatrix<Type>&,
580 const tmp<DimensionedField<Type, volMesh> >&
584 tmp<fvMatrix<Type> > operator==
586 const fvMatrix<Type>&,
587 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
591 tmp<fvMatrix<Type> > operator==
593 const tmp<fvMatrix<Type> >&,
594 const DimensionedField<Type, volMesh>&
598 tmp<fvMatrix<Type> > operator==
600 const tmp<fvMatrix<Type> >&,
601 const tmp<DimensionedField<Type, volMesh> >&
605 tmp<fvMatrix<Type> > operator==
607 const tmp<fvMatrix<Type> >&,
608 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
612 tmp<fvMatrix<Type> > operator==
614 const fvMatrix<Type>&,
615 const dimensioned<Type>&
619 tmp<fvMatrix<Type> > operator==
621 const tmp<fvMatrix<Type> >&,
622 const dimensioned<Type>&
627 tmp<fvMatrix<Type> > operator==
629 const fvMatrix<Type>&,
634 tmp<fvMatrix<Type> > operator==
636 const tmp<fvMatrix<Type> >&,
642 tmp<fvMatrix<Type> > operator-
644 const fvMatrix<Type>&
648 tmp<fvMatrix<Type> > operator-
650 const tmp<fvMatrix<Type> >&
655 tmp<fvMatrix<Type> > operator+
657 const fvMatrix<Type>&,
658 const fvMatrix<Type>&
662 tmp<fvMatrix<Type> > operator+
664 const tmp<fvMatrix<Type> >&,
665 const fvMatrix<Type>&
669 tmp<fvMatrix<Type> > operator+
671 const fvMatrix<Type>&,
672 const tmp<fvMatrix<Type> >&
676 tmp<fvMatrix<Type> > operator+
678 const tmp<fvMatrix<Type> >&,
679 const tmp<fvMatrix<Type> >&
684 tmp<fvMatrix<Type> > operator+
686 const fvMatrix<Type>&,
687 const DimensionedField<Type, volMesh>&
691 tmp<fvMatrix<Type> > operator+
693 const fvMatrix<Type>&,
694 const tmp<DimensionedField<Type, volMesh> >&
698 tmp<fvMatrix<Type> > operator+
700 const fvMatrix<Type>&,
701 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
705 tmp<fvMatrix<Type> > operator+
707 const tmp<fvMatrix<Type> >&,
708 const DimensionedField<Type, volMesh>&
712 tmp<fvMatrix<Type> > operator+
714 const tmp<fvMatrix<Type> >&,
715 const tmp<DimensionedField<Type, volMesh> >&
719 tmp<fvMatrix<Type> > operator+
721 const tmp<fvMatrix<Type> >&,
722 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
726 tmp<fvMatrix<Type> > operator+
728 const DimensionedField<Type, volMesh>&,
729 const fvMatrix<Type>&
733 tmp<fvMatrix<Type> > operator+
735 const tmp<DimensionedField<Type, volMesh> >&,
736 const fvMatrix<Type>&
740 tmp<fvMatrix<Type> > operator+
742 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
743 const fvMatrix<Type>&
747 tmp<fvMatrix<Type> > operator+
749 const DimensionedField<Type, volMesh>&,
750 const tmp<fvMatrix<Type> >&
754 tmp<fvMatrix<Type> > operator+
756 const tmp<DimensionedField<Type, volMesh> >&,
757 const tmp<fvMatrix<Type> >&
761 tmp<fvMatrix<Type> > operator+
763 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
764 const tmp<fvMatrix<Type> >&
769 tmp<fvMatrix<Type> > operator+
771 const fvMatrix<Type>&,
772 const dimensioned<Type>&
776 tmp<fvMatrix<Type> > operator+
778 const tmp<fvMatrix<Type> >&,
779 const dimensioned<Type>&
783 tmp<fvMatrix<Type> > operator+
785 const dimensioned<Type>&,
786 const fvMatrix<Type>&
790 tmp<fvMatrix<Type> > operator+
792 const dimensioned<Type>&,
793 const tmp<fvMatrix<Type> >&
798 tmp<fvMatrix<Type> > operator-
800 const fvMatrix<Type>&,
801 const fvMatrix<Type>&
805 tmp<fvMatrix<Type> > operator-
807 const tmp<fvMatrix<Type> >&,
808 const fvMatrix<Type>&
812 tmp<fvMatrix<Type> > operator-
814 const fvMatrix<Type>&,
815 const tmp<fvMatrix<Type> >&
819 tmp<fvMatrix<Type> > operator-
821 const tmp<fvMatrix<Type> >&,
822 const tmp<fvMatrix<Type> >&
827 tmp<fvMatrix<Type> > operator-
829 const fvMatrix<Type>&,
830 const DimensionedField<Type, volMesh>&
834 tmp<fvMatrix<Type> > operator-
836 const fvMatrix<Type>&,
837 const tmp<DimensionedField<Type, volMesh> >&
841 tmp<fvMatrix<Type> > operator-
843 const fvMatrix<Type>&,
844 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
848 tmp<fvMatrix<Type> > operator-
850 const tmp<fvMatrix<Type> >&,
851 const DimensionedField<Type, volMesh>&
855 tmp<fvMatrix<Type> > operator-
857 const tmp<fvMatrix<Type> >&,
858 const tmp<DimensionedField<Type, volMesh> >&
862 tmp<fvMatrix<Type> > operator-
864 const tmp<fvMatrix<Type> >&,
865 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
869 tmp<fvMatrix<Type> > operator-
871 const DimensionedField<Type, volMesh>&,
872 const fvMatrix<Type>&
876 tmp<fvMatrix<Type> > operator-
878 const tmp<DimensionedField<Type, volMesh> >&,
879 const fvMatrix<Type>&
883 tmp<fvMatrix<Type> > operator-
885 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
886 const fvMatrix<Type>&
890 tmp<fvMatrix<Type> > operator-
892 const DimensionedField<Type, volMesh>&,
893 const tmp<fvMatrix<Type> >&
897 tmp<fvMatrix<Type> > operator-
899 const tmp<DimensionedField<Type, volMesh> >&,
900 const tmp<fvMatrix<Type> >&
904 tmp<fvMatrix<Type> > operator-
906 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
907 const tmp<fvMatrix<Type> >&
912 tmp<fvMatrix<Type> > operator-
914 const fvMatrix<Type>&,
915 const dimensioned<Type>&
919 tmp<fvMatrix<Type> > operator-
921 const tmp<fvMatrix<Type> >&,
922 const dimensioned<Type>&
926 tmp<fvMatrix<Type> > operator-
928 const dimensioned<Type>&,
929 const fvMatrix<Type>&
933 tmp<fvMatrix<Type> > operator-
935 const dimensioned<Type>&,
936 const tmp<fvMatrix<Type> >&
941 tmp<fvMatrix<Type> > operator*
943 const DimensionedField<scalar, volMesh>&,
944 const fvMatrix<Type>&
948 tmp<fvMatrix<Type> > operator*
950 const tmp<DimensionedField<scalar, volMesh> >&,
951 const fvMatrix<Type>&
955 tmp<fvMatrix<Type> > operator*
957 const tmp<volScalarField>&,
958 const fvMatrix<Type>&
962 tmp<fvMatrix<Type> > operator*
964 const DimensionedField<scalar, volMesh>&,
965 const tmp<fvMatrix<Type> >&
969 tmp<fvMatrix<Type> > operator*
971 const tmp<DimensionedField<scalar, volMesh> >&,
972 const tmp<fvMatrix<Type> >&
976 tmp<fvMatrix<Type> > operator*
978 const tmp<volScalarField>&,
979 const tmp<fvMatrix<Type> >&
984 tmp<fvMatrix<Type> > operator*
986 const dimensioned<scalar>&,
987 const fvMatrix<Type>&
991 tmp<fvMatrix<Type> > operator*
993 const dimensioned<scalar>&,
994 const tmp<fvMatrix<Type> >&
998 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1000 } // End namespace Foam
1002 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1005 # include "fvMatrix.C"
1008 // Specialisation for scalars
1009 #include "fvScalarMatrix.H"
1011 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1015 // ************************************************************************* //