ddtScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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  Info<< "ddtScheme<Type>::New(const fvMesh&, Istream&) : "
54  "constructing ddtScheme<Type>"
55  << endl;
56  }
57 
58  if (schemeData.eof())
59  {
61  (
62  "ddtScheme<Type>::New(const fvMesh&, Istream&)",
63  schemeData
64  ) << "Ddt scheme not specified" << endl << endl
65  << "Valid ddt schemes are :" << endl
66  << IstreamConstructorTablePtr_->sortedToc()
67  << exit(FatalIOError);
68  }
69 
70  const word schemeName(schemeData);
71 
72  typename IstreamConstructorTable::iterator cstrIter =
73  IstreamConstructorTablePtr_->find(schemeName);
74 
75  if (cstrIter == IstreamConstructorTablePtr_->end())
76  {
78  (
79  "ddtScheme<Type>::New(const fvMesh&, Istream&)",
80  schemeData
81  ) << "Unknown ddt scheme " << schemeName << nl << nl
82  << "Valid ddt schemes are :" << endl
83  << IstreamConstructorTablePtr_->sortedToc()
84  << exit(FatalIOError);
85  }
86 
87  return cstrIter()(mesh, schemeData);
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
93 template<class Type>
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100 template<class Type>
102 (
103  const volScalarField& alpha,
104  const volScalarField& rho,
106 )
107 {
108  notImplemented("fvmDdt(alpha, rho, psi");
109 
111  (
113  );
114 }
115 
116 
117 template<class Type>
119 (
120  const volScalarField& alpha,
121  const volScalarField& rho,
123 )
124 {
125  notImplemented("fvmDdt(alpha, rho, psi");
126 
127  return tmp<fvMatrix<Type> >
128  (
129  new fvMatrix<Type>
130  (
131  vf,
132  alpha.dimensions()*rho.dimensions()
133  *vf.dimensions()*dimVol/dimTime
134  )
135  );
136 }
137 
138 
139 template<class Type>
141 (
143  const fluxFieldType& phi,
144  const fluxFieldType& phiCorr
145 )
146 {
147  tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
148  - min
149  (
150  mag(phiCorr)
151  /(mag(phi) + dimensionedScalar("small", phi.dimensions(), SMALL)),
152  scalar(1)
153  );
154 
155  surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
156 
158  {
159  if
160  (
161  U.boundaryField()[patchi].fixesValue()
162  || isA<cyclicAMIFvPatch>(mesh().boundary()[patchi])
163  )
164  {
165  ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
166  }
167  }
168 
169  if (debug > 1)
170  {
171  Info<< "ddtCouplingCoeff mean max min = "
172  << gAverage(ddtCouplingCoeff.internalField())
173  << " " << gMax(ddtCouplingCoeff.internalField())
174  << " " << gMin(ddtCouplingCoeff.internalField())
175  << endl;
176  }
177 
178  return tddtCouplingCoeff;
179 }
180 
181 
182 template<class Type>
184 (
186  const fluxFieldType& phi
187 )
188 {
189  return fvcDdtPhiCoeff(U, phi, phi - (mesh().Sf() & fvc::interpolate(U)));
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 } // End namespace fv
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace Foam
200 
201 // ************************************************************************* //
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
virtual ~ddtScheme()
Destructor.
Definition: ddtScheme.C:94
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensioned< scalar > mag(const dimensioned< Type > &)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
A class for handling words, derived from string.
Definition: word.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Surface Interpolation.
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:68
messageStream Info
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:46
dynamicFvMesh & mesh
Namespace for OpenFOAM.
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
static const char nl
Definition: Ostream.H:260
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
#define forAll(list, i)
Definition: UList.H:421
label patchi
const dimensionSet & dimensions() const
Return dimensions.
Type gMin(const FieldField< Field, Type > &f)
surfaceScalarField & phi
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Type gAverage(const FieldField< Field, Type > &f)
labelList fv(nPoints)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:356
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:141
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
U
Definition: pEqn.H:82