cyclicACMIFvPatchField.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) 2013-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 "cyclicACMIFvPatchField.H"
27 #include "transformField.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const fvPatch& p,
36 )
37 :
40  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p))
41 {}
42 
43 
44 template<class Type>
46 (
48  const fvPatch& p,
50  const fvPatchFieldMapper& mapper
51 )
52 :
54  coupledFvPatchField<Type>(ptf, p, iF, mapper),
55  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p))
56 {
57  if (!isA<cyclicACMIFvPatch>(this->patch()))
58  {
60  << "' not constraint type '" << typeName << "'"
61  << "\n for patch " << p.name()
62  << " of field " << this->dimensionedInternalField().name()
63  << " in file " << this->dimensionedInternalField().objectPath()
64  << exit(FatalIOError);
65  }
66 }
67 
68 
69 template<class Type>
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
79  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p))
80 {
81  if (!isA<cyclicACMIFvPatch>(p))
82  {
84  (
85  dict
86  ) << " patch type '" << p.type()
87  << "' not constraint type '" << typeName << "'"
88  << "\n for patch " << p.name()
89  << " of field " << this->dimensionedInternalField().name()
90  << " in file " << this->dimensionedInternalField().objectPath()
91  << exit(FatalIOError);
92  }
93 
94  if (!dict.found("value") && this->coupled())
95  {
96  this->evaluate(Pstream::blocking);
97  }
98 }
99 
100 
101 template<class Type>
103 (
105 )
106 :
109  cyclicACMIPatch_(ptf.cyclicACMIPatch_)
110 {}
111 
112 
113 template<class Type>
115 (
116  const cyclicACMIFvPatchField<Type>& ptf,
118 )
119 :
121  coupledFvPatchField<Type>(ptf, iF),
122  cyclicACMIPatch_(ptf.cyclicACMIPatch_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
128 template<class Type>
130 {
131  return cyclicACMIPatch_.coupled();
132 }
133 
134 
135 template<class Type>
138 {
139  const Field<Type>& iField = this->internalField();
140  const labelUList& nbrFaceCellsCoupled =
141  cyclicACMIPatch_.cyclicACMIPatch().neighbPatch().faceCells();
142  const labelUList& faceCellsNonOverlap =
143  cyclicACMIPatch_.cyclicACMIPatch().nonOverlapPatch().faceCells();
144 
145  Field<Type> pnfCoupled(iField, nbrFaceCellsCoupled);
146  Field<Type> pfNonOverlap(iField, faceCellsNonOverlap);
147 
148  tmp<Field<Type> > tpnf
149  (
150  new Field<Type>
151  (
152  cyclicACMIPatch_.interpolate
153  (
154  pnfCoupled,
155  pfNonOverlap
156  )
157  )
158  );
159 
160  if (doTransform())
161  {
162  tpnf() = transform(forwardT(), tpnf());
163  }
164 
165  return tpnf;
166 }
167 
168 
169 template<class Type>
172 {
175  (
176  this->internalField()
177  );
178 
179  return refCast<const cyclicACMIFvPatchField<Type> >
180  (
181  fld.boundaryField()[cyclicACMIPatch_.neighbPatchID()]
182  );
183 }
184 
185 
186 template<class Type>
189 {
192  (
193  this->internalField()
194  );
195 
196  return fld.boundaryField()[cyclicACMIPatch_.nonOverlapPatchID()];
197 }
198 
199 
200 template<class Type>
202 (
203  scalarField& result,
204  const scalarField& psiInternal,
205  const scalarField& coeffs,
206  const direction cmpt,
207  const Pstream::commsTypes
208 ) const
209 {
210  // note: only applying coupled contribution
211 
212  const labelUList& nbrFaceCellsCoupled =
213  cyclicACMIPatch_.cyclicACMIPatch().neighbPatch().faceCells();
214 
215  scalarField pnf(psiInternal, nbrFaceCellsCoupled);
216 
217  // Transform according to the transformation tensors
218  transformCoupleField(pnf, cmpt);
219 
220  const labelUList& faceCells = cyclicACMIPatch_.faceCells();
221 
222  pnf = cyclicACMIPatch_.interpolate(pnf);
223 
224  forAll(faceCells, elemI)
225  {
226  result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
227  }
228 }
229 
230 
231 template<class Type>
233 (
234  Field<Type>& result,
235  const Field<Type>& psiInternal,
236  const scalarField& coeffs,
237  const Pstream::commsTypes
238 ) const
239 {
240  // note: only applying coupled contribution
241 
242  const labelUList& nbrFaceCellsCoupled =
243  cyclicACMIPatch_.cyclicACMIPatch().neighbPatch().faceCells();
244 
245  Field<Type> pnf(psiInternal, nbrFaceCellsCoupled);
246 
247  // Transform according to the transformation tensors
248  transformCoupleField(pnf);
249 
250  const labelUList& faceCells = cyclicACMIPatch_.faceCells();
251 
252  pnf = cyclicACMIPatch_.interpolate(pnf);
253 
254  forAll(faceCells, elemI)
255  {
256  result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
257  }
258 }
259 
260 
261 template<class Type>
263 (
264  const scalarField& deltaCoeffs
265 ) const
266 {
267  // note: only applying coupled contribution
268  return coupledFvPatchField<Type>::snGrad(deltaCoeffs);
269 }
270 
271 
272 template<class Type>
274 {
275  // update non-overlap patch - some will implement updateCoeffs, and
276  // others will implement evaluate
277 
278  // scale neighbour field by (1 - mask)
279 
280  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
281  const fvPatchField<Type>& npf = nonOverlapPatchField();
282  const_cast<fvPatchField<Type>&>(npf).updateCoeffs(1.0 - mask);
283 }
284 
285 
286 template<class Type>
288 (
289  const Pstream::commsTypes comms
290 )
291 {
292  // update non-overlap patch (if not already updated by updateCoeffs)
293 
294  // scale neighbour field by (1 - mask)
295 
296  fvPatchField<Type>& npf =
297  const_cast<fvPatchField<Type>&>(nonOverlapPatchField());
298 
299  if (!npf.updated())
300  {
301  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
302 
303  npf.evaluate(comms);
304 
305  npf *= 1.0 - mask;
306  }
307 }
308 
309 
310 template<class Type>
312 (
313  const Pstream::commsTypes comms
314 )
315 {
316  // blend contributions from the coupled and non-overlap patches
317 
318  // neighbour patch field is updated via updateCoeffs or initEvaluate
319  // and is already scaled by (1 - mask)
320  const fvPatchField<Type>& npf = nonOverlapPatchField();
321 
323  const Field<Type>& cpf = *this;
324 
325  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
326  Field<Type>::operator=(mask*cpf + npf);
327 
329 }
330 
331 
332 template<class Type>
335 (
336  const tmp<scalarField>& w
337 ) const
338 {
339  // note: do not blend based on mask field
340  // - when applied this is scaled by the areas which are already scaled
342 }
343 
344 
345 template<class Type>
348 (
349  const tmp<scalarField>& w
350 ) const
351 {
352  // note: do not blend based on mask field
353  // - when applied this is scaled by the areas which are already scaled
355 }
356 
357 
358 template<class Type>
361 (
362  const scalarField& deltaCoeffs
363 ) const
364 {
365  // note: do not blend based on mask field
366  // - when applied this is scaled by the areas which are already scaled
368 }
369 
370 
371 template<class Type>
374 {
375  // note: do not blend based on mask field
376  // - when applied this is scaled by the areas which are already scaled
378 }
379 
380 
381 template<class Type>
384 (
385  const scalarField& deltaCoeffs
386 ) const
387 {
388  // note: do not blend based on mask field
389  // - when applied this is scaled by the areas which are already scaled
391 }
392 
393 
394 template<class Type>
397 {
398  // note: do not blend based on mask field
399  // - when applied this is scaled by the areas which are already scaled
401 }
402 
403 
404 template<class Type>
406 (
407  fvMatrix<Type>& matrix
408 )
409 {
410  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
411 
412  // nothing to be done by the AMI, but re-direct to non-overlap patch
413  // with non-overlap patch weights
414  const fvPatchField<Type>& npf = nonOverlapPatchField();
415 
416  const_cast<fvPatchField<Type>&>(npf).manipulateMatrix(matrix, 1.0 - mask);
417 }
418 
419 
420 template<class Type>
422 {
424  this->writeEntry("value", os);
425 }
426 
427 
428 // ************************************************************************* //
Foam::fvPatchField< Type >
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::cyclicACMIFvPatchField::valueInternalCoeffs
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:335
Foam::fvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:318
Foam::cyclicACMIFvPatchField::initEvaluate
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
Definition: cyclicACMIFvPatchField.C:288
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
p
p
Definition: pEqn.H:62
Foam::cyclicACMIFvPatchField::nonOverlapPatchField
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
Definition: cyclicACMIFvPatchField.C:188
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::cyclicACMIFvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
Definition: cyclicACMIFvPatchField.C:312
Foam::cyclicACMIFvPatchField
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
Definition: cyclicACMIFvPatchField.H:69
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::FatalIOError
IOerror FatalIOError
Foam::cyclicACMIFvPatchField::gradientInternalCoeffs
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:373
Foam::cyclicACMILduInterfaceField
Abstract base class for cyclic ACMI coupled interfaces.
Definition: cyclicACMILduInterfaceField.H:48
Foam::fvPatchField::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:348
transformField.H
Spatial transformation functions for primitive fields.
Foam::transform
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
Foam::cyclicACMIFvPatchField::updateCoeffs
void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: cyclicACMIFvPatchField.C:273
Foam::cyclicACMIFvPatchField::neighbourPatchField
const cyclicACMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
Definition: cyclicACMIFvPatchField.C:171
Foam::cyclicACMIFvPatchField::valueBoundaryCoeffs
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:348
Foam::cyclicACMIFvPatchField::updateInterfaceMatrix
virtual void updateInterfaceMatrix(scalarField &result, const scalarField &psiInternal, const scalarField &coeffs, const direction cmpt, const Pstream::commsTypes commsType) const
Update result field based on interface functionality.
Definition: cyclicACMIFvPatchField.C:202
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field< Type >
Foam::cyclicACMIFvPatchField::coupled
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
Definition: cyclicACMIFvPatchField.C:129
Foam::cyclicACMIFvPatchField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: cyclicACMIFvPatchField.C:406
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::cyclicACMIFvPatchField::cyclicACMIPatch_
const cyclicACMIFvPatch & cyclicACMIPatch_
Local reference cast into the cyclic patch.
Definition: cyclicACMIFvPatchField.H:77
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
internalField
conserve internalField()+
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::coupledFvPatchField
Abstract base class for coupled patches.
Definition: coupledFvPatchField.H:54
Foam::cyclicACMIFvPatchField::write
virtual void write(Ostream &os) const
Write.
Definition: cyclicACMIFvPatchField.C:421
Foam::cyclicACMIFvPatchField::cyclicACMIFvPatchField
cyclicACMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: cyclicACMIFvPatchField.C:33
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::cyclicACMIFvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
Definition: cyclicACMIFvPatchField.C:137
Foam::fvMatrix< Type >
Foam::direction
unsigned char direction
Definition: direction.H:43
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::coupledFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: coupledFvPatchField.H:147
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::cyclicACMIFvPatchField::gradientBoundaryCoeffs
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:396
cyclicACMIFvPatchField.H
write
Tcoeff write()
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51