cyclicLagrangianPatch.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 
26 #include "cyclicLagrangianPatch.H"
27 #include "LagrangianMesh.H"
28 #include "SubField.H"
29 #include "tracking.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 
39  (
42  polyPatch
43  );
44 }
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const polyPatch& patch,
51  const LagrangianBoundaryMesh& boundaryMesh
52 )
53 :
54  LagrangianPatch(patch, boundaryMesh),
55  cyclicPatch_(refCast<const cyclicPolyPatch>(patch)),
56  isNbrPatchMesh_(false)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
67 
69 {
70  return
71  boundaryMesh()
72  [
73  isNbrPatchMesh_
74  ? cyclicPatch_.nbrPatchIndex()
75  : patch().index()
76  ].LagrangianPatch::mesh();
77 }
78 
79 
81 (
85 ) const
86 {
87  const LagrangianSubMesh& patchMesh = this->mesh();
88 
89  // Sub-set the geometry and topology of the elements
90  SubField<barycentric> patchCoordinates = patchMesh.sub(mesh.coordinates());
91  SubField<label> patchCelli = patchMesh.sub(mesh.celli());
92  SubField<label> patchFacei = patchMesh.sub(mesh.facei());
93  SubField<label> patchFaceTrii = patchMesh.sub(mesh.faceTrii());
94 
95  // Cross between the cyclic patches
96  forAll(patchMesh, i)
97  {
99  (
100  cyclicPatch_,
101  patchCoordinates[i],
102  patchCelli[i],
103  patchFacei[i],
104  patchFaceTrii[i]
105  );
106  }
107 
108  // Set the elements to be in the adjacent cell
109  patchMesh.sub(mesh.states()) = LagrangianState::inCell;
110 
111  // The elements are now on the neighbour patch
112  isNbrPatchMesh_ = true;
113 }
114 
115 
117 {
119 
120  // The elements are back on this patch
121  isNbrPatchMesh_ = false;
122 }
123 
124 
125 // ************************************************************************* //
#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.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Pre-declare related SubField type.
Definition: SubField.H:63
Cyclic Lagrangian patch. This is used for the patches that represent translated or rotated repetition...
cyclicLagrangianPatch(const polyPatch &, const LagrangianBoundaryMesh &)
Construct from a patch and a boundary mesh.
virtual ~cyclicLagrangianPatch()
Destructor.
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.
Cyclic plane patch.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
void crossCyclic(const cyclicPolyPatch &inPatch, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Cross a cyclic patch.
Definition: tracking.C:1736
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
defineTypeNameAndDebug(combustionModel, 0)