nonConformalProcessorCyclicLagrangianPatch.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 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 
39  (
42  polyPatch
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
51 (
52  const polyPatch& patch,
53  const LagrangianBoundaryMesh& boundaryMesh
54 )
55 :
56  processorCyclicLagrangianPatch(patch, boundaryMesh),
57  nonConformalProcessorCyclicPatch_
58  (
60  )
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
74 (
75  PstreamBuffers& pBufs,
78 ) const
79 {
80  const LagrangianSubMesh& patchMesh = this->mesh();
81 
82  // Sub-set the step fractions
83  SubField<scalar> sendFraction = patchMesh.sub(fraction.primitiveField());
84 
85  // Sub-set the receiving information
86  SubField<label> receivePatchFace =
87  patchMesh.sub(mesh.receivePatchFacePtr_());
88  SubField<point> receivePosition =
89  patchMesh.sub(mesh.receivePositionPtr_());
90 
91  // Send
92  UOPstream(nonConformalProcessorCyclicPatch_.neighbProcNo(), pBufs)()
93  << receivePatchFace
94  << receivePosition
95  << sendFraction;
96 
97  // Remove the sent elements
98  patchMesh.sub(mesh.states()) = LagrangianState::toBeRemoved;
99 
100  // Invalidate the now used receiving information
101  receivePatchFace = -1;
102  receivePosition = point::nan;
103 }
104 
105 
107 (
108  PstreamBuffers& pBufs,
111 ) const
112 {
113  // Receive
114  UIPstream uips(nonConformalProcessorCyclicPatch_.neighbProcNo(), pBufs);
115  labelField receivePatchFace(uips);
116  pointField receivePosition(uips);
117  scalarField receiveFraction(uips);
118 
119  // Get a reference to the receiving original patch
120  const polyPatch& receivePp =
121  nonConformalProcessorCyclicPatch_.referPatch().origPatch();
122 
123  // Search for the elements on the receiving side
124  barycentricField receiveCoordinates(receivePatchFace.size());
125  labelField receiveCelli(receivePatchFace.size());
126  labelField receiveFacei(receivePatchFace.size());
127  labelField receiveFaceTrii(receivePatchFace.size());
128  forAll(receivePatchFace, i)
129  {
130  receiveCelli[i] =
131  mesh.mesh().faceOwner()[receivePatchFace[i] + receivePp.start()];
132 
133  if
134  (
136  (
137  mesh.mesh(),
138  receivePosition[i],
139  receiveCoordinates[i],
140  receiveCelli[i],
141  receiveFacei[i],
142  receiveFaceTrii[i],
143  receiveFraction[i]
144  )
145  )
146  {
147  // The search hit a boundary. Compute the positional error.
148  const scalar positionalErrorSqr =
149  magSqr
150  (
151  receivePosition[i]
153  (
154  mesh.mesh(),
155  receiveCoordinates[i],
156  receiveCelli[i],
157  receiveFacei[i],
158  receiveFaceTrii[i],
159  receiveFraction[i]
160  )
161  );
162 
163  // If the positional error is significant, relative to the size of
164  // the receiving face, then register as a failure
165  if
166  (
167  sqr(positionalErrorSqr)
168  > sqr(small)*magSqr(receivePp.faceAreas()[receivePatchFace[i]])
169  )
170  {
171  referPatch().nPositionalErrors_ ++;
172 
173  if (positionalErrorSqr > referPatch().maxPositionalErrorSqr_)
174  {
175  referPatch().maxPositionalErrorSqr_ = positionalErrorSqr;
176  referPatch().maxPositionalErrorReceivePosition_ =
177  receivePosition[i];
178  }
179  }
180  }
181  }
182 
183  // Insert the received elements
184  receiveMeshPtr_.set
185  (
187  (
188  mesh.append
189  (
190  receiveCoordinates,
191  receiveCelli,
192  receiveFacei,
193  receiveFaceTrii
194  )
195  )
196  );
197 
198  // Set the step fractions of the elements
199  mesh.appendSpecifiedField<scalar, LagrangianInternalDynamicField>
200  (
201  receiveMeshPtr_(),
202  const_cast<LagrangianScalarInternalDynamicField&>(fraction),
203  receiveFraction
204  );
205 
206  // Set the elements to be in the adjacent cell
207  receiveMeshPtr_().sub(mesh.states()) = LagrangianState::inCell;
208 }
209 
210 
211 // ************************************************************************* //
#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...
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
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.
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.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Pre-declare related SubField type.
Definition: SubField.H:63
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:57
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:58
static const Form nan
Definition: VectorSpace.H:124
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
virtual void initEvaluate(PstreamBuffers &, LagrangianMesh &, const LagrangianScalarInternalDynamicField &fraction) const
Initialise evaluation of changes in elements that have tracked to.
nonConformalProcessorCyclicLagrangianPatch(const polyPatch &, const LagrangianBoundaryMesh &)
Construct from a patch and a boundary mesh.
virtual void evaluate(PstreamBuffers &, LagrangianMesh &, const LagrangianScalarInternalDynamicField &fraction) const
Evaluate changes in elements that have tracked to this patch.
Non-conformal processor cyclic poly patch. As nonConformalCyclicPolyPatch, but the neighbouring patch...
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
Processor-cyclic Lagrangian patch. This is used for the patches that interface between processors acr...
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
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)