pointFieldDecomposerTemplates.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) 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
31
32// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33
34template<class Type>
37(
39) const
40{
41 // Create and map the internal field values
42 Field<Type> internalField(field.primitiveField(), pointAddressing_);
43
44 // Create a list of pointers for the patchFields
45 PtrList<pointPatchField<Type>> patchFields(boundaryAddressing_.size());
46
47 // Create and map the patch field values
48 forAll(boundaryAddressing_, patchi)
49 {
50 if (patchFieldDecomposerPtrs_.set(patchi))
51 {
52 patchFields.set
53 (
54 patchi,
56 (
57 field.boundaryField()[boundaryAddressing_[patchi]],
58 procMesh_.boundary()[patchi],
60 patchFieldDecomposerPtrs_[patchi]
61 )
62 );
63 }
64 else
65 {
66 patchFields.set
67 (
68 patchi,
70 (
71 procMesh_.boundary()[patchi],
73 )
74 );
75 }
76 }
77
78 // Create the field for the processor
79 return
81 (
83 (
84 field.name(),
85 procMesh_.thisDb().time().timeName(),
86 procMesh_.thisDb(),
89 false
90 ),
91 procMesh_,
92 field.dimensions(),
93 internalField,
94 patchFields
95 );
96}
97
98
99template<class GeoField>
101(
103) const
104{
105 for (const auto& fld : fields)
106 {
107 decomposeField(fld)().write();
108 }
109}
110
111
112// ************************************************************************* //
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.
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
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
const Time & time() const noexcept
Return time registry.
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose list of fields.
tmp< GeometricField< Type, pointPatchField, pointMesh > > decomposeField(const GeometricField< Type, pointPatchField, pointMesh > &) const
Decompose point field.
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:126
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:114
Abstract base class for point-mesh patch fields.
Foam::processorPointPatchField.
A class for managing temporary objects.
Definition: tmp.H:65
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