fvMeshStitcherTemplates.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) 2022-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 Description
25  Perform mapping of finite volume fields required by stitching.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvMeshStitcher.H"
30 #include "conformedFvPatchField.H"
31 #include "conformedFvsPatchField.H"
33 #include "setSizeFieldMapper.H"
34 #include "surfaceFields.H"
35 #include "volFields.H"
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 template<class Type, template<class> class GeoField>
40 void Foam::fvMeshStitcher::resizePatchFields()
41 {
42  UPtrList<GeoField<Type>> fields(mesh_.fields<GeoField<Type>>());
43 
44  forAll(fields, i)
45  {
46  forAll(mesh_.boundary(), patchi)
47  {
48  typename GeoField<Type>::Patch& pf =
49  fields[i].boundaryFieldRefNoStoreOldTimes()[patchi];
50 
51  if (isA<nonConformalFvPatch>(pf.patch()))
52  {
53  pf.map(pf, setSizeFieldMapper(pf.patch().size()));
54  }
55  }
56  }
57 }
58 
59 
60 template<template<class> class GeoField>
61 void Foam::fvMeshStitcher::resizePatchFields()
62 {
63  #define ResizePatchFields(Type, nullArg) \
64  resizePatchFields<Type, GeoField>();
66  #undef ResizePatchFields
67 }
68 
69 
70 template<class Type>
71 void Foam::fvMeshStitcher::preConformSurfaceFields()
72 {
73  UPtrList<SurfaceField<Type>> fields(mesh_.fields<SurfaceField<Type>>());
74 
75  forAll(fields, i)
76  {
78  (
79  fields[i].boundaryFieldRefNoStoreOldTimes()
80  );
81  }
82 }
83 
84 
85 template<class Type>
86 void Foam::fvMeshStitcher::preConformVolFields()
87 {
88  UPtrList<VolField<Type>> fields(mesh_.fields<VolField<Type>>());
89 
90  forAll(fields, i)
91  {
93  (
94  fields[i].boundaryFieldRefNoStoreOldTimes()
95  );
96  }
97 }
98 
99 
100 template<class Type>
101 void Foam::fvMeshStitcher::postUnconformSurfaceFields()
102 {
103  if (mesh_.topoChanged())
104  {
105  UPtrList<SurfaceField<Type>> curFields
106  (
107  mesh_.curFields<SurfaceField<Type>>()
108  );
109 
110  forAll(curFields, i)
111  {
112  curFields[i].clearOldTimes();
113  }
114  }
115 
116  UPtrList<SurfaceField<Type>> fields(mesh_.fields<SurfaceField<Type>>());
117 
118  forAll(fields, i)
119  {
121  (
122  fields[i].boundaryFieldRefNoStoreOldTimes()
123  );
124  }
125 }
126 
127 
128 template<class Type>
129 void Foam::fvMeshStitcher::postUnconformVolFields()
130 {
131  UPtrList<VolField<Type>> fields(mesh_.fields<VolField<Type>>());
132 
133  forAll(fields, i)
134  {
136  (
137  fields[i].boundaryFieldRefNoStoreOldTimes()
138  );
139  }
140 }
141 
142 
143 template<class Type>
144 void Foam::fvMeshStitcher::postUnconformEvaluateVolFields()
145 {
146  auto evaluate = [](const typename VolField<Type>::Patch& pf)
147  {
148  return
149  (
150  isA<nonConformalFvPatch>(pf.patch())
151  && pf.type() == pf.patch().patch().type()
152  && polyPatch::constraintType(pf.patch().patch().type())
153  )
154  || isA<nonConformalErrorFvPatch>(pf.patch());
155  };
156 
157  UPtrList<VolField<Type>> fields(mesh_.fields<VolField<Type>>());
158 
159  forAll(fields, i)
160  {
161  const label nReq = Pstream::nRequests();
162 
163  forAll(mesh_.boundary(), patchi)
164  {
165  typename VolField<Type>::Patch& pf =
166  fields[i].boundaryFieldRefNoStoreOldTimes()[patchi];
167 
168  if (evaluate(pf))
169  {
170  pf.initEvaluate(Pstream::defaultCommsType);
171  }
172  }
173 
174  if
175  (
178  )
179  {
180  Pstream::waitRequests(nReq);
181  }
182 
183  forAll(mesh_.boundary(), patchi)
184  {
185  typename VolField<Type>::Patch& pf =
186  fields[i].boundaryFieldRefNoStoreOldTimes()[patchi];
187 
188  if (evaluate(pf))
189  {
190  pf.evaluate(Pstream::defaultCommsType);
191  }
192  }
193  }
194 }
195 
196 
197 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
GeoMesh::template PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static void conform(typename VolField< Type >::Boundary &bF)
Conform the given boundary field.
static void unconform(typename VolField< Type >::Boundary &bF)
Un-conform the given boundary field.
static void unconform(typename SurfaceField< Type >::Boundary &bF)
Un-conform the given boundary field.
static void conform(typename SurfaceField< Type >::Boundary &bF)
Conform the given boundary field.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
UPtrList< GeoField > fields(bool strict=false, const HashSet< word > &geometryFields=fvMesh::geometryFields) const
Return the list of fields of type GeoField.
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
#define ResizePatchFields(Type, nullArg)
label patchi
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:234
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
void evaluate(GeometricField< Type, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, GeoMesh > &x)
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
Foam::surfaceFields.