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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2022 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 "fvMeshTools.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 template<class GeoField>
36 void Foam::fvMeshTools::addPatchFields
37 (
38  fvMesh& mesh,
39  const dictionary& patchFieldDict,
40  const word& defaultPatchFieldType,
41  const typename GeoField::value_type& defaultPatchValue
42 )
43 {
44  HashTable<GeoField*> flds
45  (
46  mesh.objectRegistry::lookupClass<GeoField>()
47  );
48 
49  forAllIters(flds, iter)
50  {
51  GeoField& fld = *iter.val();
52  auto& bfld = fld.boundaryFieldRef();
53 
54  const label newPatchi = bfld.size();
55  bfld.resize(newPatchi+1);
56 
57  const dictionary* dict = patchFieldDict.findDict(fld.name());
58 
59  if (dict)
60  {
61  bfld.set
62  (
63  newPatchi,
65  (
66  mesh.boundary()[newPatchi],
67  fld.internalField(),
68  *dict
69  )
70  );
71  }
72  else
73  {
74  bfld.set
75  (
76  newPatchi,
78  (
79  defaultPatchFieldType,
80  mesh.boundary()[newPatchi],
81  fld.internalField()
82  )
83  );
84  bfld[newPatchi] == defaultPatchValue;
85  }
86  }
87 }
88 
89 
90 template<class GeoField>
91 void Foam::fvMeshTools::setPatchFields
92 (
93  fvMesh& mesh,
94  const label patchi,
95  const dictionary& patchFieldDict
96 )
97 {
98  HashTable<GeoField*> flds
99  (
100  mesh.objectRegistry::lookupClass<GeoField>()
101  );
102 
103  forAllIters(flds, iter)
104  {
105  GeoField& fld = *iter.val();
106  auto& bfld = fld.boundaryFieldRef();
107 
108  const dictionary* dict = patchFieldDict.findDict(fld.name());
109 
110  if (dict)
111  {
112  bfld.set
113  (
114  patchi,
116  (
117  mesh.boundary()[patchi],
118  fld.internalField(),
119  *dict
120  )
121  );
122  }
123  }
124 }
125 
126 
127 template<class GeoField>
128 void Foam::fvMeshTools::setPatchFields
129 (
130  fvMesh& mesh,
131  const label patchi,
132  const typename GeoField::value_type& value
133 )
134 {
135  HashTable<GeoField*> flds
136  (
137  mesh.objectRegistry::lookupClass<GeoField>()
138  );
139 
140  forAllIters(flds, iter)
141  {
142  GeoField& fld = *iter.val();
143  auto& bfld = fld.boundaryFieldRef();
144 
145  bfld[patchi] == value;
146  }
147 }
148 
149 
150 // Remove last patch field
151 template<class GeoField>
152 void Foam::fvMeshTools::trimPatchFields(fvMesh& mesh, const label nPatches)
153 {
154  HashTable<GeoField*> flds
155  (
156  mesh.objectRegistry::lookupClass<GeoField>()
157  );
158 
159  forAllIters(flds, iter)
160  {
161  GeoField& fld = *iter.val();
162  fld.boundaryFieldRef().resize(nPatches);
163  }
164 }
165 
166 
167 // Reorder patch field
168 template<class GeoField>
169 void Foam::fvMeshTools::reorderPatchFields
170 (
171  fvMesh& mesh,
172  const labelList& oldToNew
173 )
174 {
175  HashTable<GeoField*> flds
176  (
177  mesh.objectRegistry::lookupClass<GeoField>()
178  );
179 
180  forAllIters(flds, iter)
181  {
182  GeoField& fld = *iter.val();
183  auto& bfld = fld.boundaryFieldRef();
184 
185  bfld.reorder(oldToNew);
186  }
187 }
188 
189 
190 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
Foam::surfaceFields.
dictionary dict
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dynamicFvMesh & mesh
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:329
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:397
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:777
List< label > labelList
A List of labels.
Definition: List.H:62
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:120