nonConformalCyclicLagrangianPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 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 
27 #include "LagrangianFields.H"
28 #include "RemoteData.H"
29 #include "tracking.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  template<class Type>
37  struct maxFirstOp
38  {
39  const Type& operator()(const Type& a, const Type& b) const
40  {
41  return a.first() > b.first() ? a : b;
42  }
43  };
44 }
45 
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
52 
54  (
57  polyPatch
58  );
59 }
60 
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 (
67  const polyPatch& patch,
68  const LagrangianBoundaryMesh& boundaryMesh
69 )
70 :
71  LagrangianPatch(patch, boundaryMesh),
72  nonConformalCyclicPatch_(refCast<const nonConformalCyclicPolyPatch>(patch)),
73  isNbrPatchMesh_(false),
74  nPositionalErrors_(0),
75  maxPositionalErrorSqr_(-vGreat),
76  maxPositionalErrorReceivePosition_(point::uniform(NaN))
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
87 
90 {
91  return
92  boundaryMesh()
93  [
94  isNbrPatchMesh_
95  ? nonConformalCyclicPatch_.nbrPatchIndex()
96  : patch().index()
97  ].LagrangianPatch::mesh();
98 }
99 
100 
102 (
106 ) const
107 {
108  const LagrangianSubMesh& patchMesh = this->mesh();
109 
110  // Sub-set the geometry and topology of the elements
111  SubField<barycentric> patchCoordinates = patchMesh.sub(mesh.coordinates());
112  SubField<label> patchCelli = patchMesh.sub(mesh.celli());
113  SubField<label> patchFacei = patchMesh.sub(mesh.facei());
114  SubField<label> patchFaceTrii = patchMesh.sub(mesh.faceTrii());
115 
116  // Sub-set the receiving information
117  SubField<label> receivePatchFace =
118  patchMesh.sub(mesh.receivePatchFacePtr_());
119  SubField<vector> receivePosition =
120  patchMesh.sub(mesh.receivePositionPtr_());
121 
122  // Get a reference to the receiving original patch
123  const polyPatch& receivePp =
124  nonConformalCyclicPatch_.nbrPatch().origPatch();
125 
126  // Search for the elements on the receiving side
127  forAll(patchMesh, i)
128  {
129  patchCelli[i] =
130  mesh.mesh().faceOwner()[receivePatchFace[i] + receivePp.start()];
131 
132  if
133  (
135  (
136  mesh.mesh(),
137  receivePosition[i],
138  patchCoordinates[i],
139  patchCelli[i],
140  patchFacei[i],
141  patchFaceTrii[i],
142  fraction[i + patchMesh.start()]
143  )
144  )
145  {
146  // The search hit a boundary. Compute the positional error.
147  const scalar positionalErrorSqr =
148  magSqr
149  (
150  receivePosition[i]
152  (
153  mesh.mesh(),
154  patchCoordinates[i],
155  patchCelli[i],
156  patchFacei[i],
157  patchFaceTrii[i],
158  fraction[i + patchMesh.start()]
159  )
160  );
161 
162  // If the positional error is significant, relative to the size of
163  // the receiving face, then register as a failure
164  if
165  (
166  sqr(positionalErrorSqr)
167  > sqr(small)*magSqr(receivePp.faceAreas()[receivePatchFace[i]])
168  )
169  {
170  nPositionalErrors_ ++;
171 
172  if (positionalErrorSqr > maxPositionalErrorSqr_)
173  {
174  maxPositionalErrorSqr_ = positionalErrorSqr;
175  maxPositionalErrorReceivePosition_ = receivePosition[i];
176  }
177  }
178  }
179  }
180 
181  // Set the elements to be in the adjacent cell
182  patchMesh.sub(mesh.states()) = LagrangianState::inCell;
183 
184  // The elements are now on the neighbour patch
185  isNbrPatchMesh_ = true;
186 
187  // Invalidate the now used receiving information
188  receivePatchFace = -1;
189  receivePosition = point::nan;
190 }
191 
192 
194 {
196 
197  // The elements are back on this patch
198  isNbrPatchMesh_ = false;
199 
200  // Check for positional errors and report if any are found
201  const label nPositionalErrors =
202  returnReduce(nPositionalErrors_, sumOp<label>());
203 
204  if (nPositionalErrors != 0)
205  {
206  const label nTransfers =
207  returnReduce(this->mesh().size(), sumOp<label>());
208 
209  const Tuple2<scalar, vector> maxPositionalErrorInfo =
211  (
213  (
214  maxPositionalErrorSqr_,
215  maxPositionalErrorReceivePosition_
216  ),
218  );
219 
220  maxPositionalErrorSqr_ = maxPositionalErrorInfo.first();
221  maxPositionalErrorReceivePosition_ = maxPositionalErrorInfo.second();
222 
224  << nPositionalErrors << "/" << nTransfers << " elements "
225  << "transferring to patch "
226  << nonConformalCyclicPatch_.nbrPatch().name()
227  << " were not accurately located. The largest positional error "
228  << "was " << sqrt(maxPositionalErrorSqr_)
229  << " from " << maxPositionalErrorReceivePosition_ << "."
230  << nl;
231  }
232 
233  nPositionalErrors_ = 0;
234  maxPositionalErrorSqr_ = -vGreat;
235  maxPositionalErrorReceivePosition_ = point::uniform(NaN);
236 }
237 
238 
239 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Boundary part of a Lagrangian mesh. Just a list of Lagrangian patches with some added convenience fun...
Class containing Lagrangian geometry and topology.
Base class for Lagrangian patches.
virtual void partition() const
Update for mesh changes.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
SubList< Type > sub(const List< Type > &list) const
Return a sub-list corresponding to this sub-mesh.
label start() const
Return start.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Pre-declare related SubField type.
Definition: SubField.H:63
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
static const Form nan
Definition: VectorSpace.H:124
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
nonConformalCyclicLagrangianPatch(const polyPatch &, const LagrangianBoundaryMesh &)
Construct from a patch and a boundary mesh.
virtual void partition() const
Update following partitioning of the mesh.
virtual void evaluate(PstreamBuffers &, LagrangianMesh &, const LagrangianScalarInternalDynamicField &fraction) const
Evaluate changes in elements that have tracked to this patch.
virtual const LagrangianSubMesh & mesh() const
Return the sub-mesh associated with this patch.
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const vectorField::subField faceAreas() const
Return face areas.
Definition: polyPatch.C:282
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
volScalarField & b
Definition: createFields.H:27
#define WarningInFunction
Report a warning using Foam::Warning.
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
bool locate(const polyMesh &mesh, const point &position, barycentric &coordinates, label &celli, label &facei, label &faceTrii, const scalar stepFraction, const string &debugPrefix=NullObjectRef< string >())
Initialise the location at the given position. Returns whether or not a.
Definition: tracking.C:1593
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
defineTypeNameAndDebug(combustionModel, 0)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
static const char nl
Definition: Ostream.H:267
const Type & operator()(const Type &a, const Type &b) const