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