nearWallFieldsTemplates.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-2022 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 "nearWallFields.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  PtrList<VolField<Type>>& sflds
34 ) const
35 {
37 
38  forAllConstIter(typename HashTable<const VolField<Type>*>, flds, iter)
39  {
40  const VolField<Type>& fld = *iter();
41 
42  if (fieldMap_.found(fld.name()))
43  {
44  const word& sampleFldName = fieldMap_[fld.name()];
45 
46  if (obr_.found(sampleFldName))
47  {
48  Log << " a field " << sampleFldName
49  << " already exists on the mesh."
50  << endl;
51  }
52  else
53  {
54  label sz = sflds.size();
55  sflds.setSize(sz+1);
56 
57  sflds.set
58  (
59  sz,
60  new VolField<Type>
61  (
62  IOobject
63  (
64  sampleFldName,
65  time_.name(),
66  mesh_
67  ),
68  fld,
69  calculatedFvPatchScalarField::typeName
70  )
71  );
72 
73  Log << " created " << sflds[sz].name()
74  << " to sample " << fld.name() << endl;
75  }
76  }
77  }
78 }
79 
80 
81 template<class Type>
83 (
84  const interpolationCellPoint<Type>& interpolator,
86 ) const
87 {
88  // Construct flat fields for all patch faces to be sampled
89  Field<Type> sampledValues(getPatchDataMapPtr_().constructSize());
90 
91  forAll(cellToWalls_, celli)
92  {
93  const labelList& cData = cellToWalls_[celli];
94 
95  forAll(cData, i)
96  {
97  const point& samplePt = cellToSamples_[celli][i];
98  sampledValues[cData[i]] = interpolator.interpolate(samplePt, celli);
99  }
100  }
101 
102  // Send back sampled values to patch faces
103  getPatchDataMapPtr_().reverseDistribute
104  (
105  getPatchDataMapPtr_().constructSize(),
106  sampledValues
107  );
108 
109  typename VolField<Type>::
110  Boundary& fldBf = fld.boundaryFieldRef();
111 
112  // Pick up data
113  label nPatchFaces = 0;
114  forAllConstIter(labelHashSet, patchSet_, iter)
115  {
116  label patchi = iter.key();
117 
118  fvPatchField<Type>& pfld = fldBf[patchi];
119 
120  Field<Type> newFld(pfld.size());
121  forAll(pfld, i)
122  {
123  newFld[i] = sampledValues[nPatchFaces++];
124  }
125 
126  pfld == newFld;
127  }
128 }
129 
130 
131 template<class Type>
133 (
134  PtrList<VolField<Type>>& sflds
135 ) const
136 {
137  forAll(sflds, i)
138  {
139  const word& fldName = reverseFieldMap_[sflds[i].name()];
140  const VolField<Type>& fld = obr_.lookupObject<VolField<Type>>(fldName);
141 
142  // Take over internal and boundary values
143  sflds[i] == fld;
144 
145  // Construct interpolation method
146  interpolationCellPoint<Type> interpolator(fld);
147 
148  // Override sampled values
149  sampleBoundaryField(interpolator, sflds[i]);
150  }
151 }
152 
153 
154 // ************************************************************************* //
#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
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricField class.
An STL-conforming hash table.
Definition: HashTable.H:127
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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
const word & name() const
Return const reference to name.
const Time & time_
Reference to time.
const fvMesh & mesh_
Reference to the fvMesh.
void sampleFields(PtrList< VolField< Type >> &) const
void createFields(PtrList< VolField< Type >> &) const
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, VolField< Type > &fld) const
Override boundary fields with sampled values.
HashTable< word > fieldMap_
From original field to sampled result.
const objectRegistry & obr_
Reference to the region objectRegistry.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
A class for handling words, derived from string.
Definition: word.H:62
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))
#define Log
Report write to Foam::Info if the local log switch is true.
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 & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251