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 | Copyright (C) 2011-2015 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 Description
25  Set values on a selected set of cells/patchfaces through a dictionary.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "argList.H"
30 #include "Time.H"
31 #include "fvMesh.H"
32 #include "topoSetSource.H"
33 #include "cellSet.H"
34 #include "faceSet.H"
35 #include "volFields.H"
36 
37 using namespace Foam;
38 
39 template<class Type>
40 bool setCellFieldType
41 (
42  const word& fieldTypeDesc,
43  const fvMesh& mesh,
44  const labelList& selectedCells,
45  Istream& fieldValueStream
46 )
47 {
49 
50  if (fieldTypeDesc != fieldType::typeName + "Value")
51  {
52  return false;
53  }
54 
55  word fieldName(fieldValueStream);
56 
57  // Check the current time directory
58  IOobject fieldHeader
59  (
60  fieldName,
61  mesh.time().timeName(),
62  mesh,
64  );
65 
66  // Check the "constant" directory
67  if (!fieldHeader.headerOk())
68  {
69  fieldHeader = IOobject
70  (
71  fieldName,
72  mesh.time().constant(),
73  mesh,
75  );
76  }
77 
78  // Check field exists
79  if (fieldHeader.headerOk())
80  {
81  Info<< " Setting internal values of "
82  << fieldHeader.headerClassName()
83  << " " << fieldName << endl;
84 
85  fieldType field(fieldHeader, mesh);
86 
87  const Type& value = pTraits<Type>(fieldValueStream);
88 
89  if (selectedCells.size() == field.size())
90  {
91  field.internalField() = value;
92  }
93  else
94  {
95  forAll(selectedCells, celli)
96  {
97  field[selectedCells[celli]] = value;
98  }
99  }
100 
101  forAll(field.boundaryField(), patchi)
102  {
103  field.boundaryField()[patchi] =
104  field.boundaryField()[patchi].patchInternalField();
105  }
106 
107  if (!field.write())
108  {
110  << "Failed writing field " << fieldName << endl;
111  }
112  }
113  else
114  {
116  << "Field " << fieldName << " not found" << endl;
117 
118  // Consume value
119  (void)pTraits<Type>(fieldValueStream);
120  }
121 
122  return true;
123 }
124 
125 
126 class setCellField
127 {
128 
129 public:
130 
131  setCellField()
132  {}
133 
134  autoPtr<setCellField> clone() const
135  {
136  return autoPtr<setCellField>(new setCellField());
137  }
138 
139  class iNew
140  {
141  const fvMesh& mesh_;
142  const labelList& selectedCells_;
143 
144  public:
145 
146  iNew(const fvMesh& mesh, const labelList& selectedCells)
147  :
148  mesh_(mesh),
149  selectedCells_(selectedCells)
150  {}
151 
152  autoPtr<setCellField> operator()(Istream& fieldValues) const
153  {
154  word fieldType(fieldValues);
155 
156  if
157  (
158  !(
159  setCellFieldType<scalar>
160  (fieldType, mesh_, selectedCells_, fieldValues)
161  || setCellFieldType<vector>
162  (fieldType, mesh_, selectedCells_, fieldValues)
163  || setCellFieldType<sphericalTensor>
164  (fieldType, mesh_, selectedCells_, fieldValues)
165  || setCellFieldType<symmTensor>
166  (fieldType, mesh_, selectedCells_, fieldValues)
167  || setCellFieldType<tensor>
168  (fieldType, mesh_, selectedCells_, fieldValues)
169  )
170  )
171  {
173  << "field type " << fieldType << " not currently supported"
174  << endl;
175  }
176 
177  return autoPtr<setCellField>(new setCellField());
178  }
179  };
180 };
181 
182 
183 template<class Type>
184 bool setFaceFieldType
185 (
186  const word& fieldTypeDesc,
187  const fvMesh& mesh,
188  const labelList& selectedFaces,
189  Istream& fieldValueStream
190 )
191 {
193 
194  if (fieldTypeDesc != fieldType::typeName + "Value")
195  {
196  return false;
197  }
198 
199  word fieldName(fieldValueStream);
200 
201  // Check the current time directory
202  IOobject fieldHeader
203  (
204  fieldName,
205  mesh.time().timeName(),
206  mesh,
208  );
209 
210  // Check the "constant" directory
211  if (!fieldHeader.headerOk())
212  {
213  fieldHeader = IOobject
214  (
215  fieldName,
216  mesh.time().constant(),
217  mesh,
219  );
220  }
221 
222  // Check field exists
223  if (fieldHeader.headerOk())
224  {
225  Info<< " Setting patchField values of "
226  << fieldHeader.headerClassName()
227  << " " << fieldName << endl;
228 
229  fieldType field(fieldHeader, mesh);
230 
231  const Type& value = pTraits<Type>(fieldValueStream);
232 
233  // Create flat list of selected faces and their value.
234  Field<Type> allBoundaryValues(mesh.nFaces()-mesh.nInternalFaces());
235  forAll(field.boundaryField(), patchi)
236  {
238  (
239  allBoundaryValues,
240  field.boundaryField()[patchi].size(),
241  field.boundaryField()[patchi].patch().start()
242  - mesh.nInternalFaces()
243  ).assign(field.boundaryField()[patchi]);
244  }
245 
246  // Override
247  bool hasWarned = false;
248  labelList nChanged
249  (
250  returnReduce(field.boundaryField().size(), maxOp<label>()),
251  0
252  );
253  forAll(selectedFaces, i)
254  {
255  label facei = selectedFaces[i];
256  if (mesh.isInternalFace(facei))
257  {
258  if (!hasWarned)
259  {
260  hasWarned = true;
262  << "Ignoring internal face " << facei
263  << ". Suppressing further warnings." << endl;
264  }
265  }
266  else
267  {
268  label bFaceI = facei-mesh.nInternalFaces();
269  allBoundaryValues[bFaceI] = value;
270  nChanged[mesh.boundaryMesh().patchID()[bFaceI]]++;
271  }
272  }
273 
275  Pstream::listCombineScatter(nChanged);
276 
277  // Reassign.
278  forAll(field.boundaryField(), patchi)
279  {
280  if (nChanged[patchi] > 0)
281  {
282  Info<< " On patch "
283  << field.boundaryField()[patchi].patch().name()
284  << " set " << nChanged[patchi] << " values" << endl;
285  field.boundaryField()[patchi] == SubField<Type>
286  (
287  allBoundaryValues,
288  field.boundaryField()[patchi].size(),
289  field.boundaryField()[patchi].patch().start()
290  - mesh.nInternalFaces()
291  );
292  }
293  }
294 
295  if (!field.write())
296  {
298  << "Failed writing field " << field.name() << exit(FatalError);
299  }
300  }
301  else
302  {
304  << "Field " << fieldName << " not found" << endl;
305 
306  // Consume value
307  (void)pTraits<Type>(fieldValueStream);
308  }
309 
310  return true;
311 }
312 
313 
314 class setFaceField
315 {
316 
317 public:
318 
319  setFaceField()
320  {}
321 
322  autoPtr<setFaceField> clone() const
323  {
324  return autoPtr<setFaceField>(new setFaceField());
325  }
326 
327  class iNew
328  {
329  const fvMesh& mesh_;
330  const labelList& selectedFaces_;
331 
332  public:
333 
334  iNew(const fvMesh& mesh, const labelList& selectedFaces)
335  :
336  mesh_(mesh),
337  selectedFaces_(selectedFaces)
338  {}
339 
340  autoPtr<setFaceField> operator()(Istream& fieldValues) const
341  {
342  word fieldType(fieldValues);
343 
344  if
345  (
346  !(
347  setFaceFieldType<scalar>
348  (fieldType, mesh_, selectedFaces_, fieldValues)
349  || setFaceFieldType<vector>
350  (fieldType, mesh_, selectedFaces_, fieldValues)
351  || setFaceFieldType<sphericalTensor>
352  (fieldType, mesh_, selectedFaces_, fieldValues)
353  || setFaceFieldType<symmTensor>
354  (fieldType, mesh_, selectedFaces_, fieldValues)
355  || setFaceFieldType<tensor>
356  (fieldType, mesh_, selectedFaces_, fieldValues)
357  )
358  )
359  {
361  << "field type " << fieldType << " not currently supported"
362  << endl;
363  }
364 
365  return autoPtr<setFaceField>(new setFaceField());
366  }
367  };
368 };
369 
370 
371 
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373 
374 int main(int argc, char *argv[])
375 {
376  #include "addRegionOption.H"
377  #include "setRootCase.H"
378  #include "createTime.H"
379  #include "createNamedMesh.H"
380 
381  Info<< "Reading setFieldsDict\n" << endl;
382 
383  IOdictionary setFieldsDict
384  (
385  IOobject
386  (
387  "setFieldsDict",
388  runTime.system(),
389  mesh,
392  )
393  );
394 
395  if (setFieldsDict.found("defaultFieldValues"))
396  {
397  Info<< "Setting field default values" << endl;
398  PtrList<setCellField> defaultFieldValues
399  (
400  setFieldsDict.lookup("defaultFieldValues"),
401  setCellField::iNew(mesh, labelList(mesh.nCells()))
402  );
403  Info<< endl;
404  }
405 
406 
407  Info<< "Setting field region values" << endl;
408 
409  PtrList<entry> regions(setFieldsDict.lookup("regions"));
410 
411  forAll(regions, regionI)
412  {
413  const entry& region = regions[regionI];
414 
415  autoPtr<topoSetSource> source =
416  topoSetSource::New(region.keyword(), mesh, region.dict());
417 
418  if (source().setType() == topoSetSource::CELLSETSOURCE)
419  {
420  cellSet selectedCellSet
421  (
422  mesh,
423  "cellSet",
424  mesh.nCells()/10+1 // Reasonable size estimate.
425  );
426 
427  source->applyToSet
428  (
430  selectedCellSet
431  );
432 
433  PtrList<setCellField> fieldValues
434  (
435  region.dict().lookup("fieldValues"),
436  setCellField::iNew(mesh, selectedCellSet.toc())
437  );
438  }
439  else if (source().setType() == topoSetSource::FACESETSOURCE)
440  {
441  faceSet selectedFaceSet
442  (
443  mesh,
444  "faceSet",
445  (mesh.nFaces()-mesh.nInternalFaces())/10+1
446  );
447 
448  source->applyToSet
449  (
451  selectedFaceSet
452  );
453 
454  PtrList<setFaceField> fieldValues
455  (
456  region.dict().lookup("fieldValues"),
457  setFaceField::iNew(mesh, selectedFaceSet.toc())
458  );
459  }
460  }
461 
462 
463  Info<< "\nEnd\n" << endl;
464 
465  return 0;
466 }
467 
468 
469 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:65
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::topoSetSource::FACESETSOURCE
@ FACESETSOURCE
Definition: topoSetSource.H:73
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::entry::keyword
const keyType & keyword() const
Return keyword.
Definition: entry.H:120
Foam::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::topoSetSource::NEW
@ NEW
Definition: topoSetSource.H:85
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::SubField< Type >
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
Foam::polyBoundaryMesh::patchID
const labelList & patchID() const
Per boundary face label the patch index.
Definition: polyBoundaryMesh.C:399
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::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
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::FatalError
error FatalError
Foam::topoSetSource::CELLSETSOURCE
@ CELLSETSOURCE
Definition: topoSetSource.H:72
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::plusEqOp
Definition: ops.H:71
Foam::topoSetSource::New
static autoPtr< topoSetSource > New(const word &topoSetSourceType, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected topoSetSource.
Definition: topoSetSource.C:74
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:48
Foam::entry::dict
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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:282
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
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::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:70
topoSetSource.H
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
patchi
label patchi
Definition: getPatchFieldScalar.H:1
createTime.H
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
cellSet.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:417
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259