Go to the documentation of this file.
61 template<
class Type>
class fvMatrix;
62 template<
class T>
class UIndirectList;
65 Ostream&
operator<<(Ostream&,
const fvMatrix<Type>&);
69 tmp<GeometricField<Type, fvPatchField, volMesh>>
operator&
71 const fvMatrix<Type>&,
72 const DimensionedField<Type, volMesh>&
76 tmp<GeometricField<Type, fvPatchField, volMesh>>
operator&
78 const fvMatrix<Type>&,
79 const tmp<DimensionedField<Type, volMesh>>&
83 tmp<GeometricField<Type, fvPatchField, volMesh>>
operator&
85 const fvMatrix<Type>&,
86 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
90 tmp<GeometricField<Type, fvPatchField, volMesh>>
operator&
92 const tmp<fvMatrix<Type>>&,
93 const DimensionedField<Type, volMesh>&
97 tmp<GeometricField<Type, fvPatchField, volMesh>>
operator&
99 const tmp<fvMatrix<Type>>&,
100 const tmp<DimensionedField<Type, volMesh>>&
104 tmp<GeometricField<Type, fvPatchField, volMesh>>
operator&
106 const tmp<fvMatrix<Type>>&,
107 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
125 const GeometricField<Type, fvPatchField, volMesh>& psi_;
128 dimensionSet dimensions_;
135 FieldField<Field, Type> internalCoeffs_;
139 FieldField<Field, Type> boundaryCoeffs_;
142 mutable GeometricField<Type, fvsPatchField, surfaceMesh>
143 *faceFluxCorrectionPtr_;
154 template<
class Type2>
162 template<
class Type2>
171 template<
class Type2>
179 template<
class Type2>
201 const bool couples=
true
208 template<
template<
class>
class ListType>
212 const ListType<Type>&
values
235 solver_(std::move(sol))
312 return internalCoeffs_;
319 return internalCoeffs_;
326 return boundaryCoeffs_;
333 return boundaryCoeffs_;
344 return faceFluxCorrectionPtr_;
379 const bool forceReference =
false
387 const bool forceReference =
false
395 const bool forceReference =
false
565 friend Ostream& operator<< <Type>
578 const fvMatrix<Type>&,
579 const fvMatrix<Type>&,
586 const fvMatrix<Type>&,
587 const DimensionedField<Type, volMesh>&,
594 const fvMatrix<Type>&,
595 const dimensioned<Type>&,
603 SolverPerformance<Type>
solve(fvMatrix<Type>&,
const dictionary&);
610 SolverPerformance<Type>
solve
612 const tmp<fvMatrix<Type>>&,
620 SolverPerformance<Type>
solve(fvMatrix<Type>&);
627 SolverPerformance<Type>
solve(
const tmp<fvMatrix<Type>>&);
633 tmp<fvMatrix<Type>>
correction(
const fvMatrix<Type>&);
639 tmp<fvMatrix<Type>>
correction(
const tmp<fvMatrix<Type>>&);
645 tmp<fvMatrix<Type>>
operator==
647 const fvMatrix<Type>&,
648 const fvMatrix<Type>&
652 tmp<fvMatrix<Type>>
operator==
654 const tmp<fvMatrix<Type>>&,
655 const fvMatrix<Type>&
659 tmp<fvMatrix<Type>>
operator==
661 const fvMatrix<Type>&,
662 const tmp<fvMatrix<Type>>&
666 tmp<fvMatrix<Type>>
operator==
668 const tmp<fvMatrix<Type>>&,
669 const tmp<fvMatrix<Type>>&
674 tmp<fvMatrix<Type>>
operator==
676 const fvMatrix<Type>&,
677 const DimensionedField<Type, volMesh>&
681 tmp<fvMatrix<Type>>
operator==
683 const fvMatrix<Type>&,
684 const tmp<DimensionedField<Type, volMesh>>&
688 tmp<fvMatrix<Type>>
operator==
690 const fvMatrix<Type>&,
691 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
695 tmp<fvMatrix<Type>>
operator==
697 const tmp<fvMatrix<Type>>&,
698 const DimensionedField<Type, volMesh>&
702 tmp<fvMatrix<Type>>
operator==
704 const tmp<fvMatrix<Type>>&,
705 const tmp<DimensionedField<Type, volMesh>>&
709 tmp<fvMatrix<Type>>
operator==
711 const tmp<fvMatrix<Type>>&,
712 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
716 tmp<fvMatrix<Type>>
operator==
718 const fvMatrix<Type>&,
719 const dimensioned<Type>&
723 tmp<fvMatrix<Type>>
operator==
725 const tmp<fvMatrix<Type>>&,
726 const dimensioned<Type>&
731 tmp<fvMatrix<Type>>
operator==
733 const fvMatrix<Type>&,
738 tmp<fvMatrix<Type>>
operator==
740 const tmp<fvMatrix<Type>>&,
746 tmp<fvMatrix<Type>>
operator-
748 const fvMatrix<Type>&
752 tmp<fvMatrix<Type>>
operator-
754 const tmp<fvMatrix<Type>>&
759 tmp<fvMatrix<Type>>
operator+
761 const fvMatrix<Type>&,
762 const fvMatrix<Type>&
766 tmp<fvMatrix<Type>>
operator+
768 const tmp<fvMatrix<Type>>&,
769 const fvMatrix<Type>&
773 tmp<fvMatrix<Type>>
operator+
775 const fvMatrix<Type>&,
776 const tmp<fvMatrix<Type>>&
780 tmp<fvMatrix<Type>>
operator+
782 const tmp<fvMatrix<Type>>&,
783 const tmp<fvMatrix<Type>>&
788 tmp<fvMatrix<Type>>
operator+
790 const fvMatrix<Type>&,
791 const DimensionedField<Type, volMesh>&
795 tmp<fvMatrix<Type>>
operator+
797 const fvMatrix<Type>&,
798 const tmp<DimensionedField<Type, volMesh>>&
802 tmp<fvMatrix<Type>>
operator+
804 const fvMatrix<Type>&,
805 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
809 tmp<fvMatrix<Type>>
operator+
811 const tmp<fvMatrix<Type>>&,
812 const DimensionedField<Type, volMesh>&
816 tmp<fvMatrix<Type>>
operator+
818 const tmp<fvMatrix<Type>>&,
819 const tmp<DimensionedField<Type, volMesh>>&
823 tmp<fvMatrix<Type>>
operator+
825 const tmp<fvMatrix<Type>>&,
826 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
830 tmp<fvMatrix<Type>>
operator+
832 const DimensionedField<Type, volMesh>&,
833 const fvMatrix<Type>&
837 tmp<fvMatrix<Type>>
operator+
839 const tmp<DimensionedField<Type, volMesh>>&,
840 const fvMatrix<Type>&
844 tmp<fvMatrix<Type>>
operator+
846 const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
847 const fvMatrix<Type>&
851 tmp<fvMatrix<Type>>
operator+
853 const DimensionedField<Type, volMesh>&,
854 const tmp<fvMatrix<Type>>&
858 tmp<fvMatrix<Type>>
operator+
860 const tmp<DimensionedField<Type, volMesh>>&,
861 const tmp<fvMatrix<Type>>&
865 tmp<fvMatrix<Type>>
operator+
867 const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
868 const tmp<fvMatrix<Type>>&
873 tmp<fvMatrix<Type>>
operator+
875 const fvMatrix<Type>&,
876 const dimensioned<Type>&
880 tmp<fvMatrix<Type>>
operator+
882 const tmp<fvMatrix<Type>>&,
883 const dimensioned<Type>&
887 tmp<fvMatrix<Type>>
operator+
889 const dimensioned<Type>&,
890 const fvMatrix<Type>&
894 tmp<fvMatrix<Type>>
operator+
896 const dimensioned<Type>&,
897 const tmp<fvMatrix<Type>>&
902 tmp<fvMatrix<Type>>
operator-
904 const fvMatrix<Type>&,
905 const fvMatrix<Type>&
909 tmp<fvMatrix<Type>>
operator-
911 const tmp<fvMatrix<Type>>&,
912 const fvMatrix<Type>&
916 tmp<fvMatrix<Type>>
operator-
918 const fvMatrix<Type>&,
919 const tmp<fvMatrix<Type>>&
923 tmp<fvMatrix<Type>>
operator-
925 const tmp<fvMatrix<Type>>&,
926 const tmp<fvMatrix<Type>>&
931 tmp<fvMatrix<Type>>
operator-
933 const fvMatrix<Type>&,
934 const DimensionedField<Type, volMesh>&
938 tmp<fvMatrix<Type>>
operator-
940 const fvMatrix<Type>&,
941 const tmp<DimensionedField<Type, volMesh>>&
945 tmp<fvMatrix<Type>>
operator-
947 const fvMatrix<Type>&,
948 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
952 tmp<fvMatrix<Type>>
operator-
954 const tmp<fvMatrix<Type>>&,
955 const DimensionedField<Type, volMesh>&
959 tmp<fvMatrix<Type>>
operator-
961 const tmp<fvMatrix<Type>>&,
962 const tmp<DimensionedField<Type, volMesh>>&
966 tmp<fvMatrix<Type>>
operator-
968 const tmp<fvMatrix<Type>>&,
969 const tmp<GeometricField<Type, fvPatchField, volMesh>>&
973 tmp<fvMatrix<Type>>
operator-
975 const DimensionedField<Type, volMesh>&,
976 const fvMatrix<Type>&
980 tmp<fvMatrix<Type>>
operator-
982 const tmp<DimensionedField<Type, volMesh>>&,
983 const fvMatrix<Type>&
987 tmp<fvMatrix<Type>>
operator-
989 const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
990 const fvMatrix<Type>&
994 tmp<fvMatrix<Type>>
operator-
996 const DimensionedField<Type, volMesh>&,
997 const tmp<fvMatrix<Type>>&
1000 template<
class Type>
1001 tmp<fvMatrix<Type>>
operator-
1003 const tmp<DimensionedField<Type, volMesh>>&,
1004 const tmp<fvMatrix<Type>>&
1007 template<
class Type>
1008 tmp<fvMatrix<Type>>
operator-
1010 const tmp<GeometricField<Type, fvPatchField, volMesh>>&,
1011 const tmp<fvMatrix<Type>>&
1015 template<
class Type>
1016 tmp<fvMatrix<Type>>
operator-
1018 const fvMatrix<Type>&,
1019 const dimensioned<Type>&
1022 template<
class Type>
1023 tmp<fvMatrix<Type>>
operator-
1025 const tmp<fvMatrix<Type>>&,
1026 const dimensioned<Type>&
1029 template<
class Type>
1030 tmp<fvMatrix<Type>>
operator-
1032 const dimensioned<Type>&,
1033 const fvMatrix<Type>&
1036 template<
class Type>
1037 tmp<fvMatrix<Type>>
operator-
1039 const dimensioned<Type>&,
1040 const tmp<fvMatrix<Type>>&
1044 template<
class Type>
1045 tmp<fvMatrix<Type>>
operator*
1048 const fvMatrix<Type>&
1051 template<
class Type>
1052 tmp<fvMatrix<Type>>
operator*
1054 const tmp<volScalarField::Internal>&,
1055 const fvMatrix<Type>&
1058 template<
class Type>
1059 tmp<fvMatrix<Type>>
operator*
1061 const tmp<volScalarField>&,
1062 const fvMatrix<Type>&
1065 template<
class Type>
1066 tmp<fvMatrix<Type>>
operator*
1069 const tmp<fvMatrix<Type>>&
1072 template<
class Type>
1073 tmp<fvMatrix<Type>>
operator*
1075 const tmp<volScalarField::Internal>&,
1076 const tmp<fvMatrix<Type>>&
1079 template<
class Type>
1080 tmp<fvMatrix<Type>>
operator*
1082 const tmp<volScalarField>&,
1083 const tmp<fvMatrix<Type>>&
1087 template<
class Type>
1088 tmp<fvMatrix<Type>>
operator*
1090 const dimensioned<scalar>&,
1091 const fvMatrix<Type>&
1094 template<
class Type>
1095 tmp<fvMatrix<Type>>
operator*
1097 const dimensioned<scalar>&,
1098 const tmp<fvMatrix<Type>>&
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
tmp< Field< Type > > residual() const
Return the matrix residual.
SolverPerformance< Type > solve()
Solve returning the solution statistics.
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
void addBoundarySource(Field< Type > &source, const bool couples=true) const
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
void addCmptAvBoundaryDiag(scalarField &diag) const
A field of fields is a PtrList of fields with reference counting.
A class for managing temporary objects.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
tmp< scalarField > D() const
Return the matrix scalar diagonal.
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void operator=(const fvMatrix< Type > &)
SolverPerformance< Type > solve()
Solve returning the solution statistics.
const dimensionSet & dimensions() const
void relax()
Relax matrix (for steady-state solution).
FieldField< Field, Type > & internalCoeffs()
const FieldField< Field, Type > & boundaryCoeffs() const
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
void setValuesFromList(const labelUList &cellLabels, const ListType< Type > &values)
Set solution in given cells to the specified values.
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
fvSolver(fvMatrix< Type > &fvMat, autoPtr< lduMatrix::solver > &&sol)
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Dimension set for the base types.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
surfaceTypeFieldPtr & faceFluxCorrectionPtr()
Return pointer to face-flux non-orthogonal correction field.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
SolverPerformance< Type > solveCoupled(const dictionary &)
Solve coupled returning the solution statistics.
void operator-=(const fvMatrix< Type > &)
Generic templated field type.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void operator*=(const volScalarField::Internal &)
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
tmp< volScalarField > H1() const
Return H(1)
Macro definitions for declaring ClassName(), NamespaceName(), etc.
virtual ~fvMatrix()
Destructor.
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
void setReferences(const labelUList &cellLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A scalar instance of fvMatrix.
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Generic dimensioned Type class.
tmp< volScalarField > A() const
Return the central coefficient.
GeometricField< Type, fvsPatchField, surfaceMesh > * surfaceTypeFieldPtr
Declare return type of the faceFluxCorrectionPtr() function.
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
autoPtr< fvSolver > solver()
Construct and return the solver.
void setValues(const labelUList &cellLabels, const Type &value)
FieldField< Field, Type > & boundaryCoeffs()
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
const Field< Type > & source() const
void operator+=(const fvMatrix< Type > &)
const dictionary & solverDict() const
Return the solver dictionary taking into account finalIteration.
tmp< fvMatrix< Type > > clone() const
Clone.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
A List with indirect addressing.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
const FieldField< Field, Type > & internalCoeffs() const
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...