conformedFvsPatchField.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) 2011-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 "conformedFvsPatchField.H"
27 #include "fvMeshStitcherTools.H"
28 #include "nonConformalBoundary.H"
29 #include "nonConformalFvPatch.H"
31 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * Private Constructors * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvPatch& p,
39  const DimensionedField<Type, surfaceMesh>& iF,
40  autoPtr<fvsPatchField<Type>>&& origFieldPtr,
41  autoPtr<calculatedFvsPatchField<Type>>&& ncFieldPtr
42 )
43 :
44  fvsPatchField<Type>(p, iF),
45  origFieldPtr_(origFieldPtr),
46  ncFieldPtr_(ncFieldPtr)
47 {}
48 
49 
50 // * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * * //
51 
52 template<class Type>
54 (
55  typename SurfaceField<Type>::Boundary& bF
56 )
57 {
58  const DimensionedField<Type, surfaceMesh>& iF = bF[0].internalField();
59 
60  const fvBoundaryMesh& fvbm = iF.mesh().boundary();
61 
62  const labelList origPatchIndices =
64 
65  // Evaluate the conformed orig and non-conformal boundary fields
66  const typename SurfaceField<Type>::Boundary origBf
67  (
70  );
71  const typename SurfaceField<Type>::Boundary ncBf
72  (
75  );
76 
77  // Replace every original patch field with a conformed patch field
78  // containing the conformed orig and non-conformal fields
79  forAll(origPatchIndices, i)
80  {
81  const label origPatchi = origPatchIndices[i];
82  const fvPatch& origFvp = fvbm[origPatchi];
83 
85  (
87  (
88  origFvp,
89  iF,
90  bF.set(origPatchi, nullptr),
92  (
93  new calculatedFvsPatchField<Type>(origFvp, iF)
94  )
95  )
96  );
97 
98  pF() == origBf[origPatchi];
99  pF->origFieldPtr_() == origBf[origPatchi];
100  pF->ncFieldPtr_() == ncBf[origPatchi];
101 
102  bF.set(origPatchi, pF.ptr());
103  }
104 }
105 
106 
107 template<class Type>
109 (
110  typename SurfaceField<Type>::Boundary& bF
111 )
112 {
113  const DimensionedField<Type, surfaceMesh>& iF = bF[0].internalField();
114 
115  const fvBoundaryMesh& fvbm = iF.mesh().boundary();
116 
117  const labelList origPatchIndices =
119 
120  // If this field does not contain conformed patch fields then it was
121  // created during the mesh change process, between un-stitch and stitch.
122  // The only field known for which this can happen is the Crank-Nicolson
123  // mesh flux. What we create here doesn't really matter all that much as
124  // subsequent steps will update this flux with the actual mesh flux, and
125  // there's plenty of handling elsewhere that maps that correctly. So,
126  // to get something vaguely sensible for the first time-step in which this
127  // flux exists, transfer values from the original faces to the
128  // non-conformal faces in proportion with their area ratios. This, at
129  // least, ensures that the total flux is preserved.
130  forAll(origPatchIndices, i)
131  {
132  const label origPatchi = origPatchIndices[i];
133 
134  if (!isA<conformedFvsPatchField<Type>>(bF[origPatchi]))
135  {
137  return;
138  }
139  }
140 
141  // Extract the conformed orig and non-conformal boundary fields from
142  // the stored conformed patch fields
143  PtrList<fvsPatchField<Type>> origPFs(fvbm.size());
144  PtrList<fvsPatchField<Type>> ncPFs(fvbm.size());
145  forAll(origPatchIndices, i)
146  {
147  const label origPatchi = origPatchIndices[i];
148 
150  refCast<conformedFvsPatchField<Type>>(bF[origPatchi]);
151 
152  // If the mesh has topo-changed then maintained surface fields should
153  // have been mapped or re-interpolated. So, copy the value from the
154  // base field into the original field.
155  if (iF.mesh().topoChanged())
156  {
157  cpF.origFieldPtr_() = bF[origPatchi];
158  }
159 
160  origPFs.set(origPatchi, cpF.origFieldPtr_.ptr());
161  ncPFs.set(origPatchi, cpF.ncFieldPtr_.ptr());
162  }
163  forAll(origPFs, patchi)
164  {
165  if (origPFs.set(patchi)) continue;
166 
167  origPFs.set(patchi, bF.set(patchi, nullptr));
168  ncPFs.set
169  (
170  patchi,
172  (
174  fvbm[patchi],
175  iF
176  )
177  );
178  }
179 
180  // If the mesh has topo-changed then just use the original parts and leave
181  // the non-conformal parts unset
182  if (iF.mesh().topoChanged())
183  {
184  bF.transfer(origPFs);
185  }
186  // If the mesh has not topo-changed, then combine the conformed boundary
187  // fields to create the non-conformal boundary field
188  else
189  {
190  typename SurfaceField<Type>::Boundary origAndNcBf
191  (
192  iF,
194  (
195  typename SurfaceField<Type>::Boundary(fvbm, iF, ncPFs),
196  typename SurfaceField<Type>::Boundary(fvbm, iF, origPFs)
197  )
198  );
199 
200  bF.transfer(origAndNcBf);
201 
203  }
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
208 
209 template<class Type>
211 (
212  const fvPatch& p,
214 )
215 :
216  fvsPatchField<Type>(p, iF)
217 {
219 }
220 
221 
222 template<class Type>
224 (
225  const fvPatch& p,
227  const dictionary& dict
228 )
229 :
230  fvsPatchField<Type>(p, iF, dict),
231  origFieldPtr_
232  (
233  fvsPatchField<Type>::New(p, iF, dict.subDict("origField")).ptr()
234  ),
235  ncFieldPtr_
236  (
237  new calculatedFvsPatchField<Type>(p, iF, dict.subDict("ncField"))
238  )
239 {}
240 
241 
242 template<class Type>
244 (
245  const conformedFvsPatchField<Type>& ptf,
246  const fvPatch& p,
248  const fieldMapper& mapper
249 )
250 :
251  fvsPatchField<Type>(ptf, p, iF, mapper),
252  origFieldPtr_
253  (
254  fvsPatchField<Type>::New(ptf.origFieldPtr_(), p, iF, mapper).ptr()
255  ),
256  ncFieldPtr_
257  (
258  new calculatedFvsPatchField<Type>
259  (
260  ptf.ncFieldPtr_(),
261  p,
262  iF,
263  mapper
264  )
265  )
266 {}
267 
268 
269 template<class Type>
271 (
272  const conformedFvsPatchField<Type>& ptf,
274 )
275 :
276  fvsPatchField<Type>(ptf, iF),
277  origFieldPtr_(ptf.origFieldPtr_->clone(iF).ptr()),
278  ncFieldPtr_(new calculatedFvsPatchField<Type>(ptf.ncFieldPtr_(), iF))
279 {}
280 
281 
282 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
283 
284 template<class Type>
286 (
287  const fvsPatchField<Type>& ptf,
288  const fieldMapper& mapper
289 )
290 {
291  fvsPatchField<Type>::map(ptf, mapper);
292 
294  {
295  const conformedFvsPatchField<Type>& cptf =
296  refCast<const conformedFvsPatchField<Type>>(ptf);
297 
298  origFieldPtr_->map(cptf.origFieldPtr_(), mapper);
299  ncFieldPtr_->map(cptf.ncFieldPtr_(), mapper);
300  }
301  else
302  {
303  origFieldPtr_->reset(ptf);
304  ncFieldPtr_() == origFieldPtr_();
305  }
306 }
307 
308 
309 template<class Type>
311 {
313 
315  {
316  const conformedFvsPatchField<Type>& cptf =
317  refCast<const conformedFvsPatchField<Type>>(ptf);
318 
319  origFieldPtr_->reset(cptf.origFieldPtr_());
320  ncFieldPtr_->reset(cptf.ncFieldPtr_());
321  }
322  else
323  {
324  origFieldPtr_->reset(ptf);
325  ncFieldPtr_() == origFieldPtr_();
326  }
327 }
328 
329 
330 template<class Type>
332 {
334  writeEntry(os, "value", *this);
335 
336  writeKeyword(os, "origField") << nl;
337  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
338  origFieldPtr_->write(os);
339  os << decrIndent << indent << token::END_BLOCK << nl;
340 
341  writeKeyword(os, "ncField") << nl;
342  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
343  ncFieldPtr_->write(os);
344  os << decrIndent << indent << token::END_BLOCK << nl;
345 }
346 
347 
348 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
static nonConformalBoundary & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
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
void transfer(PtrList< T > &)
Transfer the contents of the argument PtrList into this PtrList.
Definition: PtrList.C:213
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
T * ptr()
Return object pointer for reuse.
Definition: autoPtrI.H:90
Foam::calculatedFvsPatchField.
This surface field boundary condition holds data from both the original faces and any associated non-...
virtual void write(Ostream &) const
Write.
static void unconform(typename SurfaceField< Type >::Boundary &bF)
Un-conform the given boundary field.
virtual void reset(const fvsPatchField< Type > &)
Reset the fvsPatchField to the given fvsPatchField.
virtual void map(const fvsPatchField< Type > &, const fieldMapper &)
Map the given fvsPatchField onto this fvsPatchField.
static void conform(typename SurfaceField< Type >::Boundary &bF)
Conform the given boundary field.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Foam::fvBoundaryMesh.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:86
virtual void write(Ostream &) const
Write.
virtual void reset(const fvsPatchField< Type > &)
Reset the fvsPatchField to the given fvsPatchField.
virtual void map(const fvsPatchField< Type > &, const fieldMapper &)
Map the given fvsPatchField onto this fvsPatchField.
labelList allOrigPatchIndices() const
Return a list of the orig patch indices.
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
label patchi
tmp< SurfaceFieldBoundary< Type > > conformedNcBoundaryField(const SurfaceFieldBoundary< Type > &fieldb)
Extract the non-conformal parts of the boundary field and store it on the.
tmp< SurfaceFieldBoundary< Type > > synchronisedBoundaryField(const SurfaceFieldBoundary< Type > &fieldb, const bool flip, const scalar ownerWeight, const scalar neighbourWeight)
Synchronise the boundary field by combining corresponding.
tmp< SurfaceFieldBoundary< Type > > unconformedBoundaryField(const SurfaceFieldBoundary< Type > &ncFieldb, const SurfaceFieldBoundary< Type > &origFieldb)
Combine non-conformal and original parts of the boundary field from the.
tmp< SurfaceFieldBoundary< Type > > conformedOrigBoundaryField(const SurfaceFieldBoundary< Type > &fieldb)
Extract the original parts of the boundary field and store it.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
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
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:171
T clone(const T &t)
Definition: List.H:55
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
static const char nl
Definition: Ostream.H:267
dictionary dict
volScalarField & p
Foam::surfaceFields.