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-2017 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 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
53  coupledFvPatchField<Type>(p, iF, dict, dict.found("value")),
54  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p))
55 {
56  if (!isA<cyclicACMIFvPatch>(p))
57  {
59  (
60  dict
61  ) << " patch type '" << p.type()
62  << "' not constraint type '" << typeName << "'"
63  << "\n for patch " << p.name()
64  << " of field " << this->internalField().name()
65  << " in file " << this->internalField().objectPath()
66  << exit(FatalIOError);
67  }
68 
69  if (!dict.found("value") && this->coupled())
70  {
71  this->evaluate(Pstream::commsTypes::blocking);
72  }
73 }
74 
75 
76 template<class Type>
78 (
80  const fvPatch& p,
82  const fvPatchFieldMapper& mapper
83 )
84 :
86  coupledFvPatchField<Type>(ptf, p, iF, mapper),
87  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p))
88 {
89  if (!isA<cyclicACMIFvPatch>(this->patch()))
90  {
92  << "' not constraint type '" << typeName << "'"
93  << "\n for patch " << p.name()
94  << " of field " << this->internalField().name()
95  << " in file " << this->internalField().objectPath()
96  << exit(FatalIOError);
97  }
98 }
99 
100 
101 
102 template<class Type>
104 (
106 )
107 :
110  cyclicACMIPatch_(ptf.cyclicACMIPatch_)
111 {}
112 
113 
114 template<class Type>
116 (
117  const cyclicACMIFvPatchField<Type>& ptf,
119 )
120 :
122  coupledFvPatchField<Type>(ptf, iF),
123  cyclicACMIPatch_(ptf.cyclicACMIPatch_)
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
129 template<class Type>
131 {
132  return cyclicACMIPatch_.coupled();
133 }
134 
135 
136 template<class Type>
139 {
140  const Field<Type>& iField = this->primitiveField();
141  const cyclicACMIPolyPatch& cpp = cyclicACMIPatch_.cyclicACMIPatch();
142  tmp<Field<Type>> tpnf
143  (
144  cyclicACMIPatch_.interpolate
145  (
147  (
148  iField,
149  cpp.neighbPatch().faceCells()
150  )
151  )
152  );
153 
154  if (doTransform())
155  {
156  tpnf.ref() = transform(forwardT(), tpnf());
157  }
158 
159  return tpnf;
160 }
161 
162 
163 template<class Type>
166 {
169  (
170  this->primitiveField()
171  );
172 
173  return refCast<const cyclicACMIFvPatchField<Type>>
174  (
175  fld.boundaryField()[cyclicACMIPatch_.neighbPatchID()]
176  );
177 }
178 
179 
180 template<class Type>
183 {
186  (
187  this->primitiveField()
188  );
189 
190  return fld.boundaryField()[cyclicACMIPatch_.nonOverlapPatchID()];
191 }
192 
193 
194 template<class Type>
196 (
197  scalarField& result,
198  const scalarField& psiInternal,
199  const scalarField& coeffs,
200  const direction cmpt,
201  const Pstream::commsTypes
202 ) const
203 {
204  const cyclicACMIPolyPatch& cpp = cyclicACMIPatch_.cyclicACMIPatch();
205 
206  // note: only applying coupled contribution
207 
208  const labelUList& nbrFaceCellsCoupled =
209  cpp.neighbPatch().faceCells();
210 
211  scalarField pnf(psiInternal, nbrFaceCellsCoupled);
212 
213  // Transform according to the transformation tensors
214  transformCoupleField(pnf, cmpt);
215 
216  const labelUList& faceCells = cyclicACMIPatch_.faceCells();
217 
218  pnf = cyclicACMIPatch_.interpolate(pnf);
219 
220  forAll(faceCells, elemI)
221  {
222  result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
223  }
224 }
225 
226 
227 template<class Type>
229 (
230  Field<Type>& result,
231  const Field<Type>& psiInternal,
232  const scalarField& coeffs,
233  const Pstream::commsTypes
234 ) const
235 {
236  const cyclicACMIPolyPatch& cpp = cyclicACMIPatch_.cyclicACMIPatch();
237 
238  // note: only applying coupled contribution
239 
240  const labelUList& nbrFaceCellsCoupled = cpp.neighbPatch().faceCells();
241 
242  Field<Type> pnf(psiInternal, nbrFaceCellsCoupled);
243 
244  // Transform according to the transformation tensors
245  transformCoupleField(pnf);
246 
247  const labelUList& faceCells = cyclicACMIPatch_.faceCells();
248 
249  pnf = cyclicACMIPatch_.interpolate(pnf);
250 
251  forAll(faceCells, elemI)
252  {
253  result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
254  }
255 }
256 
257 
258 template<class Type>
260 (
261  fvMatrix<Type>& matrix
262 )
263 {
264  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
265 
266  // nothing to be done by the AMI, but re-direct to non-overlap patch
267  // with non-overlap patch weights
268  const fvPatchField<Type>& npf = nonOverlapPatchField();
269 
270  const_cast<fvPatchField<Type>&>(npf).manipulateMatrix(matrix, 1.0 - mask);
271 }
272 
273 
274 template<class Type>
276 {
277  // Update non-overlap patch - some will implement updateCoeffs, and
278  // others will implement evaluate
279 
280  // Pass in (1 - mask) to give non-overlap patch the chance to do
281  // manipulation of non-face based data
282 
283  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
284  const fvPatchField<Type>& npf = nonOverlapPatchField();
285  const_cast<fvPatchField<Type>&>(npf).updateWeightedCoeffs(1.0 - mask);
286 }
287 
288 
289 template<class Type>
291 {
293  this->writeEntry("value", os);
294 }
295 
296 
297 // ************************************************************************* //
dictionary dict
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual const cyclicACMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
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
const Boundary & boundaryField() const
Return const-reference to the boundary field.
commsTypes
Types of communications.
Definition: UPstream.H:64
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
uint8_t direction
Definition: direction.H:45
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:66
Generic GeometricField class.
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
Spatial transformation functions for primitive fields.
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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 labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
virtual void write(Ostream &os) const
Write.
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:72
const cyclicACMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
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:61
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.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
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.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
IOerror FatalIOError