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-2023 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"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
34 void Foam::fvMeshToFvMesh::evaluateConstraintTypes(VolField<Type>& fld)
35 {
36  typename VolField<Type>::Boundary& fldBf = fld.boundaryFieldRef();
37 
38  if
39  (
42  )
43  {
44  label nReq = Pstream::nRequests();
45 
46  forAll(fldBf, patchi)
47  {
48  fvPatchField<Type>& tgtField = fldBf[patchi];
49 
50  if
51  (
52  tgtField.type() == tgtField.patch().patch().type()
53  && polyPatch::constraintType(tgtField.patch().patch().type())
54  )
55  {
56  tgtField.initEvaluate(Pstream::defaultCommsType);
57  }
58  }
59 
60  // Block for any outstanding requests
61  if
62  (
65  )
66  {
68  }
69 
70  forAll(fldBf, patchi)
71  {
72  fvPatchField<Type>& tgtField = fldBf[patchi];
73 
74  if
75  (
76  tgtField.type() == tgtField.patch().patch().type()
77  && polyPatch::constraintType(tgtField.patch().patch().type())
78  )
79  {
80  tgtField.evaluate(Pstream::defaultCommsType);
81  }
82  }
83  }
85  {
86  const lduSchedule& patchSchedule =
87  fld.mesh().globalData().patchSchedule();
88 
89  forAll(patchSchedule, patchEvali)
90  {
91  label patchi = patchSchedule[patchEvali].patch;
92  fvPatchField<Type>& tgtField = fldBf[patchi];
93 
94  if
95  (
96  tgtField.type() == tgtField.patch().patch().type()
97  && polyPatch::constraintType(tgtField.patch().patch().type())
98  )
99  {
100  if (patchSchedule[patchEvali].init)
101  {
102  tgtField.initEvaluate(Pstream::commsTypes::scheduled);
103  }
104  else
105  {
106  tgtField.evaluate(Pstream::commsTypes::scheduled);
107  }
108  }
109  }
110  }
111 }
112 
113 
114 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
115 
116 template<class Type>
118 (
119  const VolField<Type>& srcFld
120 ) const
121 {
122  const fvMesh& tgtMesh = static_cast<const fvMesh&>(meshToMesh::tgtMesh());
123 
124  // Construct target patch fields as copies of source patch fields, but do
125  // not map values yet
126  PtrList<fvPatchField<Type>> tgtPatchFields(tgtMesh.boundary().size());
127  forAll(patchIDs(), i)
128  {
129  const label srcPatchi = patchIDs()[i].first();
130  const label tgtPatchi = patchIDs()[i].second();
131 
132  if (!tgtPatchFields.set(tgtPatchi))
133  {
134  tgtPatchFields.set
135  (
136  tgtPatchi,
138  (
139  srcFld.boundaryField()[srcPatchi],
140  tgtMesh.boundary()[tgtPatchi],
143  (
144  labelList(tgtMesh.boundary()[tgtPatchi].size(), -1)
145  )
146  )
147  );
148  }
149  }
150 
151  // Create calculated patch fields for any unset target patches. Use
152  // fvPatchField<Type>::New to construct these, rather than using the
153  // calculated constructor directly, so that constraints are maintained.
154  labelList tgtPatchFieldIsUnMapped(tgtPatchFields.size(), false);
155  forAll(tgtPatchFields, tgtPatchi)
156  {
157  if (!tgtPatchFields.set(tgtPatchi))
158  {
159  tgtPatchFields.set
160  (
161  tgtPatchi,
163  (
165  tgtMesh.boundary()[tgtPatchi],
167  )
168  );
169 
170  tgtPatchFieldIsUnMapped[tgtPatchi] =
172  (
173  tgtMesh.boundary()[tgtPatchi].patch().type()
174  );
175  }
176  }
177 
178  // Construct the result field
179  tmp<VolField<Type>> ttgtFld =
181  (
182  typedName("interpolate(" + srcFld.name() + ")"),
183  srcToTgt<Type>(srcFld.internalField())(),
184  tgtPatchFields
185  );
186  typename VolField<Type>::Boundary& tgtBfld =
187  ttgtFld.ref().boundaryFieldRef();
188 
189  // Mapped patches
190  forAll(patchIDs(), i)
191  {
192  const label srcPatchi = patchIDs()[i].first();
193  const label tgtPatchi = patchIDs()[i].second();
194 
195  tgtBfld[tgtPatchi].map
196  (
197  srcFld.boundaryField()[srcPatchi],
199  (
200  patchInterpolation(i),
201  tgtPatchStabilisation(i)
202  )
203  );
204  }
205 
206  // Un-mapped patches. Set values to that of the internal cell field.
207  forAll(tgtBfld, patchi)
208  {
209  if (tgtPatchFieldIsUnMapped[patchi])
210  {
211  fvPatchField<Type>& tgtPfld = tgtBfld[patchi];
212  tgtPfld == tgtPfld.patchInternalField();
213  }
214  }
215 
216  // Evaluate constraints
217  evaluateConstraintTypes(ttgtFld.ref());
218 
219  return ttgtFld;
220 }
221 
222 
223 template<class Type>
225 (
226  const VolField<Type>& srcFld,
227  const VolField<Type>& leftOverTgtFld,
228  const UList<wordRe>& tgtCuttingPatchNames
229 ) const
230 {
231  // Construct the result field
232  tmp<VolField<Type>> ttgtFld =
234  (
235  typedName("interpolate(" + srcFld.name() + ")"),
236  srcToTgt<Type>(srcFld.v(), leftOverTgtFld.v())(),
237  leftOverTgtFld.boundaryField()
238  );
239  typename VolField<Type>::Boundary& tgtBfld =
240  ttgtFld.ref().boundaryFieldRef();
241 
242  // Mapped patches
243  forAll(patchIDs(), i)
244  {
245  const label srcPatchi = patchIDs()[i].first();
246  const label tgtPatchi = patchIDs()[i].second();
247 
248  tgtBfld[tgtPatchi].map
249  (
250  leftOverTgtFld.boundaryField()[tgtPatchi],
252  );
253  tgtBfld[tgtPatchi].map
254  (
255  srcFld.boundaryField()[srcPatchi],
256  patchToPatchLeftOverFvPatchFieldMapper(patchInterpolation(i))
257  );
258  }
259 
260  // Cutting patches. Set values to that of the internal cell field.
261  const labelHashSet tgtCuttingPatchIDs =
262  leftOverTgtFld.mesh().boundaryMesh().patchSet(tgtCuttingPatchNames);
263  forAllConstIter(labelHashSet, tgtCuttingPatchIDs, iter)
264  {
265  tgtBfld[iter.key()] == tgtBfld[iter.key()].patchInternalField();
266  }
267 
268  // Evaluate constraints
269  evaluateConstraintTypes(ttgtFld.ref());
270 
271  return ttgtFld;
272 }
273 
274 
275 template<class Type>
277 (
278  const VolInternalField<Type>& srcFld
279 ) const
280 {
281  tmp<VolInternalField<Type>> ttgtFld =
283  (
284  typedName("interpolate(" + srcFld.name() + ")"),
285  static_cast<const fvMesh&>(meshToMesh::tgtMesh()),
286  srcFld.dimensions(),
287  cellsInterpolation().srcToTgt(srcFld)
288  );
289 
290  tgtCellsStabilisation().stabilise(ttgtFld.ref());
291 
292  return ttgtFld;
293 }
294 
295 
296 template<class Type>
298 (
299  const VolInternalField<Type>& srcFld,
300  const VolInternalField<Type>& leftOverTgtFld
301 ) const
302 {
303  return
305  (
306  typedName("interpolate(" + srcFld.name() + ")"),
307  static_cast<const fvMesh&>(meshToMesh::tgtMesh()),
308  leftOverTgtFld.dimensions(),
309  cellsInterpolation().srcToTgt(srcFld, leftOverTgtFld)
310  );
311 }
312 
313 
314 // ************************************************************************* //
#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.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
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.
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:65
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.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
const polyMesh & tgtMesh() const
Return const access to the target mesh.
Definition: meshToMeshI.H:133
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
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)
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:58
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149