setFields.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Application
27 setFields
28
29Group
30 grpPreProcessingUtilities
31
32Description
33 Set values on a selected set of cells/patch-faces via a dictionary.
34
35\*---------------------------------------------------------------------------*/
36
37#include "argList.H"
38#include "Time.H"
39#include "fvMesh.H"
40#include "topoSetSource.H"
41#include "cellSet.H"
42#include "faceSet.H"
43#include "volFields.H"
44
45using namespace Foam;
46
47template<class Type>
48bool setCellFieldType
49(
50 const word& fieldTypeDesc,
51 const fvMesh& mesh,
52 const labelList& selectedCells,
53 Istream& fieldValueStream
54)
55{
57
58 if (fieldTypeDesc != fieldType::typeName + "Value")
59 {
60 return false;
61 }
62
63 word fieldName(fieldValueStream);
64
65 // Check the current time directory
66 IOobject fieldHeader
67 (
68 fieldName,
69 mesh.time().timeName(),
70 mesh,
71 IOobject::MUST_READ
72 );
73
74 // Check the "constant" directory
75 if (!fieldHeader.typeHeaderOk<fieldType>(true))
76 {
77 fieldHeader = IOobject
78 (
79 fieldName,
80 mesh.time().constant(),
81 mesh,
82 IOobject::MUST_READ
83 );
84 }
85
86 // Check field exists
87 if (fieldHeader.typeHeaderOk<fieldType>(true))
88 {
89 Info<< " Setting internal values of "
90 << fieldHeader.headerClassName()
91 << " " << fieldName << endl;
92
93 fieldType field(fieldHeader, mesh, false);
94
95 const Type value = pTraits<Type>(fieldValueStream);
96
97 if (selectedCells.size() == field.size())
98 {
99 field.primitiveFieldRef() = value;
100 }
101 else
102 {
103 forAll(selectedCells, celli)
104 {
105 field[selectedCells[celli]] = value;
106 }
107 }
108
110 Boundary& fieldBf = field.boundaryFieldRef();
111
112 forAll(field.boundaryField(), patchi)
113 {
114 fieldBf[patchi] = fieldBf[patchi].patchInternalField();
115 }
116
117 if (!field.write())
118 {
120 << "Failed writing field " << fieldName << endl;
121 }
122 }
123 else
124 {
126 << "Field " << fieldName << " not found" << endl;
127
128 // Consume value
129 (void)pTraits<Type>(fieldValueStream);
130 }
131
132 return true;
133}
134
135
136class setCellField
137{
138
139public:
140
141 setCellField()
142 {}
143
144 autoPtr<setCellField> clone() const
145 {
147 }
148
149 class iNew
150 {
151 const fvMesh& mesh_;
152 const labelList& selectedCells_;
153
154 public:
155
156 iNew(const fvMesh& mesh, const labelList& selectedCells)
157 :
158 mesh_(mesh),
159 selectedCells_(selectedCells)
160 {}
161
162 iNew(const fvMesh& mesh, labelList&& selectedCells)
163 :
164 mesh_(mesh),
165 selectedCells_(std::move(selectedCells))
166 {}
167
168 autoPtr<setCellField> operator()(Istream& fieldValues) const
169 {
170 word fieldType(fieldValues);
171
172 if
173 (
174 !(
175 setCellFieldType<scalar>
176 (fieldType, mesh_, selectedCells_, fieldValues)
177 || setCellFieldType<vector>
178 (fieldType, mesh_, selectedCells_, fieldValues)
179 || setCellFieldType<sphericalTensor>
180 (fieldType, mesh_, selectedCells_, fieldValues)
181 || setCellFieldType<symmTensor>
182 (fieldType, mesh_, selectedCells_, fieldValues)
183 || setCellFieldType<tensor>
184 (fieldType, mesh_, selectedCells_, fieldValues)
185 )
186 )
187 {
189 << "field type " << fieldType << " not currently supported"
190 << endl;
191 }
192
194 }
195 };
196};
197
198
199template<class Type>
200bool setFaceFieldType
201(
202 const word& fieldTypeDesc,
203 const fvMesh& mesh,
204 const labelList& selectedFaces,
205 Istream& fieldValueStream
206)
207{
209
210 if (fieldTypeDesc != fieldType::typeName + "Value")
211 {
212 return false;
213 }
214
215 word fieldName(fieldValueStream);
216
217 // Check the current time directory
218 IOobject fieldHeader
219 (
220 fieldName,
221 mesh.time().timeName(),
222 mesh,
223 IOobject::MUST_READ
224 );
225
226 // Check the "constant" directory
227 if (!fieldHeader.typeHeaderOk<fieldType>(true))
228 {
229 fieldHeader = IOobject
230 (
231 fieldName,
232 mesh.time().constant(),
233 mesh,
234 IOobject::MUST_READ
235 );
236 }
237
238 // Check field exists
239 if (fieldHeader.typeHeaderOk<fieldType>(true))
240 {
241 Info<< " Setting patchField values of "
242 << fieldHeader.headerClassName()
243 << " " << fieldName << endl;
244
245 fieldType field(fieldHeader, mesh);
246
247 const Type value = pTraits<Type>(fieldValueStream);
248
249 // Create flat list of selected faces and their value.
250 Field<Type> allBoundaryValues(mesh.nBoundaryFaces());
251 forAll(field.boundaryField(), patchi)
252 {
254 (
255 allBoundaryValues,
256 field.boundaryField()[patchi].size(),
257 field.boundaryField()[patchi].patch().start()
258 - mesh.nInternalFaces()
259 ) = field.boundaryField()[patchi];
260 }
261
262 // Override
263 bool hasWarned = false;
264 labelList nChanged
265 (
266 returnReduce(field.boundaryField().size(), maxOp<label>()),
267 0
268 );
269 forAll(selectedFaces, i)
270 {
271 label facei = selectedFaces[i];
272 if (mesh.isInternalFace(facei))
273 {
274 if (!hasWarned)
275 {
276 hasWarned = true;
278 << "Ignoring internal face " << facei
279 << ". Suppressing further warnings." << endl;
280 }
281 }
282 else
283 {
284 label bFacei = facei-mesh.nInternalFaces();
285 allBoundaryValues[bFacei] = value;
286 nChanged[mesh.boundaryMesh().patchID()[bFacei]]++;
287 }
288 }
289
290 Pstream::listCombineAllGather(nChanged, plusEqOp<label>());
291
292 auto& fieldBf = field.boundaryFieldRef();
293
294 // Reassign.
295 forAll(field.boundaryField(), patchi)
296 {
297 if (nChanged[patchi] > 0)
298 {
299 Info<< " On patch "
300 << field.boundaryField()[patchi].patch().name()
301 << " set " << nChanged[patchi] << " values" << endl;
302 fieldBf[patchi] == SubField<Type>
303 (
304 allBoundaryValues,
305 fieldBf[patchi].size(),
306 fieldBf[patchi].patch().start()
307 - mesh.nInternalFaces()
308 );
309 }
310 }
311
312 if (!field.write())
313 {
315 << "Failed writing field " << field.name() << exit(FatalError);
316 }
317 }
318 else
319 {
321 << "Field " << fieldName << " not found" << endl;
322
323 // Consume value
324 (void)pTraits<Type>(fieldValueStream);
325 }
326
327 return true;
328}
329
330
331class setFaceField
332{
333
334public:
335
336 setFaceField()
337 {}
338
339 autoPtr<setFaceField> clone() const
340 {
342 }
343
344 class iNew
345 {
346 const fvMesh& mesh_;
347 const labelList& selectedFaces_;
348
349 public:
350
351 iNew(const fvMesh& mesh, const labelList& selectedFaces)
352 :
353 mesh_(mesh),
354 selectedFaces_(selectedFaces)
355 {}
356
357 iNew(const fvMesh& mesh, labelList&& selectedFaces)
358 :
359 mesh_(mesh),
360 selectedFaces_(std::move(selectedFaces))
361 {}
362
363 autoPtr<setFaceField> operator()(Istream& fieldValues) const
364 {
365 word fieldType(fieldValues);
366
367 if
368 (
369 !(
370 setFaceFieldType<scalar>
371 (fieldType, mesh_, selectedFaces_, fieldValues)
372 || setFaceFieldType<vector>
373 (fieldType, mesh_, selectedFaces_, fieldValues)
374 || setFaceFieldType<sphericalTensor>
375 (fieldType, mesh_, selectedFaces_, fieldValues)
376 || setFaceFieldType<symmTensor>
377 (fieldType, mesh_, selectedFaces_, fieldValues)
378 || setFaceFieldType<tensor>
379 (fieldType, mesh_, selectedFaces_, fieldValues)
380 )
381 )
382 {
384 << "field type " << fieldType << " not currently supported"
385 << endl;
386 }
387
389 }
390 };
391};
392
393
394
395// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
396
397int main(int argc, char *argv[])
398{
399 argList::addNote
400 (
401 "Set values on a selected set of cells/patch-faces via a dictionary"
402 );
403
404 argList::addOption("dict", "file", "Alternative setFieldsDict");
405
406 #include "addRegionOption.H"
407 #include "setRootCase.H"
408 #include "createTime.H"
409 #include "createNamedMesh.H"
410
411 const word dictName("setFieldsDict");
413
414 Info<< "Reading " << dictIO.name() << nl << endl;
415
416 IOdictionary setFieldsDict(dictIO);
417
418 if (setFieldsDict.found("defaultFieldValues"))
419 {
420 Info<< "Setting field default values" << endl;
421 PtrList<setCellField> defaultFieldValues
422 (
423 setFieldsDict.lookup("defaultFieldValues"),
424 setCellField::iNew(mesh, labelList(mesh.nCells()))
425 );
426 Info<< endl;
427 }
428
429
430 Info<< "Setting field region values" << endl;
431
432 PtrList<entry> regions(setFieldsDict.lookup("regions"));
433
434 forAll(regions, regionI)
435 {
436 const entry& region = regions[regionI];
437
439 topoSetSource::New(region.keyword(), mesh, region.dict());
440
441 if (source().setType() == topoSetSource::CELLSET_SOURCE)
442 {
443 cellSet selectedCellSet
444 (
445 mesh,
446 "cellSet",
447 mesh.nCells()/10+1 // Reasonable size estimate.
448 );
449
450 source->applyToSet
451 (
452 topoSetSource::NEW,
453 selectedCellSet
454 );
455
456 PtrList<setCellField> fieldValues
457 (
458 region.dict().lookup("fieldValues"),
459 setCellField::iNew(mesh, selectedCellSet.sortedToc())
460 );
461 }
462 else if (source().setType() == topoSetSource::FACESET_SOURCE)
463 {
464 faceSet selectedFaceSet
465 (
466 mesh,
467 "faceSet",
468 mesh.nBoundaryFaces()/10+1
469 );
470
471 source->applyToSet
472 (
473 topoSetSource::NEW,
474 selectedFaceSet
475 );
476
477 PtrList<setFaceField> fieldValues
478 (
479 region.dict().lookup("fieldValues"),
480 setFaceField::iNew(mesh, selectedFaceSet.sortedToc())
481 );
482 }
483 }
484
485
486 Info<< "\nEnd\n" << endl;
487
488 return 0;
489}
490
491
492// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: SubField.H:62
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A collection of cell labels.
Definition: cellSet.H:54
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
const keyType & keyword() const noexcept
Return keyword.
Definition: entry.H:195
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
A list of face labels.
Definition: faceSet.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
dynamicFvMesh & mesh
Required Variables.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const word dictName("faMeshDefinition")
#define WarningInFunction
Report a warning using Foam::Warning.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
IOobject dictIO
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333