diff --git a/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.C b/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.C index a7e5e1ee09..9544c8799b 100644 --- a/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.C +++ b/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.C @@ -53,34 +53,32 @@ void Foam::wallShearStress::calcShearStress volVectorField& shearStress ) { - forAll(shearStress.boundaryField(), patchI) + forAllConstIter(labelHashSet, patchSet_, iter) { + label patchI = iter.key(); const polyPatch& pp = mesh.boundaryMesh()[patchI]; - if (isA(pp)) + vectorField& ssp = shearStress.boundaryField()[patchI]; + const vectorField& Sfp = mesh.Sf().boundaryField()[patchI]; + const scalarField& magSfp = mesh.magSf().boundaryField()[patchI]; + const symmTensorField& Reffp = Reff.boundaryField()[patchI]; + + ssp = (-Sfp/magSfp) & Reffp; + + vector minSsp = gMin(ssp); + vector maxSsp = gMax(ssp); + + if (Pstream::master()) { - vectorField& ssp = shearStress.boundaryField()[patchI]; - const vectorField& Sfp = mesh.Sf().boundaryField()[patchI]; - const scalarField& magSfp = mesh.magSf().boundaryField()[patchI]; - const symmTensorField& Reffp = Reff.boundaryField()[patchI]; + file() << mesh.time().timeName() << token::TAB + << pp.name() << token::TAB << minSsp + << token::TAB << maxSsp << endl; + } - ssp = (-Sfp/magSfp) & Reffp; - - vector minSsp = gMin(ssp); - vector maxSsp = gMax(ssp); - - if (Pstream::master()) - { - file() << mesh.time().timeName() << token::TAB - << pp.name() << token::TAB << minSsp - << token::TAB << maxSsp << endl; - } - - if (log_) - { - Info<< " min/max(" << pp.name() << ") = " - << minSsp << ", " << maxSsp << endl; - } + if (log_) + { + Info<< " min/max(" << pp.name() << ") = " + << minSsp << ", " << maxSsp << endl; } } } @@ -101,7 +99,8 @@ Foam::wallShearStress::wallShearStress obr_(obr), active_(true), log_(false), - phiName_("phi") + phiName_("phi"), + patchSet_() { // Check if the available mesh is an fvMesh, otherwise deactivate if (!isA(obr_)) @@ -138,6 +137,54 @@ void Foam::wallShearStress::read(const dictionary& dict) { log_ = dict.lookupOrDefault("log", false); phiName_ = dict.lookupOrDefault("phiName", "phi"); + + const fvMesh& mesh = refCast(obr_); + const polyBoundaryMesh& pbm = mesh.boundaryMesh(); + + patchSet_ = + mesh.boundaryMesh().patchSet + ( + wordReList(dict.lookupOrDefault("patches", wordReList())) + ); + + Info<< type() << " output:" << nl; + + if (patchSet_.empty()) + { + forAll(pbm, patchI) + { + if (isA(pbm[patchI])) + { + patchSet_.insert(patchI); + } + } + + Info<< " processing all wall patches" << nl << endl; + } + else + { + Info<< " processing wall patches: " << nl; + labelHashSet filteredPatchSet; + forAllConstIter(labelHashSet, patchSet_, iter) + { + label patchI = iter.key(); + if (isA(pbm[patchI])) + { + filteredPatchSet.insert(patchI); + Info<< " " << pbm[patchI].name() << endl; + } + else + { + WarningIn("void wallShearStress::read(const dictionary&)") + << "Requested wall shear stress on non-wall boundary " + << "type patch: " << pbm[patchI].name() << endl; + } + } + + Info<< endl; + + patchSet_ = filteredPatchSet; + } } } diff --git a/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.H b/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.H index 80f12d0c24..ac16eca594 100644 --- a/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.H +++ b/src/postProcessing/functionObjects/utilities/wallShearStress/wallShearStress.H @@ -36,15 +36,34 @@ Description Stress = R \dot n \f] - The shear stress (symmetrical) tensor field is retrieved from the - turbulence model. - where \vartable R | stress tensor n | patch normal vector (into the domain) \endvartable + The shear stress (symmetrical) tensor field is retrieved from the + turbulence model. All wall patches are included by default; to restrict + the calculation to certain patches, use the optional 'patches' entry. + + Example of function object specification: + \verbatim + wallShearStress1 + { + type wallShearStress; + functionObjectLibs ("libutilityFunctionObjects.so"); + ... + patches (".*Wall"); + } + \endverbatim + + \heading Function object usage + \table + Property | Description | Required | Default value + type | type name: wallShearStress | yes | + patches | list of patches to process | no | all wall patches + \endtable + SourceFiles wallShearStress.C IOwallShearStress.H @@ -97,6 +116,9 @@ protected: //- Name of mass/volume flux field (optional, default = phi) word phiName_; + //- Optional list of patches to process + labelHashSet patchSet_; + // Protected Member Functions