cyclicFaPatch.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-2017 Wikki 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 "cyclicFaPatch.H"
29 #include "coupledPolyPatch.H"
31 #include "transform.H"
32 #include "faMesh.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(cyclicFaPatch, 0);
39  addToRunTimeSelectionTable(faPatch, cyclicFaPatch, dictionary);
40 }
41 
42 const Foam::scalar Foam::cyclicFaPatch::matchTol_ = 1e-3;
43 
44 
45 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46 
47 void Foam::cyclicFaPatch::calcTransforms()
48 {
49  if (size() > 0)
50  {
51  pointField half0Ctrs(size()/2);
52  pointField half1Ctrs(size()/2);
53  for (label i=0; i<size()/2; ++i)
54  {
55  half0Ctrs[i] = this->edgeCentres()[i];
56  half1Ctrs[i] = this->edgeCentres()[i+size()/2];
57  }
58 
59  vectorField half0Normals(size()/2);
60  vectorField half1Normals(size()/2);
61 
63 
64  scalar maxMatchError = 0;
65  label errorEdge = -1;
66 
67  for (label edgei = 0; edgei < size()/2; ++edgei)
68  {
69  half0Normals[edgei] = eN[edgei];
70  label nbrEdgei = edgei + size()/2;
71  half1Normals[edgei] = eN[nbrEdgei];
72 
73  scalar magLe = mag(half0Normals[edgei]);
74  scalar nbrMagLe = mag(half1Normals[edgei]);
75  scalar avLe = (magLe + nbrMagLe)/2.0;
76 
77  if (magLe < ROOTVSMALL && nbrMagLe < ROOTVSMALL)
78  {
79  // Undetermined normal. Use dummy normal to force separation
80  // check. (note use of sqrt(VSMALL) since that is how mag
81  // scales)
82  half0Normals[edgei] = point(1, 0, 0);
83  half1Normals[edgei] = half0Normals[edgei];
84  }
85  else if (mag(magLe - nbrMagLe)/avLe > matchTol_)
86  {
87  // Error in area matching. Find largest error
88  maxMatchError =
89  Foam::max(maxMatchError, mag(magLe - nbrMagLe)/avLe);
90  errorEdge = edgei;
91  }
92  else
93  {
94  half0Normals[edgei] /= magLe;
95  half1Normals[edgei] /= nbrMagLe;
96  }
97  }
98 
99  // Check for error in edge matching
100  if (maxMatchError > matchTol_)
101  {
102  label nbrEdgei = errorEdge + size()/2;
103  scalar magLe = mag(half0Normals[errorEdge]);
104  scalar nbrMagLe = mag(half1Normals[errorEdge]);
105  scalar avLe = (magLe + nbrMagLe)/2.0;
106 
108  << "edge " << errorEdge
109  << " area does not match neighbour "
110  << nbrEdgei << " by "
111  << 100*mag(magLe - nbrMagLe)/avLe
112  << "% -- possible edge ordering problem." << endl
113  << "patch:" << name()
114  << " my area:" << magLe
115  << " neighbour area:" << nbrMagLe
116  << " matching tolerance:" << matchTol_
117  << endl
118  << "Mesh edge:" << start() + errorEdge
119  << endl
120  << "Neighbour edge:" << start() + nbrEdgei
121  << endl
122  << "Other errors also exist, only the largest is reported. "
123  << "Please rerun with cyclic debug flag set"
124  << " for more information." << exit(FatalError);
125  }
126 
127  // Calculate transformation tensors
129  (
130  half0Ctrs,
131  half1Ctrs,
132  half0Normals,
133  half1Normals
134  );
135 
136  // Check transformation tensors
137  if (!parallel())
138  {
139  if (forwardT().size() > 1 || reverseT().size() > 1)
140  {
142  << "Transformation tensor is not constant for the cyclic "
143  << "patch. Please reconsider your setup and definition of "
144  << "cyclic boundaries." << endl;
145  }
146  }
147  }
148 }
149 
150 
152 {
153  const scalarField& magL = magEdgeLengths();
154 
155  const scalarField deltas(edgeNormals() & faPatch::delta());
156  label sizeby2 = deltas.size()/2;
157 
158  scalar maxMatchError = 0;
159  label errorEdge = -1;
160 
161  for (label edgei = 0; edgei < sizeby2; ++edgei)
162  {
163  scalar avL = (magL[edgei] + magL[edgei + sizeby2])/2.0;
164 
165  if
166  (
167  mag(magL[edgei] - magL[edgei + sizeby2])/avL
168  > matchTol_
169  )
170  {
171  // Found error. Look for largest matching error
172  maxMatchError =
173  Foam::max
174  (
175  maxMatchError,
176  mag(magL[edgei] - magL[edgei + sizeby2])/avL
177  );
178 
179  errorEdge = edgei;
180  }
181 
182  scalar di = deltas[edgei];
183  scalar dni = deltas[edgei + sizeby2];
184 
185  w[edgei] = dni/(di + dni);
186  w[edgei + sizeby2] = 1 - w[edgei];
187  }
188 
189  // Check for error in matching
190  if (maxMatchError > matchTol_)
191  {
192  scalar avL = (magL[errorEdge] + magL[errorEdge + sizeby2])/2.0;
193 
195  << "edge " << errorEdge << " and " << errorEdge + sizeby2
196  << " areas do not match by "
197  << 100*mag(magL[errorEdge] - magL[errorEdge + sizeby2])/avL
198  << "% -- possible edge ordering problem." << nl
199  << "Cyclic area match tolerance = "
200  << matchTol_ << " patch: " << name()
201  << abort(FatalError);
202  }
203 }
204 
205 
207 {
208  const scalarField deltas(edgeNormals() & faPatch::delta());
209  label sizeby2 = deltas.size()/2;
210 
211  for (label edgei = 0; edgei < sizeby2; ++edgei)
212  {
213  scalar di = deltas[edgei];
214  scalar dni = deltas[edgei + sizeby2];
215 
216  dc[edgei] = 1.0/(di + dni);
217  dc[edgei + sizeby2] = dc[edgei];
218  }
219 }
220 
221 
223 {
225 }
226 
227 
229 {
231  calcTransforms();
232 }
233 
234 
236 {
238 }
239 
240 
242 {
244  calcTransforms();
245 }
246 
247 
249 {
250  const vectorField patchD(faPatch::delta());
251  label sizeby2 = patchD.size()/2;
252 
253  tmp<vectorField> tpdv(new vectorField(patchD.size()));
254  vectorField& pdv = tpdv.ref();
255 
256  // Do the transformation if necessary
257  if (parallel())
258  {
259  for (label edgei = 0; edgei < sizeby2; ++edgei)
260  {
261  const vector& ddi = patchD[edgei];
262  vector dni = patchD[edgei + sizeby2];
263 
264  pdv[edgei] = ddi - dni;
265  pdv[edgei + sizeby2] = -pdv[edgei];
266  }
267  }
268  else
269  {
270  for (label edgei = 0; edgei < sizeby2; ++edgei)
271  {
272  const vector& ddi = patchD[edgei];
273  vector dni = patchD[edgei + sizeby2];
274 
275  pdv[edgei] = ddi - transform(forwardT()[0], dni);
276  pdv[edgei + sizeby2] = -transform(reverseT()[0], pdv[edgei]);
277  }
278  }
279 
280  return tpdv;
281 }
282 
283 
285 (
286  const labelUList& internalData
287 ) const
288 {
289  return patchInternalField(internalData);
290 }
291 
292 
294 (
295  const Pstream::commsTypes,
296  const labelUList& interfaceData
297 ) const
298 {
299  tmp<labelField> tpnf(new labelField(this->size()));
300  labelField& pnf = tpnf.ref();
301 
302  label sizeby2 = this->size()/2;
303 
304  for (label edgei=0; edgei<sizeby2; ++edgei)
305  {
306  pnf[edgei] = interfaceData[edgei + sizeby2];
307  pnf[edgei + sizeby2] = interfaceData[edgei];
308  }
309 
310  return tpnf;
311 }
312 
313 
315 (
316  const Pstream::commsTypes commsType,
317  const labelUList& iF
318 ) const
319 {
320  const labelUList& edgeCells = this->faceCells();
321 
322  tmp<labelField> tpnf(new labelField(this->size()));
323  labelField& pnf = tpnf.ref();
324 
325  label sizeby2 = this->size()/2;
326 
327  for (label edgei=0; edgei<sizeby2; ++edgei)
328  {
329  pnf[edgei] = iF[edgeCells[edgei + sizeby2]];
330  pnf[edgei + sizeby2] = iF[edgeCells[edgei]];
331  }
332 
333  return tpnf;
334 }
335 
336 
337 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::cyclicFaPatch::initMovePoints
virtual void initMovePoints(const pointField &)
Initialise the patches for moving points.
Definition: cyclicFaPatch.C:235
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::cyclicFaPatch::calcGeometry
virtual void calcGeometry()
Calculate the patch geometry.
Definition: cyclicFaPatch.C:228
Foam::cyclicFaPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: cyclicFaPatch.C:151
Foam::coupledFaPatch::parallel
bool parallel() const
Are the cyclic planes parallel.
Definition: coupledFaPatch.H:198
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::faPatch::movePoints
virtual void movePoints(const pointField &)
Correct patch after moving points.
Definition: faPatch.C:468
Foam::cyclicFaPatch::makeDeltaCoeffs
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
Definition: cyclicFaPatch.C:206
cyclicFaPatch.H
Foam::faPatch::calcGeometry
virtual void calcGeometry()
Calculate the patch geometry.
Definition: faPatch.H:125
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
faMesh.H
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::faPatch::delta
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:428
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::faPatch::initMovePoints
virtual void initMovePoints(const pointField &)
Initialise the patches for moving points.
Definition: faPatch.H:129
coupledPolyPatch.H
Foam::Field< scalar >
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::faPatch::edgeCentres
const vectorField & edgeCentres() const
Return edge centres.
Definition: faPatch.C:380
SeriousErrorInFunction
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
Definition: messageStream.H:281
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::cyclicFaPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to.
Definition: cyclicFaPatch.C:285
Foam::cyclicFaPatch::matchTol_
static const scalar matchTol_
Relative tolerance (for geometric matching). Is factor of.
Definition: cyclicFaPatch.H:72
Foam::cyclicFaPatch::initGeometry
virtual void initGeometry()
Initialise the calculation of the patch geometry.
Definition: cyclicFaPatch.C:222
Foam::FatalError
error FatalError
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::cyclicFaPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: cyclicFaPatch.H:139
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::cyclicFaPatch::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: cyclicFaPatch.H:145
Foam::faPatch::start
label start() const
Patch start in edge list.
Definition: faPatch.C:132
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::coupledFaPatch::calcTransformTensors
void calcTransformTensors(const vector &Cf, const vector &Cr, const vector &nf, const vector &nr) const
Calculate the uniform transformation tensors.
Definition: coupledFaPatch.C:41
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::faPatch::magEdgeLengths
const scalarField & magEdgeLengths() const
Return edge length magnitudes.
Definition: faPatch.C:392
Foam::Vector< scalar >
Foam::cyclicFaPatch::transfer
virtual tmp< labelField > transfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Transfer and return neighbour field.
Definition: cyclicFaPatch.C:294
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::cyclicFaPatch::movePoints
virtual void movePoints(const pointField &)
Correct patches after moving points.
Definition: cyclicFaPatch.C:241
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::faPatch::edgeNormals
tmp< vectorField > edgeNormals() const
Return edge normals.
Definition: faPatch.C:398
Foam::cyclicFaPatch::delta
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
Definition: cyclicFaPatch.C:248
Foam::faPatch::initGeometry
virtual void initGeometry()
Initialise the calculation of the patch geometry.
Definition: faPatch.H:121
transform.H
3D tensor transformation operations.
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::cyclicFaPatch::internalFieldTransfer
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
Definition: cyclicFaPatch.C:315
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
The patch name.
Definition: patchIdentifier.H:134
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::faPatch::size
virtual label size() const
Patch size.
Definition: faPatch.H:249