fvMeshToFvMeshTemplates.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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMeshToFvMesh.H"
27 #include "setSizeFieldMapper.H"
28 #include "identityFieldMapper.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class Type>
35 void Foam::fvMeshToFvMesh::evaluateConstraintTypes(VolField<Type>& fld)
36 {
37  typename VolField<Type>::Boundary& fldBf = fld.boundaryFieldRef();
38 
39  if
40  (
43  )
44  {
45  label nReq = Pstream::nRequests();
46 
47  forAll(fldBf, patchi)
48  {
49  fvPatchField<Type>& tgtField = fldBf[patchi];
50 
51  if
52  (
53  tgtField.type() == tgtField.patch().patch().type()
54  && polyPatch::constraintType(tgtField.patch().patch().type())
55  )
56  {
57  tgtField.initEvaluate(Pstream::defaultCommsType);
58  }
59  }
60 
61  // Block for any outstanding requests
62  if
63  (
66  )
67  {
69  }
70 
71  forAll(fldBf, patchi)
72  {
73  fvPatchField<Type>& tgtField = fldBf[patchi];
74 
75  if
76  (
77  tgtField.type() == tgtField.patch().patch().type()
78  && polyPatch::constraintType(tgtField.patch().patch().type())
79  )
80  {
81  tgtField.evaluate(Pstream::defaultCommsType);
82  }
83  }
84  }
86  {
87  const lduSchedule& patchSchedule =
88  fld.mesh().globalData().patchSchedule();
89 
90  forAll(patchSchedule, patchEvali)
91  {
92  label patchi = patchSchedule[patchEvali].patch;
93  fvPatchField<Type>& tgtField = fldBf[patchi];
94 
95  if
96  (
97  tgtField.type() == tgtField.patch().patch().type()
98  && polyPatch::constraintType(tgtField.patch().patch().type())
99  )
100  {
101  if (patchSchedule[patchEvali].init)
102  {
103  tgtField.initEvaluate(Pstream::commsTypes::scheduled);
104  }
105  else
106  {
107  tgtField.evaluate(Pstream::commsTypes::scheduled);
108  }
109  }
110  }
111  }
112 }
113 
114 
115 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
116 
117 template<class Type>
119 (
120  const VolField<Type>& srcFld
121 ) const
122 {
123  // Construct target patch fields as copies of source patch fields, but do
124  // not map values yet
125  PtrList<fvPatchField<Type>> tgtPatchFields(tgtMesh_.boundary().size());
126  forAll(patchIndices(), i)
127  {
128  const label srcPatchi = patchIndices()[i].first();
129  const label tgtPatchi = patchIndices()[i].second();
130 
131  if (!tgtPatchFields.set(tgtPatchi))
132  {
133  tgtPatchFields.set
134  (
135  tgtPatchi,
137  (
138  srcFld.boundaryField()[srcPatchi],
139  tgtMesh_.boundary()[tgtPatchi],
141  setSizeFieldMapper(tgtMesh_.boundary()[tgtPatchi].size())
142  )
143  );
144  }
145  }
146 
147  // Create calculated patch fields for any unset target patches. Use
148  // fvPatchField<Type>::New to construct these, rather than using the
149  // calculated constructor directly, so that constraints are maintained.
150  labelList tgtPatchFieldIsUnMapped(tgtPatchFields.size(), false);
151  forAll(tgtPatchFields, tgtPatchi)
152  {
153  if (!tgtPatchFields.set(tgtPatchi))
154  {
155  tgtPatchFields.set
156  (
157  tgtPatchi,
159  (
161  tgtMesh_.boundary()[tgtPatchi],
163  )
164  );
165 
166  tgtPatchFieldIsUnMapped[tgtPatchi] =
168  (
169  tgtMesh_.boundary()[tgtPatchi].patch().type()
170  );
171  }
172  }
173 
174  // Construct the result field
175  tmp<VolField<Type>> ttgtFld =
177  (
178  typedName("interpolate(" + srcFld.name() + ")"),
179  srcToTgt<Type>(srcFld.internalField())(),
180  tgtPatchFields
181  );
182  typename VolField<Type>::Boundary& tgtBfld =
183  ttgtFld.ref().boundaryFieldRef();
184 
185  // Mapped patches
186  forAll(patchIndices(), i)
187  {
188  const label srcPatchi = patchIndices()[i].first();
189  const label tgtPatchi = patchIndices()[i].second();
190 
191  tgtBfld[tgtPatchi].map
192  (
193  srcFld.boundaryField()[srcPatchi],
195  (
196  patchInterpolation(i),
197  tgtPatchStabilisation(i)
198  )
199  );
200  }
201 
202  // Un-mapped patches. Set values to that of the internal cell field.
203  forAll(tgtBfld, patchi)
204  {
205  if (tgtPatchFieldIsUnMapped[patchi])
206  {
207  fvPatchField<Type>& tgtPfld = tgtBfld[patchi];
208  tgtPfld == tgtPfld.patchInternalField();
209  }
210  }
211 
212  // Evaluate constraints
213  evaluateConstraintTypes(ttgtFld.ref());
214 
215  return ttgtFld;
216 }
217 
218 
219 template<class Type>
221 (
222  const VolField<Type>& srcFld,
223  const VolField<Type>& leftOverTgtFld,
224  const UList<wordRe>& tgtCuttingPatchNames
225 ) const
226 {
227  // Construct the result field
228  tmp<VolField<Type>> ttgtFld =
230  (
231  typedName("interpolate(" + srcFld.name() + ")"),
232  srcToTgt<Type>(srcFld.v(), leftOverTgtFld.v())(),
233  leftOverTgtFld.boundaryField()
234  );
235  typename VolField<Type>::Boundary& tgtBfld =
236  ttgtFld.ref().boundaryFieldRef();
237 
238  // Mapped patches
239  forAll(patchIndices(), i)
240  {
241  const label srcPatchi = patchIndices()[i].first();
242  const label tgtPatchi = patchIndices()[i].second();
243 
244  tgtBfld[tgtPatchi].map
245  (
246  leftOverTgtFld.boundaryField()[tgtPatchi],
248  );
249  tgtBfld[tgtPatchi].map
250  (
251  srcFld.boundaryField()[srcPatchi],
252  patchToPatchLeftOverFieldMapper(patchInterpolation(i))
253  );
254  }
255 
256  // Cutting patches. Set values to that of the internal cell field.
257  const labelHashSet tgtCuttingPatchIDs =
258  leftOverTgtFld.mesh().boundaryMesh().patchSet(tgtCuttingPatchNames);
259  forAllConstIter(labelHashSet, tgtCuttingPatchIDs, iter)
260  {
261  tgtBfld[iter.key()] == tgtBfld[iter.key()].patchInternalField();
262  }
263 
264  // Evaluate constraints
265  evaluateConstraintTypes(ttgtFld.ref());
266 
267  return ttgtFld;
268 }
269 
270 
271 template<class Type>
273 (
274  const VolInternalField<Type>& srcFld
275 ) const
276 {
277  tmp<VolInternalField<Type>> ttgtFld =
279  (
280  typedName("interpolate(" + srcFld.name() + ")"),
281  tgtMesh_,
282  srcFld.dimensions(),
283  cellsInterpolation().srcToTgt(srcFld)
284  );
285 
286  tgtCellsStabilisation().stabilise(ttgtFld.ref());
287 
288  return ttgtFld;
289 }
290 
291 
292 template<class Type>
294 (
295  const VolInternalField<Type>& srcFld,
296  const VolInternalField<Type>& leftOverTgtFld
297 ) const
298 {
299  return
301  (
302  typedName("interpolate(" + srcFld.name() + ")"),
303  tgtMesh_,
304  leftOverTgtFld.dimensions(),
305  cellsInterpolation().srcToTgt(srcFld, leftOverTgtFld)
306  );
307 }
308 
309 
310 template<class Type>
313 (
314  const SurfaceFieldBoundary<Type>& srcBfld
315 ) const
316 {
317  // Map all patch fields
318  PtrList<fvsPatchField<Type>> tgtPatchFields(tgtMesh_.boundary().size());
319  forAll(patchIndices(), i)
320  {
321  const label srcPatchi = patchIndices()[i].first();
322  const label tgtPatchi = patchIndices()[i].second();
323 
324  if (!tgtPatchFields.set(tgtPatchi))
325  {
326  tgtPatchFields.set
327  (
328  tgtPatchi,
330  (
331  srcBfld[srcPatchi],
332  tgtMesh_.boundary()[tgtPatchi],
335  (
336  patchInterpolation(i),
337  tgtPatchStabilisation(i)
338  )
339  )
340  );
341  }
342  }
343 
344  // Create any patch fields not explicitly mapped; e.g., constraints
345  forAll(tgtPatchFields, tgtPatchi)
346  {
347  if (!tgtPatchFields.set(tgtPatchi))
348  {
349  tgtPatchFields.set
350  (
351  tgtPatchi,
353  (
355  tgtMesh_.boundary()[tgtPatchi],
357  )
358  );
359  }
360  }
361 
362  return
364  (
366  (
367  tgtMesh_.boundary(),
369  tgtPatchFields
370  )
371  );
372 }
373 
374 
375 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:31
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
const word & name() const
Return name.
Definition: IOobject.H:310
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
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
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
tmp< VolField< Type > > srcToTgt(const VolField< Type > &srcFld) const
Interpolate a source vol field to the target with no left.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:170
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
Identity field mapper.
Patch-to-patch fieldMapper which retains values in the target field for parts of the patch that do no...
Patch-to-patch fieldMapper which fills values for non-overlapping parts of the target patch by normal...
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
Mapper which sets the field size. It does not actually map values.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
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))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:74
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:61
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176