faFieldDecomposerFields.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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "processorFaPatchField.H"
31 #include "processorFaePatchField.H"
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
38 (
40 ) const
41 {
42  // Create and map the internal field values
43  Field<Type> internalField(field.internalField(), faceAddressing_);
44 
45  // Create and map the patch field values
46  PtrList<faPatchField<Type>> patchFields(boundaryAddressing_.size());
47 
48  forAll(boundaryAddressing_, patchi)
49  {
50  const label oldPatchi = boundaryAddressing_[patchi];
51 
52  if (oldPatchi >= 0)
53  {
54  patchFields.set
55  (
56  patchi,
58  (
59  field.boundaryField()[oldPatchi],
60  procMesh_.boundary()[patchi],
62  patchFieldDecomposerPtrs_[patchi]
63  )
64  );
65  }
66  else
67  {
68  patchFields.set
69  (
70  patchi,
72  (
73  procMesh_.boundary()[patchi],
76  (
77  field.internalField(),
78  processorAreaPatchFieldDecomposerPtrs_[patchi]
79  )
80  )
81  );
82  }
83  }
84 
85  // Create the field for the processor
86  return
88  (
89  IOobject
90  (
91  field.name(),
92  procMesh_.time().timeName(),
93  procMesh_(),
94  IOobject::NO_READ,
95  IOobject::NO_WRITE
96  ),
97  procMesh_,
98  field.dimensions(),
99  internalField,
100  patchFields
101  );
102 }
103 
104 
105 template<class Type>
108 (
110 ) const
111 {
112  labelList mapAddr
113  (
115  (
116  edgeAddressing_,
117  procMesh_.nInternalEdges()
118  )
119  );
120  forAll(mapAddr, i)
121  {
122  mapAddr[i] -= 1;
123  }
124 
125  // Create and map the internal field values
126  Field<Type> internalField
127  (
128  field.internalField(),
129  mapAddr
130  );
131 
132  // Problem with addressing when a processor patch picks up both internal
133  // edges and edges from cyclic boundaries. This is a bit of a hack, but
134  // I cannot find a better solution without making the internal storage
135  // mechanism for edgeFields correspond to the one of edges in polyMesh
136  // (i.e. using slices)
137  Field<Type> allEdgeField(field.mesh().nEdges());
138 
139  forAll(field.internalField(), i)
140  {
141  allEdgeField[i] = field.internalField()[i];
142  }
143 
144  forAll(field.boundaryField(), patchi)
145  {
146  const Field<Type>& p = field.boundaryField()[patchi];
147 
148  const label patchStart = field.mesh().boundary()[patchi].start();
149 
150  forAll(p, i)
151  {
152  allEdgeField[patchStart + i] = p[i];
153  }
154  }
155 
156  // Create and map the patch field values
157  PtrList<faePatchField<Type>> patchFields(boundaryAddressing_.size());
158 
159  forAll(boundaryAddressing_, patchi)
160  {
161  const label oldPatchi = boundaryAddressing_[patchi];
162 
163  if (oldPatchi >= 0)
164  {
165  patchFields.set
166  (
167  patchi,
169  (
170  field.boundaryField()[oldPatchi],
171  procMesh_.boundary()[patchi],
173  patchFieldDecomposerPtrs_[patchi]
174  )
175  );
176  }
177  else
178  {
179  patchFields.set
180  (
181  patchi,
183  (
184  procMesh_.boundary()[patchi],
187  (
188  allEdgeField,
189  processorEdgePatchFieldDecomposerPtrs_[patchi]
190  )
191  )
192  );
193  }
194  }
195 
196  // Create the field for the processor
197  return
199  (
200  IOobject
201  (
202  field.name(),
203  procMesh_.time().timeName(),
204  procMesh_(),
205  IOobject::NO_READ,
206  IOobject::NO_WRITE
207  ),
208  procMesh_,
209  field.dimensions(),
210  internalField,
211  patchFields
212  );
213 }
214 
215 
216 template<class GeoField>
218 (
220 ) const
221 {
222  forAll(fields, fieldi)
223  {
224  decomposeField(fields[fieldi])().write();
225  }
226 }
227 
228 
229 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::faPatchField
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: areaFieldsFwd.H:50
Foam::faePatchField
faePatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cove...
Definition: edgeFieldsFwd.H:49
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
processorFaPatchField.H
Foam::processorFaePatchField
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
Definition: processorFaePatchField.H:55
processorFaePatchField.H
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::faFieldDecomposer::decomposeField
tmp< GeometricField< Type, faPatchField, areaMesh > > decomposeField(const GeometricField< Type, faPatchField, areaMesh > &field) const
Decompose area field.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::processorFaPatchField
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
Definition: processorFaPatchField.H:58
Foam::Field
Generic templated field type.
Definition: Field.H:63
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::List< label >
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::faFieldDecomposer::decomposeFields
void decomposeFields(const PtrList< GeoField > &fields) const
Definition: faFieldDecomposerFields.C:218
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
faFieldDecomposer.H