nearWallFieldsTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 (
34 ) const
35 {
37 
38  HashTable<const vfType*> flds(obr_.lookupClass<vfType>());
39 
40  forAllConstIter(typename HashTable<const vfType*>, flds, iter)
41  {
42  const vfType& fld = *iter();
43 
44  if (fieldMap_.found(fld.name()))
45  {
46  const word& sampleFldName = fieldMap_[fld.name()];
47 
48  if (obr_.found(sampleFldName))
49  {
51  << " a field named " << sampleFldName
52  << " already exists on the mesh"
53  << endl;
54  }
55  else
56  {
57  label sz = sflds.size();
58  sflds.setSize(sz+1);
59 
60  IOobject io(fld);
61  io.readOpt() = IOobject::NO_READ;
62  io.writeOpt() = IOobject::NO_WRITE;
63  io.rename(sampleFldName);
64 
65  sflds.set(sz, new vfType(io, fld));
66 
67  if (log_) Info
68  << " created " << sflds[sz].name() << " to sample "
69  << fld.name() << endl;
70  }
71  }
72  }
73 }
74 
75 
76 template<class Type>
78 (
79  const interpolationCellPoint<Type>& interpolator,
81 ) const
82 {
83  // Construct flat fields for all patch faces to be sampled
84  Field<Type> sampledValues(getPatchDataMapPtr_().constructSize());
85 
86  forAll(cellToWalls_, cellI)
87  {
88  const labelList& cData = cellToWalls_[cellI];
89 
90  forAll(cData, i)
91  {
92  const point& samplePt = cellToSamples_[cellI][i];
93  sampledValues[cData[i]] = interpolator.interpolate(samplePt, cellI);
94  }
95  }
96 
97  // Send back sampled values to patch faces
98  getPatchDataMapPtr_().reverseDistribute
99  (
100  getPatchDataMapPtr_().constructSize(),
101  sampledValues
102  );
103 
104  // Pick up data
105  label nPatchFaces = 0;
106  forAllConstIter(labelHashSet, patchSet_, iter)
107  {
108  label patchI = iter.key();
109 
110  fvPatchField<Type>& pfld = fld.boundaryField()[patchI];
111 
112  Field<Type> newFld(pfld.size());
113  forAll(pfld, i)
114  {
115  newFld[i] = sampledValues[nPatchFaces++];
116  }
117 
118  pfld == newFld;
119  }
120 }
121 
122 
123 template<class Type>
125 (
127 ) const
128 {
130 
131  forAll(sflds, i)
132  {
133  const word& fldName = reverseFieldMap_[sflds[i].name()];
134  const vfType& fld = obr_.lookupObject<vfType>(fldName);
135 
136  // Take over internal and boundary values
137  sflds[i] == fld;
138 
139  // Construct interpolation method
140  interpolationCellPoint<Type> interpolator(fld);
141 
142  // Override sampled values
143  sampleBoundaryField(interpolator, sflds[i]);
144  }
145 }
146 
147 
148 // ************************************************************************* //
Foam::fvPatchField< Type >
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::IOobject::writeOpt
writeOption writeOpt() const
Definition: IOobject.H:327
Foam::HashSet< label, Hash< label > >
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::nearWallFields::createFields
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
Definition: nearWallFieldsTemplates.C:32
nearWallFields.H
Foam::label
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
Foam::nearWallFields::sampleFields
void sampleFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
Definition: nearWallFieldsTemplates.C:125
Foam::Field< Type >
Foam::IOobject::readOpt
readOption readOpt() const
Definition: IOobject.H:317
Foam::interpolationCellPoint
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Definition: interpolationCellPoint.H:48
Foam::Info
messageStream Info
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::interpolationCellPoint::interpolate
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
Definition: interpolationCellPointI.H:30
Foam::IOobject::rename
virtual void rename(const word &newName)
Rename.
Definition: IOobject.H:297
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::nearWallFields::sampleBoundaryField
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, GeometricField< Type, fvPatchField, volMesh > &fld) const
Override boundary fields with sampled values.
Definition: nearWallFieldsTemplates.C:78
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259