ENH: changeDictionary: understand patchGroups

This commit is contained in:
mattijs
2012-08-03 17:36:48 +01:00
parent e316a4d647
commit 62c5dd128a

View File

@ -58,11 +58,15 @@ Usage
- changeDictionary [OPTION] - changeDictionary [OPTION]
\param -literalRE \n \param -literalRE \n
Do not interpret regular expressions; treat them as any other keyword. Do not interpret regular expressions or patchGroups;
treat them as any other keyword.
\param -enableFunctionEntries \n \param -enableFunctionEntries \n
By default all dictionary preprocessing of fields is disabled By default all dictionary preprocessing of fields is disabled
\param -disablePatchGroups \n
By default all keys are also checked for being patchGroups
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "argList.H" #include "argList.H"
@ -84,7 +88,47 @@ namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
bool merge(dictionary&, const dictionary&, const bool); // Extract groupPatch (= shortcut) info from boundary file info
HashTable<wordList, word> extractPatchGroups(const dictionary& boundaryDict)
{
HashTable<wordList, word> groupToPatch;
forAllConstIter(dictionary, boundaryDict, iter)
{
const word& patchName = iter().keyword();
const dictionary& patchDict = iter().dict();
wordList groups;
if (patchDict.readIfPresent("inGroups", groups))
{
forAll(groups, i)
{
HashTable<wordList, word>::iterator fndGroup = groupToPatch.find
(
groups[i]
);
if (fndGroup == groupToPatch.end())
{
groupToPatch.insert(groups[i], wordList(1, patchName));
}
else
{
fndGroup().append(patchName);
}
}
}
}
return groupToPatch;
}
bool merge
(
dictionary&,
const dictionary&,
const bool,
const HashTable<wordList, word>&
);
// Add thisEntry to dictionary thisDict. // Add thisEntry to dictionary thisDict.
@ -93,7 +137,8 @@ bool addEntry
dictionary& thisDict, dictionary& thisDict,
entry& thisEntry, entry& thisEntry,
const entry& mergeEntry, const entry& mergeEntry,
const bool literalRE const bool literalRE,
const HashTable<wordList, word>& shortcuts
) )
{ {
bool changed = false; bool changed = false;
@ -108,7 +153,8 @@ bool addEntry
( (
const_cast<dictionary&>(thisEntry.dict()), const_cast<dictionary&>(thisEntry.dict()),
mergeEntry.dict(), mergeEntry.dict(),
literalRE literalRE,
shortcuts
) )
) )
{ {
@ -135,7 +181,8 @@ bool merge
( (
dictionary& thisDict, dictionary& thisDict,
const dictionary& mergeDict, const dictionary& mergeDict,
const bool literalRE const bool literalRE,
const HashTable<wordList, word>& shortcuts
) )
{ {
bool changed = false; bool changed = false;
@ -156,7 +203,7 @@ bool merge
{ {
const keyType& key = mergeIter().keyword(); const keyType& key = mergeIter().keyword();
if (literalRE || !key.isPattern()) if (literalRE || !(key.isPattern() || shortcuts.found(key)))
{ {
entry* entryPtr = thisDict.lookupEntryPtr entry* entryPtr = thisDict.lookupEntryPtr
( (
@ -179,7 +226,8 @@ bool merge
thisDict, thisDict,
*entryPtr, *entryPtr,
mergeIter(), mergeIter(),
literalRE literalRE,
shortcuts
) )
) )
{ {
@ -196,29 +244,54 @@ bool merge
} }
// Pass 2. Wildcard matches (if any) on any non-match keys. // Pass 2. Wildcard or shortcut matches (if any) on any non-match keys.
if (!literalRE && thisKeysSet.size() > 0) if (!literalRE && thisKeysSet.size() > 0)
{ {
// Pick up remaining dictionary entries
wordList thisKeys(thisKeysSet.toc()); wordList thisKeys(thisKeysSet.toc());
forAllConstIter(IDLList<entry>, mergeDict, mergeIter) forAllConstIter(IDLList<entry>, mergeDict, mergeIter)
{ {
const keyType& key = mergeIter().keyword(); const keyType& key = mergeIter().keyword();
// List of indices into thisKeys
labelList matches;
if (key.isPattern()) if (key.isPattern())
{ {
// Find all matching entries in the original thisDict // Wildcard match
matches = findStrings(key, thisKeys);
labelList matches = findStrings(key, thisKeys); }
else if (shortcuts.size())
{
// See if patchGroups expand to valid thisKeys
const wordList shortcutNames = shortcuts.toc();
labelList indices = findStrings(key, shortcutNames);
forAll(indices, i)
{
const word& name = shortcutNames[indices[i]];
const wordList& keys = shortcuts[name];
forAll(keys, j)
{
label index = findIndex(thisKeys, keys[j]);
if (index != -1)
{
matches.append(index);
}
}
}
}
// Add all matches
forAll(matches, i) forAll(matches, i)
{ {
label matchI = matches[i]; const word& thisKey = thisKeys[matches[i]];
entry& thisEntry = const_cast<entry&> entry& thisEntry = const_cast<entry&>
( (
thisDict.lookupEntry(thisKeys[matchI], false, false) thisDict.lookupEntry(thisKey, false, false)
); );
if if
@ -228,7 +301,9 @@ bool merge
thisDict, thisDict,
thisEntry, thisEntry,
mergeIter(), mergeIter(),
literalRE literalRE,
HashTable<wordList, word>(0) // no shortcuts
// at deeper levels
) )
) )
{ {
@ -237,7 +312,6 @@ bool merge
} }
} }
} }
}
return changed; return changed;
} }
@ -273,6 +347,12 @@ int main(int argc, char *argv[])
"enableFunctionEntries", "enableFunctionEntries",
"enable expansion of dictionary directives - #include, #codeStream etc" "enable expansion of dictionary directives - #include, #codeStream etc"
); );
argList::addBoolOption
(
"disablePatchGroups",
"disable matching keys to patch groups"
);
#include "addRegionOption.H" #include "addRegionOption.H"
#include "setRootCase.H" #include "setRootCase.H"
@ -304,7 +384,6 @@ int main(int argc, char *argv[])
} }
const bool literalRE = args.optionFound("literalRE"); const bool literalRE = args.optionFound("literalRE");
if (literalRE) if (literalRE)
{ {
Info<< "Not interpreting any regular expressions (RE)" Info<< "Not interpreting any regular expressions (RE)"
@ -328,6 +407,15 @@ int main(int argc, char *argv[])
} }
const bool disablePatchGroups = args.optionFound("disablePatchGroups");
if (disablePatchGroups)
{
Info<< "Not interpreting any keys in the changeDictionary"
<< " as patchGroups"
<< endl;
}
fileName regionPrefix = ""; fileName regionPrefix = "";
if (regionName != fvMesh::defaultRegion) if (regionName != fvMesh::defaultRegion)
{ {
@ -369,20 +457,12 @@ int main(int argc, char *argv[])
<< replaceDicts.toc() << endl; << replaceDicts.toc() << endl;
// Every replacement is a dictionary name and a keyword in this
forAllConstIter(dictionary, replaceDicts, fieldIter) // Always read boundary to get patch groups
{ // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
const word& fieldName = fieldIter().keyword();
Info<< "Replacing entries in dictionary " << fieldName << endl;
// Handle 'boundary' specially: Info<< "Reading polyMesh/boundary file to extract patch names"
// - is PtrList of dictionaries << endl;
// - is in polyMesh/
if (fieldName == "boundary")
{
Info<< "Special handling of " << fieldName
<< " as polyMesh/boundary file." << endl;
// Read PtrList of dictionary as dictionary. // Read PtrList of dictionary as dictionary.
const word oldTypeName = IOPtrList<entry>::typeName; const word oldTypeName = IOPtrList<entry>::typeName;
@ -391,15 +471,16 @@ int main(int argc, char *argv[])
( (
IOobject IOobject
( (
fieldName, "boundary",
runTime.findInstance runTime.findInstance
( (
regionPrefix/polyMesh::meshSubDir, regionPrefix/polyMesh::meshSubDir,
fieldName "boundary",
IOobject::READ_IF_PRESENT
), ),
polyMesh::meshSubDir, polyMesh::meshSubDir,
mesh, mesh,
IOobject::MUST_READ, IOobject::READ_IF_PRESENT,
IOobject::NO_WRITE, IOobject::NO_WRITE,
false false
) )
@ -415,15 +496,52 @@ int main(int argc, char *argv[])
fieldDict.add(dictList[i].keyword(), dictList[i].dict()); fieldDict.add(dictList[i].keyword(), dictList[i].dict());
} }
Info<< "Loaded dictionary " << fieldName if (dictList.size())
{
Info<< "Loaded dictionary " << dictList.name()
<< " with entries " << fieldDict.toc() << endl; << " with entries " << fieldDict.toc() << endl;
}
// Extract any patchGroups information (= shortcut for set of
// patches)
HashTable<wordList, word> patchGroups;
if (!disablePatchGroups)
{
patchGroups = extractPatchGroups(fieldDict);
if (patchGroups.size())
{
Info<< "Extracted patch groups:" << endl;
wordList groups(patchGroups.sortedToc());
forAll(groups, i)
{
Info<< " group " << groups[i] << " with patches "
<< patchGroups[groups[i]] << endl;
}
}
}
// Every replacement is a dictionary name and a keyword in this
forAllConstIter(dictionary, replaceDicts, fieldIter)
{
const word& fieldName = fieldIter().keyword();
Info<< "Replacing entries in dictionary " << fieldName << endl;
// Handle 'boundary' specially:
// - is PtrList of dictionaries
// - is in polyMesh/
if (fieldName == "boundary")
{
Info<< "Special handling of " << fieldName
<< " as polyMesh/boundary file." << endl;
// Get the replacement dictionary for the field // Get the replacement dictionary for the field
const dictionary& replaceDict = fieldIter().dict(); const dictionary& replaceDict = fieldIter().dict();
Info<< "Merging entries from " << replaceDict.toc() << endl; Info<< "Merging entries from " << replaceDict.toc() << endl;
// Merge the replacements in // Merge the replacements in
merge(fieldDict, replaceDict, literalRE); merge(fieldDict, replaceDict, literalRE, patchGroups);
Info<< "fieldDict:" << fieldDict << endl; Info<< "fieldDict:" << fieldDict << endl;
@ -492,7 +610,7 @@ int main(int argc, char *argv[])
Info<< "Merging entries from " << replaceDict.toc() << endl; Info<< "Merging entries from " << replaceDict.toc() << endl;
// Merge the replacements in // Merge the replacements in
merge(fieldDict, replaceDict, literalRE); merge(fieldDict, replaceDict, literalRE, patchGroups);
Info<< "Writing modified fieldDict " << fieldName << endl; Info<< "Writing modified fieldDict " << fieldName << endl;
fieldDict.regIOobject::write(); fieldDict.regIOobject::write();