foamUpgradeCyclics.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  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  foamUpgradeCyclics
29 
30 Group
31  grpPreProcessingUtilities
32 
33 Description
34  Tool to upgrade mesh and fields for split cyclics.
35 
36 Usage
37  \b foamUpgradeCyclics [OPTION]
38 
39  Options:
40  - \par -dry-run
41  Suppress writing the updated files with split cyclics
42 
43  - \par -enableFunctionEntries
44  By default all dictionary preprocessing of fields is disabled
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "argList.H"
49 #include "Time.H"
50 #include "timeSelector.H"
51 #include "IOdictionary.H"
52 #include "polyMesh.H"
53 #include "entry.H"
54 #include "IOPtrList.H"
55 #include "cyclicPolyPatch.H"
56 #include "dictionaryEntry.H"
57 #include "IOobjectList.H"
58 #include "volFields.H"
59 #include "pointFields.H"
60 #include "surfaceFields.H"
61 
62 using namespace Foam;
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
69 }
70 
71 
72 // Read boundary file without reading mesh
73 void rewriteBoundary
74 (
75  const bool dryrun,
76  const IOobject& io,
77  const fileName& regionPrefix,
78  HashTable<word>& thisNames,
79  HashTable<word>& nbrNames
80 )
81 {
82  Info<< "Reading boundary from " << typeFilePath<IOPtrList<entry>>(io)
83  << endl;
84 
85  // Read PtrList of dictionary.
86  const word oldTypeName = IOPtrList<entry>::typeName;
89  const_cast<word&>(IOPtrList<entry>::typeName) = oldTypeName;
90  // Fake type back to what was in field
91  const_cast<word&>(patches.type()) = patches.headerClassName();
92 
93 
94  // Replace any 'cyclic'
95  label nOldCyclics = 0;
96  forAll(patches, patchi)
97  {
98  const dictionary& patchDict = patches[patchi].dict();
99 
100  if (patchDict.get<word>("type") == cyclicPolyPatch::typeName)
101  {
102  if (!patchDict.found("neighbourPatch"))
103  {
104  Info<< "Patch " << patches[patchi].keyword()
105  << " does not have 'neighbourPatch' entry; assuming it"
106  << " is of the old type." << endl;
107  nOldCyclics++;
108  }
109  }
110  }
111 
112  Info<< "Detected " << nOldCyclics << " old cyclics." << nl << endl;
113 
114 
115  // Save old patches.
116  PtrList<entry> oldPatches(patches);
117 
118  // Extend
119  label nOldPatches = patches.size();
120  patches.setSize(nOldPatches+nOldCyclics);
121 
122  // Create reordering map
123  labelList oldToNew(patches.size());
124 
125 
126  // Add new entries
127  label addedPatchi = nOldPatches;
128  label newPatchi = 0;
129  forAll(oldPatches, patchi)
130  {
131  const dictionary& patchDict = oldPatches[patchi].dict();
132 
133  if
134  (
135  patchDict.get<word>("type") == cyclicPolyPatch::typeName
136  )
137  {
138  const word& name = oldPatches[patchi].keyword();
139 
140  if (patchDict.found("neighbourPatch"))
141  {
142  patches.set(patchi, oldPatches.set(patchi, nullptr));
143  oldToNew[patchi] = newPatchi++;
144 
145  // Check if patches come from automatic conversion
146  word oldName;
147 
148  string::size_type i = name.rfind("_half0");
149  if (i != string::npos)
150  {
151  oldName = name.substr(0, i);
152  thisNames.insert(oldName, name);
153  Info<< "Detected converted cyclic patch " << name
154  << " ; assuming it originates from " << oldName
155  << endl;
156  }
157  else
158  {
159  i = name.rfind("_half1");
160  if (i != string::npos)
161  {
162  oldName = name.substr(0, i);
163  nbrNames.insert(oldName, name);
164  Info<< "Detected converted cyclic patch " << name
165  << " ; assuming it originates from " << oldName
166  << endl;
167  }
168  }
169  }
170  else
171  {
172  label nFaces = patchDict.get<label>("nFaces");
173  label startFace = patchDict.get<label>("startFace");
174 
175  Info<< "Detected old style " << patchDict.get<word>("type")
176  << " patch " << name << " with" << nl
177  << " nFaces : " << nFaces << nl
178  << " startFace : " << startFace << endl;
179 
180  word thisName = name + "_half0";
181  word nbrName = name + "_half1";
182 
183  thisNames.insert(name, thisName);
184  nbrNames.insert(name, nbrName);
185 
186  // Save current dictionary
187  const dictionary patchDict(patches[patchi].dict());
188 
189  // Change entry on this side
190  patches.set(patchi, oldPatches.set(patchi, nullptr));
191  oldToNew[patchi] = newPatchi++;
192  dictionary& thisPatchDict = patches[patchi].dict();
193  thisPatchDict.add("neighbourPatch", nbrName);
194  thisPatchDict.set("nFaces", nFaces/2);
195  patches[patchi].keyword() = thisName;
196 
197  // Add entry on other side
198  patches.set
199  (
200  addedPatchi,
201  new dictionaryEntry
202  (
203  nbrName,
205  patchDict
206  )
207  );
208  oldToNew[addedPatchi] = newPatchi++;
209  dictionary& nbrPatchDict = patches[addedPatchi].dict();
210  nbrPatchDict.set("neighbourPatch", thisName);
211  nbrPatchDict.set("nFaces", nFaces/2);
212  nbrPatchDict.set("startFace", startFace+nFaces/2);
213  patches[addedPatchi].keyword() = nbrName;
214 
215  Info<< "Replaced with patches" << nl
216  << patches[patchi].keyword() << " with" << nl
217  << " nFaces : "
218  << thisPatchDict.get<label>("nFaces") << nl
219  << " startFace : "
220  << thisPatchDict.get<label>("startFace") << nl
221  << patches[addedPatchi].keyword() << " with" << nl
222  << " nFaces : "
223  << nbrPatchDict.get<label>("nFaces") << nl
224  << " startFace : "
225  << nbrPatchDict.get<label>("startFace") << nl
226  << endl;
227 
228  addedPatchi++;
229  }
230  }
231  else
232  {
233  patches.set(patchi, oldPatches.set(patchi, nullptr));
234  oldToNew[patchi] = newPatchi++;
235  }
236  }
237 
238  patches.reorder(oldToNew);
239 
240  if (returnReduceOr(nOldCyclics))
241  {
242  if (dryrun)
243  {
244  //Info<< "-dry-run option: no changes made" << nl << endl;
245  }
246  else
247  {
248  if (mvBak(patches.objectPath(), "old"))
249  {
250  Info<< "Backup to "
251  << (patches.objectPath() + ".old") << nl;
252  }
253 
254  Info<< "Write to "
255  << patches.objectPath() << nl << endl;
256  patches.write();
257  }
258  }
259  else
260  {
261  Info<< "No changes made to boundary file." << nl << endl;
262  }
263 }
264 
265 
266 void rewriteField
267 (
268  const bool dryrun,
269  const Time& runTime,
270  const word& fieldName,
271  const HashTable<word>& thisNames,
272  const HashTable<word>& nbrNames
273 )
274 {
275  // Read dictionary. (disable class type checking so we can load
276  // field)
277  Info<< "Loading field " << fieldName << endl;
278  const word oldTypeName = IOdictionary::typeName;
279  const_cast<word&>(IOdictionary::typeName) = word::null;
280 
281  IOdictionary fieldDict
282  (
283  IOobject
284  (
285  fieldName,
286  runTime.timeName(),
287  runTime,
291  )
292  );
293  const_cast<word&>(IOdictionary::typeName) = oldTypeName;
294  // Fake type back to what was in field
295  const_cast<word&>(fieldDict.type()) = fieldDict.headerClassName();
296 
297 
298 
299  dictionary& boundaryField = fieldDict.subDict("boundaryField");
300 
301  bool hasChange = false;
302 
303  forAllConstIters(thisNames, iter)
304  {
305  const word& patchName = iter.key();
306  const word& newName = iter.val();
307 
308  Info<< "Looking for entry for patch " << patchName << endl;
309 
310  // Find old patch name either direct or through wildcards
311  // Find new patch name direct only
312 
313  if
314  (
315  boundaryField.found(patchName)
316  && !boundaryField.found(newName, keyType::LITERAL)
317  )
318  {
319  Info<< " Changing entry " << patchName << " to " << newName
320  << endl;
321 
322  dictionary& patchDict = boundaryField.subDict(patchName);
323 
324  if (patchDict.found("value"))
325  {
326  // Remove any value field since wrong size.
327  patchDict.remove("value");
328  }
329 
330 
331  boundaryField.changeKeyword(patchName, newName);
332  boundaryField.add
333  (
334  nbrNames[patchName],
335  patchDict
336  );
337  Info<< " Adding entry " << nbrNames[patchName] << endl;
338 
339  hasChange = true;
340  }
341  }
342 
343  //Info<< "New boundaryField:" << boundaryField << endl;
344 
345  if (returnReduceOr(hasChange))
346  {
347  if (dryrun)
348  {
349  //Info<< "-test option: no changes made" << endl;
350  }
351  else
352  {
353  if (mvBak(fieldDict.objectPath(), "old"))
354  {
355  Info<< "Backup to "
356  << (fieldDict.objectPath() + ".old") << nl;
357  }
358 
359  Info<< "Write to "
360  << fieldDict.objectPath() << endl;
361  fieldDict.regIOobject::write();
362  }
363  }
364  else
365  {
366  Info<< "No changes made to field " << fieldName << endl;
367  }
368  Info<< endl;
369 }
370 
371 
372 void rewriteFields
373 (
374  const bool dryrun,
375  const Time& runTime,
376  const wordList& fieldNames,
377  const HashTable<word>& thisNames,
378  const HashTable<word>& nbrNames
379 )
380 {
381  for (const word& fieldName : fieldNames)
382  {
383  rewriteField
384  (
385  dryrun,
386  runTime,
387  fieldName,
388  thisNames,
389  nbrNames
390  );
391  }
392 }
393 
394 
395 
396 int main(int argc, char *argv[])
397 {
399  (
400  "Tool to upgrade mesh and fields for split cyclics"
401  );
402 
404 
406  (
407  "Test only do not change any files"
408  );
410  (
411  "enableFunctionEntries",
412  "Enable expansion of dictionary directives - #include, #codeStream etc"
413  );
414  #include "addRegionOption.H"
415 
416  #include "setRootCase.H"
417  #include "createTime.H"
418 
419 
420  // Make sure we do not use the master-only reading since we read
421  // fields (different per processor) as dictionaries.
423 
424 
426 
427  const bool dryrun = args.dryRun();
428  if (dryrun)
429  {
430  Info<< "-dry-run option: no changes made" << nl << endl;
431  }
432  const bool enableEntries = args.found("enableFunctionEntries");
433 
434  const word regionName =
436 
437  fileName regionPrefix;
439  {
440  regionPrefix = regionName;
441  }
442 
443 
444  // Per cyclic patch the new name for this side and the other side
445  HashTable<word> thisNames;
446  HashTable<word> nbrNames;
447 
448  // Rewrite constant boundary file. Return any patches that have been split.
449  IOobject io
450  (
451  "boundary",
452  runTime.constant(),
454  runTime,
458  );
459 
460  if (io.typeHeaderOk<IOPtrList<entry>>(false))
461  {
462  rewriteBoundary
463  (
464  dryrun,
465  io,
466  regionPrefix,
467  thisNames,
468  nbrNames
469  );
470  }
471 
472 
473 
474  // Convert any fields
475 
476  forAll(timeDirs, timeI)
477  {
478  runTime.setTime(timeDirs[timeI], timeI);
479 
480  Info<< "Time: " << runTime.timeName() << endl;
481 
482  // See if mesh in time directory
483  IOobject io
484  (
485  "boundary",
486  runTime.timeName(),
488  runTime,
492  );
493 
494  if (io.typeHeaderOk<IOPtrList<entry>>(false))
495  {
496  rewriteBoundary
497  (
498  dryrun,
499  io,
500  regionPrefix,
501  thisNames,
502  nbrNames
503  );
504  }
505 
506 
507  IOobjectList objects(runTime, runTime.timeName());
508 
509 
510  const int oldFlag = entry::disableFunctionEntries;
511  if (!enableEntries)
512  {
513  // By default disable dictionary expansion for fields
515  }
516 
517  // volFields
518  // ~~~~~~~~~
519 
520  rewriteFields
521  (
522  dryrun,
523  runTime,
525  thisNames,
526  nbrNames
527  );
528  rewriteFields
529  (
530  dryrun,
531  runTime,
533  thisNames,
534  nbrNames
535  );
536  rewriteFields
537  (
538  dryrun,
539  runTime,
541  thisNames,
542  nbrNames
543  );
544  rewriteFields
545  (
546  dryrun,
547  runTime,
549  thisNames,
550  nbrNames
551  );
552  rewriteFields
553  (
554  dryrun,
555  runTime,
557  thisNames,
558  nbrNames
559  );
560 
561 
562  // pointFields
563  // ~~~~~~~~~~~
564 
565  rewriteFields
566  (
567  dryrun,
568  runTime,
570  thisNames,
571  nbrNames
572  );
573  rewriteFields
574  (
575  dryrun,
576  runTime,
578  thisNames,
579  nbrNames
580  );
581  rewriteFields
582  (
583  dryrun,
584  runTime,
586  thisNames,
587  nbrNames
588  );
589  rewriteFields
590  (
591  dryrun,
592  runTime,
594  thisNames,
595  nbrNames
596  );
597  rewriteFields
598  (
599  dryrun,
600  runTime,
602  thisNames,
603  nbrNames
604  );
605 
606 
607  // surfaceFields
608  // ~~~~~~~~~~~
609 
610  rewriteFields
611  (
612  dryrun,
613  runTime,
615  thisNames,
616  nbrNames
617  );
618  rewriteFields
619  (
620  dryrun,
621  runTime,
623  thisNames,
624  nbrNames
625  );
626  rewriteFields
627  (
628  dryrun,
629  runTime,
631  thisNames,
632  nbrNames
633  );
634  rewriteFields
635  (
636  dryrun,
637  runTime,
639  thisNames,
640  nbrNames
641  );
642  rewriteFields
643  (
644  dryrun,
645  runTime,
647  thisNames,
648  nbrNames
649  );
650 
652  }
653 
654  return 0;
655 }
656 
657 
658 // ************************************************************************* //
Foam::surfaceFields.
bool mvBak(const fileName &src, const std::string &ext="bak")
Rename to a corresponding backup file.
Definition: POSIX.C:1357
dictionary dict
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
A class for handling file names.
Definition: fileName.H:71
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
wordList names() const
The unsorted names of all objects.
void reorder(const labelUList &oldToNew, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:637
Ignore writing from objectRegistry::writeObject()
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
A keyword and a list of tokens is a &#39;dictionaryEntry&#39;.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:453
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:251
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:100
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
bool remove(const word &keyword)
Remove an entry specified by keyword.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
Foam::word regionName(Foam::polyMesh::defaultRegion)
static int disableFunctionEntries
Enable or disable use of function entries and variable expansions.
Definition: entry.H:139
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
A class for handling words, derived from Foam::string.
Definition: word.H:63
static void addDryRunOption(const string &usage, bool advanced=false)
Enable a &#39;dry-run&#39; bool option, with usage information.
Definition: argList.C:504
int dryRun() const noexcept
Return the dry-run flag.
Definition: argListI.H:109
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:397
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:113
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:465
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static const word null
An empty word.
Definition: word.H:84
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:997
String literal.
Definition: keyType.H:82
A HashTable similar to std::unordered_map.
Definition: HashTable.H:102
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:68
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:326
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
Definition: timeSelector.C:234
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:316
const word typeName("volScalarField::Internal")
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:49
const word & headerClassName() const noexcept
Return name of the class name read from header.
Definition: IOobjectI.H:180
const polyBoundaryMesh & patches
messageStream Info
Information stream (stdout output on master, null elsewhere)
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:777
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
bool changeKeyword(const keyType &oldKeyword, const keyType &newKeyword, bool overwrite=false)
Change the keyword for an entry,.
Foam::argList args(argc, argv)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
Do not request registration (bool: false)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28