Loading...
Searching...
No Matches
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-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Application
28 foamUpgradeCyclics
29
30Group
31 grpPreProcessingUtilities
32
33Description
34 Tool to upgrade mesh and fields for split cyclics.
35
36Usage
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
62using namespace Foam;
63
64// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65
66namespace Foam
67{
69}
70
71
72// Read boundary file without reading mesh
73void 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 "
83 << io.typeFilePath<IOPtrList<entry>>() << nl;
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.release(patchi));
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.release(patchi));
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,
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.release(patchi));
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
266void 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 (
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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373
374int main(int argc, char *argv[])
375{
377 (
378 "Tool to upgrade mesh and fields for split cyclics"
379 );
380
382
384 (
385 "Test only do not change any files"
386 );
388 (
389 "enableFunctionEntries",
390 "Enable expansion of dictionary directives - #include, #codeStream etc"
391 );
392 #include "addRegionOption.H"
393
394 #include "setRootCase.H"
395 #include "createTime.H"
396
397
398 // Make sure we do not use the master-only reading since we read
399 // fields (different per processor) as dictionaries.
401
402
404
405 const bool dryrun = args.dryRun();
406 if (dryrun)
407 {
408 Info<< "-dry-run option: no changes made" << nl << endl;
409 }
410 const bool enableEntries = args.found("enableFunctionEntries");
411
412 fileName regionPrefix;
413 {
414 // Specified region or default region
415 #include "getRegionOption.H"
416
417 regionPrefix = polyMesh::regionName(regionName);
418 }
419
420 // Per cyclic patch the new name for this side and the other side
421 HashTable<word> thisNames;
422 HashTable<word> nbrNames;
423
424 // Rewrite constant boundary file. Return any patches that have been split.
426 (
427 "boundary",
428 runTime.constant(),
430 runTime,
434 );
435
436 if (io.typeHeaderOk<IOPtrList<entry>>(false))
437 {
438 rewriteBoundary
439 (
440 dryrun,
441 io,
442 regionPrefix,
443 thisNames,
444 nbrNames
445 );
446 }
447
448
449
450 // Convert any fields
451
452 forAll(timeDirs, timei)
453 {
454 runTime.setTime(timeDirs[timei], timei);
455
456 Info<< "Time: " << runTime.timeName() << endl;
457
458 // See if mesh in time directory
460 (
461 "boundary",
462 runTime.timeName(),
464 runTime,
468 );
469
470 if (io.typeHeaderOk<IOPtrList<entry>>(false))
471 {
472 rewriteBoundary
473 (
474 dryrun,
475 io,
476 regionPrefix,
477 thisNames,
478 nbrNames
479 );
480 }
481
482
483 IOobjectList objects(runTime, runTime.timeName());
484
485
486 const int oldFlag = entry::disableFunctionEntries;
487 if (!enableEntries)
488 {
489 // By default disable dictionary expansion for fields
491 }
492
493 // Rewriting of fields:
494
495 #undef rewriteFields
496 #define rewriteFields(FieldType) \
497 for (const word& fieldName : objects.sortedNames<FieldType>()) \
498 { \
499 rewriteField(dryrun, runTime, fieldName, thisNames, nbrNames); \
500 }
501
502 // volFields
503 // ~~~~~~~~~
504
505 rewriteFields(volScalarField);
506 rewriteFields(volVectorField);
507 rewriteFields(volSphericalTensorField);
508 rewriteFields(volSymmTensorField);
509 rewriteFields(volTensorField);
510
511 // pointFields
512 // ~~~~~~~~~~~
513
514 rewriteFields(pointScalarField);
515 rewriteFields(pointVectorField);
516 rewriteFields(pointSphericalTensorField);
517 rewriteFields(pointSymmTensorField);
518 rewriteFields(pointTensorField);
519
520
521 // surfaceFields
522 // ~~~~~~~~~~~~~
523
524 rewriteFields(surfaceScalarField);
525 rewriteFields(surfaceVectorField);
526 rewriteFields(surfaceSphericalTensorField);
527 rewriteFields(surfaceSymmTensorField);
528 rewriteFields(surfaceTensorField);
529
530 #undef rewriteFields
531
533 }
534
535 return 0;
536}
537
538
539// ************************************************************************* //
Required Classes.
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
A PtrList of objects of type <T> with automated input and output.
Definition IOPtrList.H:53
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition IOobject.H:358
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static void addDryRunOption(const string &usage, bool advanced=false)
Enable a 'dry-run' bool option, with usage information.
Definition argList.C:519
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
A keyword and a list of tokens is a 'dictionaryEntry'.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
bool changeKeyword(const keyType &oldKeyword, const keyType &newKeyword, bool overwrite=false)
Change the keyword for an entry,.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool remove(const word &keyword)
Remove an entry specified by keyword.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
static int disableFunctionEntries
Enable or disable use of function entries and variable expansions.
Definition entry.H:139
A class for handling file names.
Definition fileName.H:75
@ LITERAL
String literal.
Definition keyType.H:82
static const word & regionName(const word &region)
The mesh region name or word::null if polyMesh::defaultRegion.
Definition polyMesh.C:796
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
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...
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
#define defineTemplateTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information for templates, useful.
Definition className.H:158
const polyBoundaryMesh & patches
engineTime & runTime
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
const auto & io
auto & name
Required Classes.
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< tensor, pointPatchField, pointMesh > pointTensorField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< sphericalTensor, pointPatchField, pointMesh > pointSphericalTensorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
List< instant > instantList
List of instants.
Definition instantList.H:41
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
GeometricField< symmTensor, pointPatchField, pointMesh > pointSymmTensorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
GeometricField< tensor, fvsPatchField, surfaceMesh > surfaceTensorField
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
GeometricField< sphericalTensor, fvsPatchField, surfaceMesh > surfaceSphericalTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
bool mvBak(const fileName &src, const std::string &ext="bak")
Rename to a corresponding backup file.
Definition POSIX.C:1359
GeometricField< symmTensor, fvsPatchField, surfaceMesh > surfaceSymmTensorField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Foam::surfaceFields.