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-2024 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 =
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::preConformVolFields()
72 {
73  UPtrList<VolField<Type>> fields(mesh_.curFields<VolField<Type>>());
74 
75  forAll(fields, i)
76  {
77  VolField<Type>& field = fields[i];
78 
79  for (label ti=0; ti<=field.nOldTimes(false); ti++)
80  {
82  (
83  boundaryFieldRefNoUpdate(field.oldTime(ti))
84  );
85  }
86  }
87 }
88 
89 
90 template<class Type>
91 void Foam::fvMeshStitcher::preConformSurfaceFields()
92 {
93  UPtrList<SurfaceField<Type>> fields(mesh_.curFields<SurfaceField<Type>>());
94 
95  forAll(fields, i)
96  {
97  SurfaceField<Type>& field = fields[i];
98 
99  for (label ti=0; ti<=field.nOldTimes(false); ti++)
100  {
102  (
103  boundaryFieldRefNoUpdate(field.oldTime(ti))
104  );
105  }
106  }
107 }
108 
109 
110 template<class Type>
111 void Foam::fvMeshStitcher::postUnconformVolFields()
112 {
113  UPtrList<VolField<Type>> fields(mesh_.curFields<VolField<Type>>());
114 
115  forAll(fields, i)
116  {
117  VolField<Type>& field = fields[i];
118 
119  for (label ti=0; ti<=field.nOldTimes(false); ti++)
120  {
122  (
123  boundaryFieldRefNoUpdate(field.oldTime(ti))
124  );
125  }
126  }
127 }
128 
129 
130 template<class Type>
131 void Foam::fvMeshStitcher::postUnconformEvaluateVolFields()
132 {
133  auto evaluate = [](const typename VolField<Type>::Patch& pf)
134  {
135  return
136  (
137  isA<nonConformalFvPatch>(pf.patch())
138  && pf.type() == pf.patch().patch().type()
139  && polyPatch::constraintType(pf.patch().patch().type())
140  )
141  || isA<nonConformalErrorFvPatch>(pf.patch());
142  };
143 
144  UPtrList<VolField<Type>> fields(mesh_.fields<VolField<Type>>());
145 
146  forAll(fields, i)
147  {
148  const label nReq = Pstream::nRequests();
149 
150  forAll(mesh_.boundary(), patchi)
151  {
152  typename VolField<Type>::Patch& pf =
153  boundaryFieldRefNoUpdate(fields[i])[patchi];
154 
155  if (evaluate(pf))
156  {
157  pf.initEvaluate(Pstream::defaultCommsType);
158  }
159  }
160 
161  if
162  (
165  )
166  {
167  Pstream::waitRequests(nReq);
168  }
169 
170  forAll(mesh_.boundary(), patchi)
171  {
172  typename VolField<Type>::Patch& pf =
173  boundaryFieldRefNoUpdate(fields[i])[patchi];
174 
175  if (evaluate(pf))
176  {
177  pf.evaluate(Pstream::defaultCommsType);
178  }
179  }
180  }
181 }
182 
183 
184 template<class Type>
185 void Foam::fvMeshStitcher::postUnconformSurfaceFields()
186 {
187  UPtrList<SurfaceField<Type>> fields(mesh_.curFields<SurfaceField<Type>>());
188 
189  forAll(fields, i)
190  {
191  SurfaceField<Type>& field = fields[i];
192 
193  for (label ti=0; ti<=field.nOldTimes(false); ti++)
194  {
196  (
197  boundaryFieldRefNoUpdate(field.oldTime(ti))
198  );
199  }
200  }
201 }
202 
203 
204 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
205 
206 template<class GeoField>
208 (
209  GeoField& fld
210 )
211 {
212  return const_cast<typename GeoField::Boundary&>(fld.boundaryField());
213 }
214 
215 
216 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
static GeoField::Boundary & boundaryFieldRefNoUpdate(GeoField &fld)
Access the boundary field reference of a field, without updating.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
UPtrList< GeoField > fields(const bool strict=false) 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
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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:228
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, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
Foam::surfaceFields.