faFieldDecomposerTemplates.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) 2016-2017 Wikki Ltd
9 Copyright (C) 2021-2022 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 "faFieldDecomposer.H"
30#include "GeometricField.H"
31#include "IOobjectList.H"
34
35// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36
37template<class Type>
40(
42) const
43{
44 // Create and map the internal field values
45 Field<Type> internalField(field.internalField(), faceAddressing_);
46
47 // Create and map the patch field values
48 PtrList<faPatchField<Type>> patchFields(boundaryAddressing_.size());
49
50 forAll(boundaryAddressing_, patchi)
51 {
52 const label oldPatchi = boundaryAddressing_[patchi];
53
54 if (oldPatchi >= 0)
55 {
56 patchFields.set
57 (
58 patchi,
60 (
61 field.boundaryField()[oldPatchi],
62 procMesh_.boundary()[patchi],
64 patchFieldDecomposerPtrs_[patchi]
65 )
66 );
67 }
68 else
69 {
70 patchFields.set
71 (
72 patchi,
74 (
75 procMesh_.boundary()[patchi],
78 (
79 field.internalField(),
80 processorAreaPatchFieldDecomposerPtrs_[patchi]
81 )
82 )
83 );
84 }
85 }
86
87 // Create the field for the processor
88 return
90 (
92 (
93 field.name(),
94 procMesh_.thisDb().time().timeName(),
95 procMesh_.thisDb(),
98 ),
99 procMesh_,
100 field.dimensions(),
101 internalField,
102 patchFields
103 );
104}
105
106
107template<class Type>
110(
112) const
113{
114 labelList mapAddr
115 (
117 (
118 edgeAddressing_,
119 procMesh_.nInternalEdges()
120 )
121 );
122 forAll(mapAddr, i)
123 {
124 mapAddr[i] -= 1;
125 }
126
127 // Create and map the internal field values
128 Field<Type> internalField
129 (
130 field.internalField(),
131 mapAddr
132 );
133
134 // Problem with addressing when a processor patch picks up both internal
135 // edges and edges from cyclic boundaries. This is a bit of a hack, but
136 // I cannot find a better solution without making the internal storage
137 // mechanism for edgeFields correspond to the one of edges in polyMesh
138 // (i.e. using slices)
139 Field<Type> allEdgeField(field.mesh().nEdges());
140
141 forAll(field.internalField(), i)
142 {
143 allEdgeField[i] = field.internalField()[i];
144 }
145
146 forAll(field.boundaryField(), patchi)
147 {
148 const Field<Type>& p = field.boundaryField()[patchi];
149
150 const label patchStart = field.mesh().boundary()[patchi].start();
151
152 forAll(p, i)
153 {
154 allEdgeField[patchStart + i] = p[i];
155 }
156 }
157
158 // Create and map the patch field values
159 PtrList<faePatchField<Type>> patchFields(boundaryAddressing_.size());
160
161 forAll(boundaryAddressing_, patchi)
162 {
163 const label oldPatchi = boundaryAddressing_[patchi];
164
165 if (oldPatchi >= 0)
166 {
167 patchFields.set
168 (
169 patchi,
171 (
172 field.boundaryField()[oldPatchi],
173 procMesh_.boundary()[patchi],
175 patchFieldDecomposerPtrs_[patchi]
176 )
177 );
178 }
179 else
180 {
181 patchFields.set
182 (
183 patchi,
185 (
186 procMesh_.boundary()[patchi],
189 (
190 allEdgeField,
191 processorEdgePatchFieldDecomposerPtrs_[patchi]
192 )
193 )
194 );
195 }
196 }
197
198 // Create the field for the processor
199 return
201 (
203 (
204 field.name(),
205 procMesh_.thisDb().time().timeName(),
206 procMesh_.thisDb(),
209 ),
210 procMesh_,
211 field.dimensions(),
212 internalField,
213 patchFields
214 );
215}
216
217
218template<class GeoField>
220(
222) const
223{
224 forAll(fields, fieldi)
225 {
226 decomposeField(fields[fieldi])().write();
227 }
228}
229
230
231// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
A List obtained as a section of another List.
Definition: SubList.H:70
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
tmp< GeometricField< Type, faPatchField, areaMesh > > decomposeField(const GeometricField< Type, faPatchField, areaMesh > &field) const
Decompose area field.
void decomposeFields(const PtrList< GeoField > &fields) const
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:38
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:697
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: faPatchField.H:82
faePatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cove...
Definition: faePatchField.H:83
const Time & time() const noexcept
Return time registry.
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
A class for managing temporary objects.
Definition: tmp.H:65
volScalarField & p
rDeltaTY field()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
runTime write()
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333