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-2016 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->internalField().name()
63  << " in file " << this->internalField().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->internalField().name()
90  << " in file " << this->internalField().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->primitiveField();
140  const cyclicACMIPolyPatch& cpp = cyclicACMIPatch_.cyclicACMIPatch();
141  tmp<Field<Type>> tpnf
142  (
143  cyclicACMIPatch_.interpolate
144  (
146  (
147  iField,
148  cpp.neighbPatch().faceCells()
149  )
150  )
151  );
152 
153  if (doTransform())
154  {
155  tpnf.ref() = transform(forwardT(), tpnf());
156  }
157 
158  return tpnf;
159 }
160 
161 
162 template<class Type>
165 {
168  (
169  this->primitiveField()
170  );
171 
172  return refCast<const cyclicACMIFvPatchField<Type>>
173  (
174  fld.boundaryField()[cyclicACMIPatch_.neighbPatchID()]
175  );
176 }
177 
178 
179 template<class Type>
182 {
185  (
186  this->primitiveField()
187  );
188 
189  return fld.boundaryField()[cyclicACMIPatch_.nonOverlapPatchID()];
190 }
191 
192 
193 template<class Type>
195 (
196  scalarField& result,
197  const scalarField& psiInternal,
198  const scalarField& coeffs,
199  const direction cmpt,
200  const Pstream::commsTypes
201 ) const
202 {
203  const cyclicACMIPolyPatch& cpp = cyclicACMIPatch_.cyclicACMIPatch();
204 
205  // note: only applying coupled contribution
206 
207  const labelUList& nbrFaceCellsCoupled =
208  cpp.neighbPatch().faceCells();
209 
210  scalarField pnf(psiInternal, nbrFaceCellsCoupled);
211 
212  // Transform according to the transformation tensors
213  transformCoupleField(pnf, cmpt);
214 
215  const labelUList& faceCells = cyclicACMIPatch_.faceCells();
216 
217  pnf = cyclicACMIPatch_.interpolate(pnf);
218 
219  forAll(faceCells, elemI)
220  {
221  result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
222  }
223 }
224 
225 
226 template<class Type>
228 (
229  Field<Type>& result,
230  const Field<Type>& psiInternal,
231  const scalarField& coeffs,
232  const Pstream::commsTypes
233 ) const
234 {
235  const cyclicACMIPolyPatch& cpp = cyclicACMIPatch_.cyclicACMIPatch();
236 
237  // note: only applying coupled contribution
238 
239  const labelUList& nbrFaceCellsCoupled = cpp.neighbPatch().faceCells();
240 
241  Field<Type> pnf(psiInternal, nbrFaceCellsCoupled);
242 
243  // Transform according to the transformation tensors
244  transformCoupleField(pnf);
245 
246  const labelUList& faceCells = cyclicACMIPatch_.faceCells();
247 
248  pnf = cyclicACMIPatch_.interpolate(pnf);
249 
250  forAll(faceCells, elemI)
251  {
252  result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
253  }
254 }
255 
256 
257 template<class Type>
259 (
260  fvMatrix<Type>& matrix
261 )
262 {
263  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
264 
265  // nothing to be done by the AMI, but re-direct to non-overlap patch
266  // with non-overlap patch weights
267  const fvPatchField<Type>& npf = nonOverlapPatchField();
268 
269  const_cast<fvPatchField<Type>&>(npf).manipulateMatrix(matrix, 1.0 - mask);
270 }
271 
272 
273 template<class Type>
275 {
276  // Update non-overlap patch - some will implement updateCoeffs, and
277  // others will implement evaluate
278 
279  // Pass in (1 - mask) to give non-overlap patch the chance to do
280  // manipulation of non-face based data
281 
282  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
283  const fvPatchField<Type>& npf = nonOverlapPatchField();
284  const_cast<fvPatchField<Type>&>(npf).updateWeightedCoeffs(1.0 - mask);
285 }
286 
287 
288 template<class Type>
290 {
292  this->writeEntry("value", os);
293 }
294 
295 
296 // ************************************************************************* //
dictionary dict
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
commsTypes
Types of communications.
Definition: UPstream.H:64
virtual const cyclicACMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Generic GeometricField class.
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
Spatial transformation functions for primitive fields.
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
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))
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvPatchFieldMapper.
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:71
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
cyclicACMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Abstract base class for coupled patches.
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
Abstract base class for cyclic ACMI coupled interfaces.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &os) const
Write.
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
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.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
runTime write()
const cyclicACMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
IOerror FatalIOError