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-2023 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 "timeSelector.H"
31 #include "Time.H"
32 #include "fvMesh.H"
33 #include "topoSetSource.H"
34 #include "cellSet.H"
35 #include "faceSet.H"
36 #include "volFields.H"
37 #include "systemDict.H"
38 #include "processorFvPatch.H"
39 #include "CompactListList.H"
40 
41 using namespace Foam;
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 template<class Type>
46 bool setCellFieldType
47 (
48  const word& fieldTypeDesc,
49  const fvMesh& mesh,
50  const labelList& selectedCells,
51  Istream& fieldValueStream
52 )
53 {
54  if (fieldTypeDesc != VolField<Type>::typeName + "Value")
55  {
56  return false;
57  }
58 
59  word fieldName(fieldValueStream);
60 
61  // Check the current time directory
62  typeIOobject<VolField<Type>> fieldHeader
63  (
64  fieldName,
65  mesh.time().name(),
66  mesh,
68  );
69 
70  // Check the "constant" directory
71  if (!fieldHeader.headerOk())
72  {
73  fieldHeader = typeIOobject<VolField<Type>>
74  (
75  fieldName,
76  mesh.time().constant(),
77  mesh,
79  );
80  }
81 
82  // Check field exists
83  if (fieldHeader.headerOk())
84  {
85  Info<< " Setting internal values of "
86  << fieldHeader.headerClassName()
87  << " " << fieldName << endl;
88 
89  VolField<Type> field(fieldHeader, mesh);
90 
91  const Type value = pTraits<Type>(fieldValueStream);
92 
93  if (selectedCells.size() == field.size())
94  {
95  field.primitiveFieldRef() = value;
96  }
97  else
98  {
99  forAll(selectedCells, celli)
100  {
101  field[selectedCells[celli]] = value;
102  }
103  }
104 
105  typename VolField<Type>::
106  Boundary& fieldBf = field.boundaryFieldRef();
107 
108  forAll(field.boundaryField(), patchi)
109  {
110  fieldBf[patchi] = fieldBf[patchi].patchInternalField();
111  }
112 
113  if (!field.write())
114  {
116  << "Failed writing field " << fieldName << endl;
117  }
118  }
119  else
120  {
122  << "Field " << fieldName << " not found" << endl;
123 
124  // Consume value
125  (void)pTraits<Type>(fieldValueStream);
126  }
127 
128  return true;
129 }
130 
131 
132 class setCellField
133 {
134 
135 public:
136 
137  setCellField()
138  {}
139 
141  {
142  return autoPtr<setCellField>(new setCellField());
143  }
144 
145  class iNew
146  {
147  const fvMesh& mesh_;
148  const labelList& selectedCells_;
149 
150  public:
151 
152  iNew(const fvMesh& mesh, const labelList& selectedCells)
153  :
154  mesh_(mesh),
155  selectedCells_(selectedCells)
156  {}
157 
158  autoPtr<setCellField> operator()(Istream& fieldValues) const
159  {
160  word fieldType(fieldValues);
161 
162  if
163  (
164  !(
165  setCellFieldType<scalar>
166  (fieldType, mesh_, selectedCells_, fieldValues)
167  || setCellFieldType<vector>
168  (fieldType, mesh_, selectedCells_, fieldValues)
169  || setCellFieldType<sphericalTensor>
170  (fieldType, mesh_, selectedCells_, fieldValues)
171  || setCellFieldType<symmTensor>
172  (fieldType, mesh_, selectedCells_, fieldValues)
173  || setCellFieldType<tensor>
174  (fieldType, mesh_, selectedCells_, fieldValues)
175  )
176  )
177  {
179  << "field type " << fieldType << " not currently supported"
180  << endl;
181  }
182 
183  return autoPtr<setCellField>(new setCellField());
184  }
185  };
186 };
187 
188 
189 template<class Type>
190 bool setFaceFieldType
191 (
192  const word& fieldTypeDesc,
193  const fvMesh& mesh,
194  const labelList& selectedFaces,
195  Istream& fieldValueStream
196 )
197 {
198  if (fieldTypeDesc != VolField<Type>::typeName + "Value")
199  {
200  return false;
201  }
202 
203  word fieldName(fieldValueStream);
204 
205  // Check the current time directory
206  typeIOobject<VolField<Type>> fieldHeader
207  (
208  fieldName,
209  mesh.time().name(),
210  mesh,
212  );
213 
214  // Check the "constant" directory
215  if (!fieldHeader.headerOk())
216  {
217  fieldHeader = typeIOobject<VolField<Type>>
218  (
219  fieldName,
220  mesh.time().constant(),
221  mesh,
223  );
224  }
225 
226  // Check field exists
227  if (fieldHeader.headerOk())
228  {
229  Info<< " Setting patchField values of "
230  << fieldHeader.headerClassName()
231  << " " << fieldName << endl;
232 
233  // Read the field
234  VolField<Type> field(fieldHeader, mesh);
235  typename VolField<Type>::Boundary& fieldBf = field.boundaryFieldRef();
236 
237  // Read the value
238  const Type value = pTraits<Type>(fieldValueStream);
239 
240  // Determine the number of non-processor patches
241  label nNonProcPatches = 0;
242  forAll(fieldBf, patchi)
243  {
244  nNonProcPatches = patchi;
245  if (isA<processorFvPatch>(mesh.boundary()[patchi]))
246  {
247  break;
248  }
249  }
250 
251  // Create a copy of the boundary field
252  typename VolField<Type>::Boundary fieldBfCopy
253  (
255  fieldBf
256  );
257 
258  // Loop selected faces and set values in the copied boundary field
259  bool haveWarnedInternal = false, haveWarnedProc = false;
260  labelList nonProcPatchNChangedFaces(nNonProcPatches, 0);
261  forAll(selectedFaces, i)
262  {
263  const label facei = selectedFaces[i];
264 
265  if (mesh.isInternalFace(facei))
266  {
267  if (!haveWarnedInternal)
268  {
270  << "Ignoring internal face " << facei
271  << ". Suppressing further warnings." << endl;
272  haveWarnedInternal = true;
273  }
274  }
275  else
276  {
277  const labelUList patches =
278  mesh.polyBFacePatches()[facei - mesh.nInternalFaces()];
279  const labelUList patchFaces =
280  mesh.polyBFacePatchFaces()[facei - mesh.nInternalFaces()];
281 
282  forAll(patches, i)
283  {
284  if (patches[i] >= nNonProcPatches)
285  {
286  if (!haveWarnedProc)
287  {
289  << "Ignoring face " << patchFaces[i]
290  << " of processor patch " << patches[i]
291  << ". Suppressing further warnings." << endl;
292  haveWarnedProc = true;
293  }
294  }
295  else
296  {
297  fieldBfCopy[patches[i]][patchFaces[i]] = value;
298  nonProcPatchNChangedFaces[patches[i]] ++;
299  }
300  }
301  }
302  }
304  (
305  nonProcPatchNChangedFaces,
307  );
309  (
310  nonProcPatchNChangedFaces
311  );
312 
313  // Reassign boundary values
314  forAll(nonProcPatchNChangedFaces, patchi)
315  {
316  if (nonProcPatchNChangedFaces[patchi] > 0)
317  {
318  Info<< " On patch "
319  << field.boundaryField()[patchi].patch().name()
320  << " set " << nonProcPatchNChangedFaces[patchi]
321  << " values" << endl;
322  fieldBf[patchi] == fieldBfCopy[patchi];
323  }
324  }
325 
326  if (!field.write())
327  {
329  << "Failed writing field " << field.name() << exit(FatalError);
330  }
331  }
332  else
333  {
335  << "Field " << fieldName << " not found" << endl;
336 
337  // Consume value
338  (void)pTraits<Type>(fieldValueStream);
339  }
340 
341  return true;
342 }
343 
344 
345 class setFaceField
346 {
347 
348 public:
349 
350  setFaceField()
351  {}
352 
354  {
355  return autoPtr<setFaceField>(new setFaceField());
356  }
357 
358  class iNew
359  {
360  const fvMesh& mesh_;
361  const labelList& selectedFaces_;
362 
363  public:
364 
365  iNew(const fvMesh& mesh, const labelList& selectedFaces)
366  :
367  mesh_(mesh),
368  selectedFaces_(selectedFaces)
369  {}
370 
371  autoPtr<setFaceField> operator()(Istream& fieldValues) const
372  {
373  word fieldType(fieldValues);
374 
375  if
376  (
377  !(
378  setFaceFieldType<scalar>
379  (fieldType, mesh_, selectedFaces_, fieldValues)
380  || setFaceFieldType<vector>
381  (fieldType, mesh_, selectedFaces_, fieldValues)
382  || setFaceFieldType<sphericalTensor>
383  (fieldType, mesh_, selectedFaces_, fieldValues)
384  || setFaceFieldType<symmTensor>
385  (fieldType, mesh_, selectedFaces_, fieldValues)
386  || setFaceFieldType<tensor>
387  (fieldType, mesh_, selectedFaces_, fieldValues)
388  )
389  )
390  {
392  << "field type " << fieldType << " not currently supported"
393  << endl;
394  }
395 
396  return autoPtr<setFaceField>(new setFaceField());
397  }
398  };
399 };
400 
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 int main(int argc, char *argv[])
405 {
407  #include "addDictOption.H"
408  #include "addRegionOption.H"
409  #include "setRootCase.H"
410  #include "createTime.H"
411  timeSelector::select0(runTime, args);
412  #include "createNamedMesh.H"
413 
414  const dictionary setFieldsDict(systemDict("setFieldsDict", args, mesh));
415 
416  if (setFieldsDict.found("defaultFieldValues"))
417  {
418  Info<< "Setting field default values" << endl;
419  PtrList<setCellField> defaultFieldValues
420  (
421  setFieldsDict.lookup("defaultFieldValues"),
422  setCellField::iNew(mesh, labelList(mesh.nCells()))
423  );
424  Info<< endl;
425  }
426 
427  Info<< "Setting field region values" << endl;
428 
429  PtrList<entry> regions(setFieldsDict.lookup("regions"));
430 
431  forAll(regions, regionI)
432  {
433  const entry& region = regions[regionI];
434 
435  autoPtr<topoSetSource> source =
436  topoSetSource::New(region.keyword(), mesh, region.dict());
437 
438  if (source().setType() == topoSetSource::CELLSETSOURCE)
439  {
440  cellSet selectedCellSet
441  (
442  mesh,
443  "cellSet",
444  mesh.nCells()/10+1 // Reasonable size estimate.
445  );
446 
447  source->applyToSet
448  (
450  selectedCellSet
451  );
452 
453  PtrList<setCellField> fieldValues
454  (
455  region.dict().lookup("fieldValues"),
456  setCellField::iNew(mesh, selectedCellSet.toc())
457  );
458  }
459  else if (source().setType() == topoSetSource::FACESETSOURCE)
460  {
461  faceSet selectedFaceSet
462  (
463  mesh,
464  "faceSet",
465  (mesh.nFaces()-mesh.nInternalFaces())/10+1
466  );
467 
468  source->applyToSet
469  (
471  selectedFaceSet
472  );
473 
474  PtrList<setFaceField> fieldValues
475  (
476  region.dict().lookup("fieldValues"),
477  setFaceField::iNew(mesh, selectedFaceSet.toc())
478  );
479  }
480  }
481 
482 
483  Info<< "\nEnd\n" << endl;
484 
485  return 0;
486 }
487 
488 
489 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
Generic GeometricField class.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A collection of cell labels.
Definition: cellSet.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const word & name() const
Return const reference to name.
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:68
const keyType & keyword() const
Return keyword.
Definition: entry.H:123
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
A list of face labels.
Definition: faceSet.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:402
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:895
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:953
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:1029
label nCells() const
label nInternalFaces() const
label nFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
static autoPtr< topoSetSource > New(const word &topoSetSourceType, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected topoSetSource.
Definition: topoSetSource.C:70
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
T clone(const T &t)
Definition: List.H:55
error FatalError
Foam::argList args(argc, argv)