1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 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
331 const bool forceReference = false
334 //- Set reference level for a component of the solution
335 // on a given patch face
336 void setComponentReference
340 const direction cmpt,
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
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;
386 tmp<volScalarField> H1() const;
388 //- Return the face-flux field from the matrix
389 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
395 void operator=(const fvMatrix<Type>&);
396 void operator=(const tmp<fvMatrix<Type> >&);
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> >&);
408 const DimensionedField<Type, volMesh>&
412 const tmp<DimensionedField<Type, volMesh> >&
416 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
421 const DimensionedField<Type, volMesh>&
425 const tmp<DimensionedField<Type, volMesh> >&
429 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
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>&);
447 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
450 const fvMatrix<Type>&,
451 const DimensionedField<Type, volMesh>&
454 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
457 const fvMatrix<Type>&,
458 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
461 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
464 const tmp<fvMatrix<Type> >&,
465 const DimensionedField<Type, volMesh>&
468 friend tmp<GeometricField<Type, fvPatchField, volMesh> >
471 const tmp<fvMatrix<Type> >&,
472 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
478 friend Ostream& operator<< <Type>
481 const fvMatrix<Type>&
486 // * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * //
491 const fvMatrix<Type>&,
492 const fvMatrix<Type>&,
499 const fvMatrix<Type>&,
500 const DimensionedField<Type, volMesh>&,
507 const fvMatrix<Type>&,
508 const dimensioned<Type>&,
513 //- Solve returning the solution statistics given convergence tolerance
514 // Solver controls read Istream
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
523 lduMatrix::solverPerformance solve(const tmp<fvMatrix<Type> >&, Istream&);
526 //- Solve returning the solution statistics given convergence tolerance
527 // Solver controls read fvSolution
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
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
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
548 tmp<fvMatrix<Type> > correction(const tmp<fvMatrix<Type> >&);
551 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
554 tmp<fvMatrix<Type> > operator==
556 const fvMatrix<Type>&,
557 const fvMatrix<Type>&
561 tmp<fvMatrix<Type> > operator==
563 const tmp<fvMatrix<Type> >&,
564 const fvMatrix<Type>&
568 tmp<fvMatrix<Type> > operator==
570 const fvMatrix<Type>&,
571 const tmp<fvMatrix<Type> >&
575 tmp<fvMatrix<Type> > operator==
577 const tmp<fvMatrix<Type> >&,
578 const tmp<fvMatrix<Type> >&
583 tmp<fvMatrix<Type> > operator==
585 const fvMatrix<Type>&,
586 const DimensionedField<Type, volMesh>&
590 tmp<fvMatrix<Type> > operator==
592 const fvMatrix<Type>&,
593 const tmp<DimensionedField<Type, volMesh> >&
597 tmp<fvMatrix<Type> > operator==
599 const fvMatrix<Type>&,
600 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
604 tmp<fvMatrix<Type> > operator==
606 const tmp<fvMatrix<Type> >&,
607 const DimensionedField<Type, volMesh>&
611 tmp<fvMatrix<Type> > operator==
613 const tmp<fvMatrix<Type> >&,
614 const tmp<DimensionedField<Type, volMesh> >&
618 tmp<fvMatrix<Type> > operator==
620 const tmp<fvMatrix<Type> >&,
621 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
625 tmp<fvMatrix<Type> > operator==
627 const fvMatrix<Type>&,
628 const dimensioned<Type>&
632 tmp<fvMatrix<Type> > operator==
634 const tmp<fvMatrix<Type> >&,
635 const dimensioned<Type>&
640 tmp<fvMatrix<Type> > operator==
642 const fvMatrix<Type>&,
647 tmp<fvMatrix<Type> > operator==
649 const tmp<fvMatrix<Type> >&,
655 tmp<fvMatrix<Type> > operator-
657 const fvMatrix<Type>&
661 tmp<fvMatrix<Type> > operator-
663 const tmp<fvMatrix<Type> >&
668 tmp<fvMatrix<Type> > operator+
670 const fvMatrix<Type>&,
671 const fvMatrix<Type>&
675 tmp<fvMatrix<Type> > operator+
677 const tmp<fvMatrix<Type> >&,
678 const fvMatrix<Type>&
682 tmp<fvMatrix<Type> > operator+
684 const fvMatrix<Type>&,
685 const tmp<fvMatrix<Type> >&
689 tmp<fvMatrix<Type> > operator+
691 const tmp<fvMatrix<Type> >&,
692 const tmp<fvMatrix<Type> >&
697 tmp<fvMatrix<Type> > operator+
699 const fvMatrix<Type>&,
700 const DimensionedField<Type, volMesh>&
704 tmp<fvMatrix<Type> > operator+
706 const fvMatrix<Type>&,
707 const tmp<DimensionedField<Type, volMesh> >&
711 tmp<fvMatrix<Type> > operator+
713 const fvMatrix<Type>&,
714 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
718 tmp<fvMatrix<Type> > operator+
720 const tmp<fvMatrix<Type> >&,
721 const DimensionedField<Type, volMesh>&
725 tmp<fvMatrix<Type> > operator+
727 const tmp<fvMatrix<Type> >&,
728 const tmp<DimensionedField<Type, volMesh> >&
732 tmp<fvMatrix<Type> > operator+
734 const tmp<fvMatrix<Type> >&,
735 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
739 tmp<fvMatrix<Type> > operator+
741 const DimensionedField<Type, volMesh>&,
742 const fvMatrix<Type>&
746 tmp<fvMatrix<Type> > operator+
748 const tmp<DimensionedField<Type, volMesh> >&,
749 const fvMatrix<Type>&
753 tmp<fvMatrix<Type> > operator+
755 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
756 const fvMatrix<Type>&
760 tmp<fvMatrix<Type> > operator+
762 const DimensionedField<Type, volMesh>&,
763 const tmp<fvMatrix<Type> >&
767 tmp<fvMatrix<Type> > operator+
769 const tmp<DimensionedField<Type, volMesh> >&,
770 const tmp<fvMatrix<Type> >&
774 tmp<fvMatrix<Type> > operator+
776 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
777 const tmp<fvMatrix<Type> >&
782 tmp<fvMatrix<Type> > operator+
784 const fvMatrix<Type>&,
785 const dimensioned<Type>&
789 tmp<fvMatrix<Type> > operator+
791 const tmp<fvMatrix<Type> >&,
792 const dimensioned<Type>&
796 tmp<fvMatrix<Type> > operator+
798 const dimensioned<Type>&,
799 const fvMatrix<Type>&
803 tmp<fvMatrix<Type> > operator+
805 const dimensioned<Type>&,
806 const tmp<fvMatrix<Type> >&
811 tmp<fvMatrix<Type> > operator-
813 const fvMatrix<Type>&,
814 const fvMatrix<Type>&
818 tmp<fvMatrix<Type> > operator-
820 const tmp<fvMatrix<Type> >&,
821 const fvMatrix<Type>&
825 tmp<fvMatrix<Type> > operator-
827 const fvMatrix<Type>&,
828 const tmp<fvMatrix<Type> >&
832 tmp<fvMatrix<Type> > operator-
834 const tmp<fvMatrix<Type> >&,
835 const tmp<fvMatrix<Type> >&
840 tmp<fvMatrix<Type> > operator-
842 const fvMatrix<Type>&,
843 const DimensionedField<Type, volMesh>&
847 tmp<fvMatrix<Type> > operator-
849 const fvMatrix<Type>&,
850 const tmp<DimensionedField<Type, volMesh> >&
854 tmp<fvMatrix<Type> > operator-
856 const fvMatrix<Type>&,
857 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
861 tmp<fvMatrix<Type> > operator-
863 const tmp<fvMatrix<Type> >&,
864 const DimensionedField<Type, volMesh>&
868 tmp<fvMatrix<Type> > operator-
870 const tmp<fvMatrix<Type> >&,
871 const tmp<DimensionedField<Type, volMesh> >&
875 tmp<fvMatrix<Type> > operator-
877 const tmp<fvMatrix<Type> >&,
878 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
882 tmp<fvMatrix<Type> > operator-
884 const DimensionedField<Type, volMesh>&,
885 const fvMatrix<Type>&
889 tmp<fvMatrix<Type> > operator-
891 const tmp<DimensionedField<Type, volMesh> >&,
892 const fvMatrix<Type>&
896 tmp<fvMatrix<Type> > operator-
898 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
899 const fvMatrix<Type>&
903 tmp<fvMatrix<Type> > operator-
905 const DimensionedField<Type, volMesh>&,
906 const tmp<fvMatrix<Type> >&
910 tmp<fvMatrix<Type> > operator-
912 const tmp<DimensionedField<Type, volMesh> >&,
913 const tmp<fvMatrix<Type> >&
917 tmp<fvMatrix<Type> > operator-
919 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
920 const tmp<fvMatrix<Type> >&
925 tmp<fvMatrix<Type> > operator-
927 const fvMatrix<Type>&,
928 const dimensioned<Type>&
932 tmp<fvMatrix<Type> > operator-
934 const tmp<fvMatrix<Type> >&,
935 const dimensioned<Type>&
939 tmp<fvMatrix<Type> > operator-
941 const dimensioned<Type>&,
942 const fvMatrix<Type>&
946 tmp<fvMatrix<Type> > operator-
948 const dimensioned<Type>&,
949 const tmp<fvMatrix<Type> >&
954 tmp<fvMatrix<Type> > operator*
956 const DimensionedField<scalar, volMesh>&,
957 const fvMatrix<Type>&
961 tmp<fvMatrix<Type> > operator*
963 const tmp<DimensionedField<scalar, volMesh> >&,
964 const fvMatrix<Type>&
968 tmp<fvMatrix<Type> > operator*
970 const tmp<volScalarField>&,
971 const fvMatrix<Type>&
975 tmp<fvMatrix<Type> > operator*
977 const DimensionedField<scalar, volMesh>&,
978 const tmp<fvMatrix<Type> >&
982 tmp<fvMatrix<Type> > operator*
984 const tmp<DimensionedField<scalar, volMesh> >&,
985 const tmp<fvMatrix<Type> >&
989 tmp<fvMatrix<Type> > operator*
991 const tmp<volScalarField>&,
992 const tmp<fvMatrix<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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1018 # include "fvMatrix.C"
1021 // Specialisation for scalars
1022 #include "fvScalarMatrix.H"
1024 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1028 // ************************************************************************* //