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 -------------------------------------------------------------------------------
10 License
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 
26 Application
27  setFields
28 
29 Group
30  grpPreProcessingUtilities
31 
32 Description
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 
45 using namespace Foam;
46 
47 template<class Type>
48 bool 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,
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,
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 
136 class setCellField
137 {
138 
139 public:
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 
199 template<class Type>
200 bool 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,
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,
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 
291  Pstream::listCombineScatter(nChanged);
292 
294  Boundary& fieldBf = field.boundaryFieldRef();
295 
296  // Reassign.
297  forAll(field.boundaryField(), patchi)
298  {
299  if (nChanged[patchi] > 0)
300  {
301  Info<< " On patch "
302  << field.boundaryField()[patchi].patch().name()
303  << " set " << nChanged[patchi] << " values" << endl;
304  fieldBf[patchi] == SubField<Type>
305  (
306  allBoundaryValues,
307  fieldBf[patchi].size(),
308  fieldBf[patchi].patch().start()
309  - mesh.nInternalFaces()
310  );
311  }
312  }
313 
314  if (!field.write())
315  {
317  << "Failed writing field " << field.name() << exit(FatalError);
318  }
319  }
320  else
321  {
323  << "Field " << fieldName << " not found" << endl;
324 
325  // Consume value
326  (void)pTraits<Type>(fieldValueStream);
327  }
328 
329  return true;
330 }
331 
332 
333 class setFaceField
334 {
335 
336 public:
337 
338  setFaceField()
339  {}
340 
341  autoPtr<setFaceField> clone() const
342  {
344  }
345 
346  class iNew
347  {
348  const fvMesh& mesh_;
349  const labelList& selectedFaces_;
350 
351  public:
352 
353  iNew(const fvMesh& mesh, const labelList& selectedFaces)
354  :
355  mesh_(mesh),
356  selectedFaces_(selectedFaces)
357  {}
358 
359  iNew(const fvMesh& mesh, labelList&& selectedFaces)
360  :
361  mesh_(mesh),
362  selectedFaces_(std::move(selectedFaces))
363  {}
364 
365  autoPtr<setFaceField> operator()(Istream& fieldValues) const
366  {
367  word fieldType(fieldValues);
368 
369  if
370  (
371  !(
372  setFaceFieldType<scalar>
373  (fieldType, mesh_, selectedFaces_, fieldValues)
374  || setFaceFieldType<vector>
375  (fieldType, mesh_, selectedFaces_, fieldValues)
376  || setFaceFieldType<sphericalTensor>
377  (fieldType, mesh_, selectedFaces_, fieldValues)
378  || setFaceFieldType<symmTensor>
379  (fieldType, mesh_, selectedFaces_, fieldValues)
380  || setFaceFieldType<tensor>
381  (fieldType, mesh_, selectedFaces_, fieldValues)
382  )
383  )
384  {
386  << "field type " << fieldType << " not currently supported"
387  << endl;
388  }
389 
391  }
392  };
393 };
394 
395 
396 
397 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398 
399 int main(int argc, char *argv[])
400 {
402  (
403  "Set values on a selected set of cells/patch-faces via a dictionary"
404  );
405 
406  argList::addOption("dict", "file", "Alternative setFieldsDict");
407 
408  #include "addRegionOption.H"
409  #include "setRootCase.H"
410  #include "createTime.H"
411  #include "createNamedMesh.H"
412 
413  const word dictName("setFieldsDict");
414  #include "setSystemMeshDictionaryIO.H"
415 
416  Info<< "Reading " << dictIO.name() << nl << endl;
417 
418  IOdictionary setFieldsDict(dictIO);
419 
420  if (setFieldsDict.found("defaultFieldValues"))
421  {
422  Info<< "Setting field default values" << endl;
423  PtrList<setCellField> defaultFieldValues
424  (
425  setFieldsDict.lookup("defaultFieldValues"),
426  setCellField::iNew(mesh, labelList(mesh.nCells()))
427  );
428  Info<< endl;
429  }
430 
431 
432  Info<< "Setting field region values" << endl;
433 
434  PtrList<entry> regions(setFieldsDict.lookup("regions"));
435 
436  forAll(regions, regionI)
437  {
438  const entry& region = regions[regionI];
439 
440  autoPtr<topoSetSource> source =
441  topoSetSource::New(region.keyword(), mesh, region.dict());
442 
443  if (source().setType() == topoSetSource::CELLSET_SOURCE)
444  {
445  cellSet selectedCellSet
446  (
447  mesh,
448  "cellSet",
449  mesh.nCells()/10+1 // Reasonable size estimate.
450  );
451 
452  source->applyToSet
453  (
455  selectedCellSet
456  );
457 
458  PtrList<setCellField> fieldValues
459  (
460  region.dict().lookup("fieldValues"),
461  setCellField::iNew(mesh, selectedCellSet.sortedToc())
462  );
463  }
464  else if (source().setType() == topoSetSource::FACESET_SOURCE)
465  {
466  faceSet selectedFaceSet
467  (
468  mesh,
469  "faceSet",
470  mesh.nBoundaryFaces()/10+1
471  );
472 
473  source->applyToSet
474  (
476  selectedFaceSet
477  );
478 
479  PtrList<setFaceField> fieldValues
480  (
481  region.dict().lookup("fieldValues"),
482  setFaceField::iNew(mesh, selectedFaceSet.sortedToc())
483  );
484  }
485  }
486 
487 
488  Info<< "\nEnd\n" << endl;
489 
490  return 0;
491 }
492 
493 
494 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:63
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Foam::maxOp
Definition: ops.H:217
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
dictName
const word dictName("faMeshDefinition")
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::faceSet
A list of face labels.
Definition: faceSet.H:47
Foam::entry::keyword
const keyType & keyword() const noexcept
Definition: entry.H:191
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::topoSetSource::NEW
@ NEW
Create a new set and ADD elements to it.
Definition: topoSetSource.H:100
setSystemMeshDictionaryIO.H
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::topoSetSource::CELLSET_SOURCE
@ CELLSET_SOURCE
Cells as set.
Definition: topoSetSource.H:78
Foam::SubField< Type >
Foam::polyBoundaryMesh::patchID
const labelList & patchID() const
Definition: polyBoundaryMesh.C:449
Foam::primitiveMesh::nCells
label nCells() const noexcept
Definition: primitiveMeshI.H:89
Foam::topoSetSource::FACESET_SOURCE
@ FACESET_SOURCE
Faces as set.
Definition: topoSetSource.H:79
Foam::Field< Type >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::Info
messageStream Info
argList.H
faceSet.H
addRegionOption.H
Foam::primitiveMesh::nBoundaryFaces
label nBoundaryFaces() const noexcept
Definition: primitiveMeshI.H:77
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:379
Foam::FatalError
error FatalError
createNamedMesh.H
Required Variables.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
dictIO
IOobject dictIO
Definition: setConstantMeshDictionaryIO.H:1
Foam::topoSetSource::New
static autoPtr< topoSetSource > New(const word &topoSetSourceType, const polyMesh &mesh, const dictionary &dict)
Definition: topoSetSource.C:103
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
fvMesh.H
Foam
Definition: atmBoundaryLayer.C:26
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:47
Foam::entry::dict
virtual const dictionary & dict() const =0
Foam::HashTable::sortedToc
List< Key > sortedToc() const
Definition: HashTable.C:129
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:284
Foam::IOobject::name
const word & name() const noexcept
Definition: IOobjectI.H:58
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Definition: primitiveMeshI.H:96
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::foamVersion::patch
const std::string patch
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Definition: primitiveMeshI.H:71
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::pTraits< Type >
topoSetSource.H
createTime.H
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::plusEqOp
Definition: ops.H:66
cellSet.H
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Definition: combineGatherScatter.C:426
Foam::TimePaths::constant
const word & constant() const
Definition: TimePathsI.H:89
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181