cyclicAMIFvPatch.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) 2011-2016 OpenFOAM Foundation
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 "cyclicAMIFvPatch.H"
30 #include "fvMesh.H"
31 #include "transform.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(cyclicAMIFvPatch, 0);
38  addToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch);
40  (
41  fvPatch,
42  cyclicAMIFvPatch,
43  polyPatch,
44  cyclicPeriodicAMI
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
50 
52 {
53  return Pstream::parRun() || (this->size() && neighbFvPatch().size());
54 }
55 
56 
58 {
59  if (coupled())
60  {
61  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
62 
63  const scalarField deltas(nf() & coupledFvPatch::delta());
64 
65  tmp<scalarField> tnbrDeltas;
66  if (applyLowWeightCorrection())
67  {
68  tnbrDeltas =
70  (
71  nbrPatch.nf() & nbrPatch.coupledFvPatch::delta(),
72  scalarField(this->size(), 1.0)
73  );
74  }
75  else
76  {
77  tnbrDeltas =
78  interpolate(nbrPatch.nf() & nbrPatch.coupledFvPatch::delta());
79  }
80 
81  const scalarField& nbrDeltas = tnbrDeltas();
82 
83  forAll(deltas, facei)
84  {
85  scalar di = deltas[facei];
86  scalar dni = nbrDeltas[facei];
87 
88  w[facei] = dni/(di + dni);
89  }
90  }
91  else
92  {
93  // Behave as uncoupled patch
95  }
96 }
97 
98 
100 {
101  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
102 
103  if (coupled())
104  {
105  const vectorField patchD(coupledFvPatch::delta());
106 
107  tmp<vectorField> tnbrPatchD;
108  if (applyLowWeightCorrection())
109  {
110  tnbrPatchD =
112  (
113  nbrPatch.coupledFvPatch::delta(),
114  vectorField(this->size(), Zero)
115  );
116  }
117  else
118  {
119  tnbrPatchD = interpolate(nbrPatch.coupledFvPatch::delta());
120  }
121 
122  const vectorField& nbrPatchD = tnbrPatchD();
123 
124  tmp<vectorField> tpdv(new vectorField(patchD.size()));
125  vectorField& pdv = tpdv.ref();
126 
127  // do the transformation if necessary
128  if (parallel())
129  {
130  forAll(patchD, facei)
131  {
132  const vector& ddi = patchD[facei];
133  const vector& dni = nbrPatchD[facei];
134 
135  pdv[facei] = ddi - dni;
136  }
137  }
138  else
139  {
140  forAll(patchD, facei)
141  {
142  const vector& ddi = patchD[facei];
143  const vector& dni = nbrPatchD[facei];
144 
145  pdv[facei] = ddi - transform(forwardT()[0], dni);
146  }
147  }
148 
149  return tpdv;
150  }
151  else
152  {
153  return coupledFvPatch::delta();
154  }
155 }
156 
157 
159 (
160  const labelUList& internalData
161 ) const
162 {
163  return patchInternalField(internalData);
164 }
165 
166 
168 (
169  const Pstream::commsTypes commsType,
170  const labelUList& iF
171 ) const
172 {
173  return neighbFvPatch().patchInternalField(iF);
174 }
175 
176 
177 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::cyclicAMIFvPatch::internalFieldTransfer
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
Definition: cyclicAMIFvPatch.C:168
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::cyclicAMIFvPatch::delta
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
Definition: cyclicAMIFvPatch.C:99
Foam::coupledFvPatch::delta
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Definition: coupledFvPatch.C:46
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::cyclicAMIFvPatch::neighbFvPatch
const cyclicAMIFvPatch & neighbFvPatch() const
Definition: cyclicAMIFvPatch.H:149
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::fvPatch::size
virtual label size() const
Return size.
Definition: fvPatch.H:163
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:120
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::Field< scalar >
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::cyclicAMIFvPatch::coupled
virtual bool coupled() const
Return true if this patch is coupled. This is equivalent.
Definition: cyclicAMIFvPatch.C:51
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::cyclicAMIFvPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to.
Definition: cyclicAMIFvPatch.C:159
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:66
Foam::Vector< scalar >
Foam::cyclicAMIFvPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: cyclicAMIFvPatch.C:57
Foam::UList< label >
Foam::cyclicAMIFvPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIFvPatch.H:53
transform.H
3D tensor transformation operations.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::fvPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:146
cyclicAMIFvPatch.H