Loading...
Searching...
No Matches
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
41 for (const VolFieldType& fld : obr_.csorted<VolFieldType>())
42 {
43 const auto fieldMapIter = fieldMap_.cfind(fld.name());
44
45 if (fieldMapIter.good())
46 {
47 const word& sampleFldName = fieldMapIter.val();
48
49 if (obr_.found(sampleFldName))
50 {
52 << " a field named " << sampleFldName
53 << " already exists on the mesh"
54 << endl;
55 }
56 else
57 {
61 io.rename(sampleFldName);
62
63 // Override bc to be calculated
64 const label newFieldi = sflds.size();
65 sflds.resize(newFieldi+1);
66
67 sflds.set
68 (
69 newFieldi,
70 new VolFieldType
71 (
72 io,
73 fld,
76 )
77 );
78
79 Log << " created " << io.name()
80 << " to sample " << fld.name() << endl;
81 }
82 }
83 }
84}
85
86
87template<class Type>
89(
90 const interpolationCellPoint<Type>& interpolator,
91 GeometricField<Type, fvPatchField, volMesh>& fld
92) const
93{
94 // Construct flat fields for all patch faces to be sampled
95 Field<Type> sampledValues(getPatchDataMapPtr_().constructSize());
96
97 forAll(cellToWalls_, celli)
98 {
99 const labelList& cData = cellToWalls_[celli];
100
101 forAll(cData, i)
102 {
103 const point& samplePt = cellToSamples_[celli][i];
104 sampledValues[cData[i]] = interpolator.interpolate(samplePt, celli);
105 }
106 }
107
108 // Send back sampled values to patch faces
109 getPatchDataMapPtr_().reverseDistribute
110 (
111 getPatchDataMapPtr_().constructSize(),
112 sampledValues
113 );
114
115 typename GeometricField<Type, fvPatchField, volMesh>::
116 Boundary& fldBf = fld.boundaryFieldRef();
117
118 // Pick up data
119 label nPatchFaces = 0;
120 for (const label patchi : patchIDs_)
121 {
122 fvPatchField<Type>& pfld = fldBf[patchi];
123
124 Field<Type> newFld(pfld.size());
125 forAll(pfld, i)
126 {
127 newFld[i] = sampledValues[nPatchFaces++];
128 }
130 pfld == newFld;
131 }
132}
133
134
135template<class Type>
137(
139) const
140{
142
143 forAll(sflds, i)
144 {
145 const word& fldName = reverseFieldMap_[sflds[i].name()];
146 const VolFieldType& fld = obr_.lookupObject<VolFieldType>(fldName);
147
148 // Take over internal and boundary values
149 sflds[i] == fld;
150
151 // Construct interpolation method
152 interpolationCellPoint<Type> interpolator(fld);
153
154 // Override sampled values
155 sampleBoundaryField(interpolator, sflds[i]);
156 }
157}
158
159
160// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
HashTable< word > reverseFieldMap_
From resulting back to original field.
autoPtr< mapDistribute > getPatchDataMapPtr_
Map from cell based data back to patch based data.
List< List< point > > cellToSamples_
From cell to tracked end point.
labelListList cellToWalls_
From cell to seed patch faces.
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
Create the local fields.
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
Sample the fields.
labelList patchIDs_
Patches to sample.
HashTable< word > fieldMap_
From original field to sampled result.
const objectRegistry & obr_
Reference to the region objectRegistry.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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.
A class for handling words, derived from Foam::string.
Definition word.H:66
const auto & io
#define WarningInFunction
Report a warning using Foam::Warning.
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299