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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "nearWallFields.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33template<class Type>
35(
37) const
38{
40
42
43 forAllConstIters(flds, iter)
44 {
45 const VolFieldType& fld = *(iter.val());
46
47 if (fieldMap_.found(fld.name()))
48 {
49 const word& sampleFldName = fieldMap_[fld.name()];
50
51 if (obr_.found(sampleFldName))
52 {
54 << " a field named " << sampleFldName
55 << " already exists on the mesh"
56 << endl;
57 }
58 else
59 {
60 label sz = sflds.size();
61 sflds.setSize(sz+1);
62
66
67 io.rename(sampleFldName);
68
69 // Override bc to be calculated
70 sflds.set
71 (
72 sz,
73 new VolFieldType
74 (
75 io,
76 fld,
77 patchSet_.toc(),
79 )
80 );
81
82 Log << " created " << sflds[sz].name()
83 << " to sample " << fld.name() << endl;
84 }
85 }
86 }
87}
88
89
90template<class Type>
92(
93 const interpolationCellPoint<Type>& interpolator,
95) const
96{
97 // Construct flat fields for all patch faces to be sampled
98 Field<Type> sampledValues(getPatchDataMapPtr_().constructSize());
99
100 forAll(cellToWalls_, celli)
101 {
102 const labelList& cData = cellToWalls_[celli];
103
104 forAll(cData, i)
105 {
106 const point& samplePt = cellToSamples_[celli][i];
107 sampledValues[cData[i]] = interpolator.interpolate(samplePt, celli);
108 }
109 }
110
111 // Send back sampled values to patch faces
112 getPatchDataMapPtr_().reverseDistribute
113 (
114 getPatchDataMapPtr_().constructSize(),
115 sampledValues
116 );
117
119 Boundary& fldBf = fld.boundaryFieldRef();
120
121 // Pick up data
122 label nPatchFaces = 0;
123 for (const label patchi : patchSet_)
124 {
125 fvPatchField<Type>& pfld = fldBf[patchi];
126
127 Field<Type> newFld(pfld.size());
128 forAll(pfld, i)
129 {
130 newFld[i] = sampledValues[nPatchFaces++];
131 }
132
133 pfld == newFld;
134 }
135}
136
137
138template<class Type>
140(
142) const
143{
145
146 forAll(sflds, i)
147 {
148 const word& fldName = reverseFieldMap_[sflds[i].name()];
149 const VolFieldType& fld = obr_.lookupObject<VolFieldType>(fldName);
150
151 // Take over internal and boundary values
152 sflds[i] == fld;
153
154 // Construct interpolation method
155 interpolationCellPoint<Type> interpolator(fld);
156
157 // Override sampled values
158 sampleBoundaryField(interpolator, sflds[i]);
159 }
160}
161
162
163// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:122
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
virtual void rename(const word &newName)
Rename the object.
Definition: IOobject.H:497
writeOption writeOpt() const noexcept
The write option.
Definition: IOobjectI.H:179
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, GeometricField< Type, fvPatchField, volMesh > &fld) const
Override boundary fields with sampled values.
void sampleFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
labelHashSet patchSet_
Patches to sample.
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:82
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
Return all objects with a class satisfying isA<Type>
bool found(const word &name, const bool recursive=false) const
Can the regIOobject object be found (by name).
A class for handling words, derived from Foam::string.
Definition: word.H:68
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
#define WarningInFunction
Report a warning using Foam::Warning.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278