oversetFvPatchField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "volFields.H"
29 #include "cellCellStencil.H"
30 #include "cellCellStencilObject.H"
31 #include "dynamicOversetFvMesh.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvPatch& p,
40 )
41 :
43  oversetPatch_(refCast<const oversetFvPatch>(p))
44 {}
45 
46 
47 template<class Type>
49 (
50  const oversetFvPatchField<Type>& ptf,
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  zeroGradientFvPatchField<Type>(ptf, p, iF, mapper),
57  oversetPatch_(refCast<const oversetFvPatch>(p))
58 {}
59 
60 
61 template<class Type>
63 (
64  const fvPatch& p,
66  const dictionary& dict
67 )
68 :
70  oversetPatch_(refCast<const oversetFvPatch>(p, dict))
71 {}
72 
73 
74 template<class Type>
76 (
77  const oversetFvPatchField<Type>& ptf
78 )
79 :
81  oversetPatch_(ptf.oversetPatch_)
82 {}
83 
84 
85 template<class Type>
87 (
88  const oversetFvPatchField<Type>& ptf,
90 )
91 :
93  oversetPatch_(ptf.oversetPatch_)
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class Type>
101 (
102  const Pstream::commsTypes commsType
103 )
104 {
105  if (this->oversetPatch_.master())
106  {
107  // Trigger interpolation
108  const fvMesh& mesh = this->internalField().mesh();
109  const dictionary& fvSchemes = mesh.schemesDict();
110  const word& fldName = this->internalField().name();
111 
112  if (&mesh.lduAddr() != &mesh.fvMesh::lduAddr())
113  {
114  // Running extended addressing. Flux correction already done
115  // in the linear solver (through the initUpdateInterfaceMatrix
116  // method below)
117  if (debug)
118  {
119  Info<< "Skipping overset interpolation for solved-for field "
120  << fldName << endl;
121  }
122  }
123  else if (!fvSchemes.found("oversetInterpolation"))
124  {
126  << "Missing required dictionary entry"
127  << " 'oversetInterpolation'"
128  << ". Skipping overset interpolation for field "
129  << fldName << endl;
130  }
131  else if (fvSchemes.found("oversetInterpolationRequired"))
132  {
133  // Backwards compatibility mode: only interpolate what is
134  // explicitly mentioned
135 
136  if (fvSchemes.found("oversetInterpolationSuppressed"))
137  {
139  << "Cannot have both dictionary entry"
140  << " 'oversetInterpolationSuppresed' and "
141  << " 'oversetInterpolationRequired' for field "
142  << fldName << exit(FatalIOError);
143  }
144  const dictionary& intDict = fvSchemes.subDict
145  (
146  "oversetInterpolationRequired"
147  );
148  if (intDict.found(fldName))
149  {
150  if (debug)
151  {
152  Info<< "Interpolating field " << fldName << endl;
153  }
154 
155  // Interpolate without bc update (since would trigger infinite
156  // recursion back to oversetFvPatchField<Type>::evaluate)
157  // The only problem is bcs that use the new cell values
158  // (e.g. zeroGradient, processor). These need to appear -after-
159  // the 'overset' bc.
160  mesh.interpolate
161  (
162  const_cast<Field<Type>&>
163  (
164  this->primitiveField()
165  )
166  );
167  }
168  else if (debug)
169  {
170  Info<< "Skipping overset interpolation for field "
171  << fldName << endl;
172  }
173  }
174  else
175  {
176  const dictionary* dictPtr
177  (
178  fvSchemes.findDict("oversetInterpolationSuppressed")
179  );
180 
181  const wordHashSet& suppress =
182  Stencil::New(mesh).nonInterpolatedFields();
183 
184  bool skipInterpolate = suppress.found(fldName);
185 
186  if (dictPtr)
187  {
188  skipInterpolate =
189  skipInterpolate
190  || dictPtr->found(fldName);
191  }
192 
193  if (skipInterpolate)
194  {
195  if (debug)
196  {
197  Info<< "Skipping suppressed overset interpolation"
198  << " for field " << fldName << endl;
199  }
200  }
201  else
202  {
203  if (debug)
204  {
205  Info<< "Interpolating non-suppressed field " << fldName
206  << endl;
207  }
208  mesh.interpolate
209  (
210  const_cast<Field<Type>&>
211  (
212  this->primitiveField()
213  )
214  );
215  }
216  }
217  }
218 }
219 
220 
221 template<class Type>
223 {
225  // Make sure to write the value for ease of postprocessing.
226  this->writeEntry("value", os);
227 }
228 
229 
230 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
volFields.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:52
dynamicOversetFvMesh.H
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::HashSet< word, Hash< word > >
Foam::oversetFvPatchField::oversetFvPatchField
oversetFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: oversetFvPatchField.C:37
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:64
Foam::oversetFvPatchField::write
virtual void write(Ostream &os) const
Write.
Definition: oversetFvPatchField.C:222
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::oversetFvPatchField::oversetPatch_
const oversetFvPatch & oversetPatch_
Local reference cast into the overset patch.
Definition: oversetFvPatchField.H:65
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::oversetFvPatchField
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
Definition: oversetFvPatchField.H:56
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
IOWarningInFunction
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Definition: messageStream.H:340
cellCellStencil.H
cellCellStencilObject.H
Foam::oversetFvPatchField::initEvaluate
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
Definition: oversetFvPatchField.C:101
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54