fvMeshToolsTemplates.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 | Copyright (C) 2012-2013 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMeshTools.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 template<class GeoField>
34 (
35  fvMesh& mesh,
36  const dictionary& patchFieldDict,
37  const word& defaultPatchFieldType,
38  const typename GeoField::value_type& defaultPatchValue
39 )
40 {
42  (
43  mesh.objectRegistry::lookupClass<GeoField>()
44  );
45 
46  forAllIter(typename HashTable<GeoField*>, flds, iter)
47  {
48  GeoField& fld = *iter();
49 
50  typename GeoField::GeometricBoundaryField& bfld =
51  fld.boundaryField();
52 
53  label sz = bfld.size();
54  bfld.setSize(sz+1);
55 
56  if (patchFieldDict.found(fld.name()))
57  {
58  bfld.set
59  (
60  sz,
62  (
63  mesh.boundary()[sz],
64  fld.dimensionedInternalField(),
65  patchFieldDict.subDict(fld.name())
66  )
67  );
68  }
69  else
70  {
71  bfld.set
72  (
73  sz,
75  (
76  defaultPatchFieldType,
77  mesh.boundary()[sz],
78  fld.dimensionedInternalField()
79  )
80  );
81  bfld[sz] == defaultPatchValue;
82  }
83  }
84 }
85 
86 
87 template<class GeoField>
89 (
90  fvMesh& mesh,
91  const label patchI,
92  const dictionary& patchFieldDict
93 )
94 {
96  (
97  mesh.objectRegistry::lookupClass<GeoField>()
98  );
99 
100  forAllIter(typename HashTable<GeoField*>, flds, iter)
101  {
102  GeoField& fld = *iter();
103 
104  typename GeoField::GeometricBoundaryField& bfld =
105  fld.boundaryField();
106 
107  if (patchFieldDict.found(fld.name()))
108  {
109  bfld.set
110  (
111  patchI,
113  (
114  mesh.boundary()[patchI],
115  fld.dimensionedInternalField(),
116  patchFieldDict.subDict(fld.name())
117  )
118  );
119  }
120  }
121 }
122 
123 
124 
125 
126 template<class GeoField>
128 (
129  fvMesh& mesh,
130  const label patchI,
131  const typename GeoField::value_type& value
132 )
133 {
135  (
136  mesh.objectRegistry::lookupClass<GeoField>()
137  );
138 
139  forAllIter(typename HashTable<GeoField*>, flds, iter)
140  {
141  GeoField& fld = *iter();
142 
143  typename GeoField::GeometricBoundaryField& bfld =
144  fld.boundaryField();
145 
146  bfld[patchI] == value;
147  }
148 }
149 
150 
151 // Remove last patch field
152 template<class GeoField>
154 {
156  (
157  mesh.objectRegistry::lookupClass<GeoField>()
158  );
159 
160  forAllIter(typename HashTable<GeoField*>, flds, iter)
161  {
162  GeoField& fld = *iter();
163  fld.boundaryField().setSize(nPatches);
164  }
165 }
166 
167 
168 // Reorder patch field
169 template<class GeoField>
171 (
172  fvMesh& mesh,
173  const labelList& oldToNew
174 )
175 {
177  (
178  mesh.objectRegistry::lookupClass<GeoField>()
179  );
180 
181  forAllIter(typename HashTable<GeoField*>, flds, iter)
182  {
183  GeoField& fld = *iter();
184 
185  typename GeoField::GeometricBoundaryField& bfld =
186  fld.boundaryField();
187 
188  bfld.reorder(oldToNew);
189  }
190 }
191 
192 
193 // ************************************************************************* //
volFields.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::fvMeshTools::addPatchFields
static void addPatchFields(fvMesh &, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const typename GeoField::value_type &defaultPatchValue)
Definition: fvMeshToolsTemplates.C:34
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
nPatches
label nPatches
Definition: readKivaGrid.H:402
surfaceFields.H
Foam::surfaceFields.
Foam::fvMeshTools::reorderPatchFields
static void reorderPatchFields(fvMesh &, const labelList &oldToNew)
Definition: fvMeshToolsTemplates.C:171
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::fvMeshTools::trimPatchFields
static void trimPatchFields(fvMesh &, const label nPatches)
Definition: fvMeshToolsTemplates.C:153
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
fvMeshTools.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::fvMeshTools::setPatchFields
static void setPatchFields(fvMesh &mesh, const label patchI, const dictionary &patchFieldDict)
Set patchFields according to dictionary.
Definition: fvMeshToolsTemplates.C:89
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631