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