ddtScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fv.H"
27 #include "HashTable.H"
28 #include "surfaceInterpolate.H"
29 #include "fvMatrix.H"
30 #include "cyclicAMIFvPatch.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fv
40 {
41 
42 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
43 
44 template<class Type>
45 tmp<ddtScheme<Type>> ddtScheme<Type>::New
46 (
47  const fvMesh& mesh,
48  Istream& schemeData
49 )
50 {
51  if (fv::debug)
52  {
53  InfoInFunction << "Constructing ddtScheme<Type>" << endl;
54  }
55 
56  if (schemeData.eof())
57  {
59  (
60  schemeData
61  ) << "Ddt scheme not specified" << endl << endl
62  << "Valid ddt schemes are :" << endl
63  << IstreamConstructorTablePtr_->sortedToc()
64  << exit(FatalIOError);
65  }
66 
67  const word schemeName(schemeData);
68 
69  typename IstreamConstructorTable::iterator cstrIter =
70  IstreamConstructorTablePtr_->find(schemeName);
71 
72  if (cstrIter == IstreamConstructorTablePtr_->end())
73  {
75  (
76  schemeData
77  ) << "Unknown ddt scheme " << schemeName << nl << nl
78  << "Valid ddt schemes are :" << endl
79  << IstreamConstructorTablePtr_->sortedToc()
80  << exit(FatalIOError);
81  }
82 
83  return cstrIter()(mesh, schemeData);
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
89 template<class Type>
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95 
96 template<class Type>
98 (
99  const volScalarField& alpha,
100  const volScalarField& rho,
102 )
103 {
105 
107  (
109  );
110 }
111 
112 
113 template<class Type>
115 (
116  const volScalarField& alpha,
117  const volScalarField& rho,
119 )
120 {
122 
123  return tmp<fvMatrix<Type>>
124  (
125  new fvMatrix<Type>
126  (
127  vf,
128  alpha.dimensions()*rho.dimensions()
129  *vf.dimensions()*dimVol/dimTime
130  )
131  );
132 }
133 
134 
135 template<class Type>
137 (
139 )
140 {
142 
144  (
146  );
147 }
148 
149 
150 
151 template<class Type>
153 (
155  const fluxFieldType& phi,
156  const fluxFieldType& phiCorr
157 )
158 {
159  // Courant number limited formulation
160  /*
161  tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
162  - min
163  (
164  mag(phiCorr)
165  *mesh().time().deltaT()*mesh().deltaCoeffs()/mesh().magSf(),
166  scalar(1)
167  );
168  */
169 
170  // Flux normalised formulation
171  tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
172  - min
173  (
174  mag(phiCorr)
175  /(mag(phi) + dimensionedScalar("small", phi.dimensions(), SMALL)),
176  scalar(1)
177  );
178 
179  surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff.ref();
180 
181  surfaceScalarField::Boundary& ccbf = ddtCouplingCoeff.boundaryFieldRef();
182 
184  {
185  if
186  (
187  !U.boundaryField()[patchi].coupled()
188  || isA<cyclicAMIFvPatch>(mesh().boundary()[patchi])
189  )
190  {
191  ccbf[patchi] = 0.0;
192  }
193  }
194 
195  if (debug > 1)
196  {
198  << "ddtCouplingCoeff mean max min = "
199  << gAverage(ddtCouplingCoeff.primitiveField())
200  << " " << gMax(ddtCouplingCoeff.primitiveField())
201  << " " << gMin(ddtCouplingCoeff.primitiveField())
202  << endl;
203  }
204 
205  return tddtCouplingCoeff;
206 }
207 
208 
209 template<class Type>
211 (
213  const fluxFieldType& phi,
214  const fluxFieldType& phiCorr,
215  const volScalarField& rho
216 )
217 {
218  return fvcDdtPhiCoeff(U, phi, phiCorr);
219 }
220 
221 
222 template<class Type>
224 (
226  const fluxFieldType& phi
227 )
228 {
229  return fvcDdtPhiCoeff(U, phi, phi - fvc::dotInterpolate(mesh().Sf(), U));
230 }
231 
232 
233 template<class Type>
235 (
237  const fluxFieldType& phi,
238  const volScalarField& rho
239 )
240 {
241  return fvcDdtPhiCoeff
242  (
243  U,
244  phi,
245  (phi - fvc::dotInterpolate(mesh().Sf(), U))
246  );
247 }
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace fv
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 } // End namespace Foam
257 
258 // ************************************************************************* //
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gMin(const FieldField< Field, Type > &f)
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static const scalar SMALL
Definition: scalar.H:115
fvMesh & mesh
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
virtual ~ddtScheme()
Destructor.
Definition: ddtScheme.C:90
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
const dimensionSet dimVol
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:336
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
virtual tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:153
label patchi
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:46
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
Namespace for OpenFOAM.
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.