Skip to content
Snippets Groups Projects
surfaceFieldValueTemplates.C 14.6 KiB
Newer Older
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2015-2019 OpenCFD Ltd.
     \\/     M anipulation  |
-------------------------------------------------------------------------------
                            | Copyright (C) 2011-2017 OpenFOAM Foundation
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

#include "surfaceFields.H"
#include "volFields.H"
#include "interpolationCell.H"
#include "interpolationCellPoint.H"

// * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //

template<class WeightType>
inline bool Foam::functionObjects::fieldValues::surfaceFieldValue::canWeight
(
    const Field<WeightType>& weightField
) const
{
    return
    (
        usesWeight()
     && returnReduce(!weightField.empty(), orOp<bool>()) // On some processor
    );
}


template<class Type>
bool Foam::functionObjects::fieldValues::surfaceFieldValue::validField
{
    typedef GeometricField<Type, fvsPatchField, surfaceMesh> sf;
    typedef GeometricField<Type, fvPatchField, volMesh> vf;
    typedef DimensionedField<Type, polySurfaceGeoMesh> smt;

    return
    (
        foundObject<smt>(fieldName)
     || foundObject<vf>(fieldName)
     || (withSurfaceFields() && foundObject<sf>(fieldName))
template<class Type>
Foam::functionObjects::fieldValues::surfaceFieldValue::getFieldValues
) const
{
    typedef GeometricField<Type, fvsPatchField, surfaceMesh> sf;
    typedef GeometricField<Type, fvPatchField, volMesh> vf;
    typedef DimensionedField<Type, polySurfaceGeoMesh> smt;
    if (foundObject<smt>(fieldName))
        return lookupObject<smt>(fieldName);
    else if (withSurfaceFields() && foundObject<sf>(fieldName))
        return filterField(lookupObject<sf>(fieldName));
    }
    else if (foundObject<vf>(fieldName))
    {
        const vf& fld = lookupObject<vf>(fieldName);
            if (sampledPtr_().interpolate())
            {
                const interpolationCellPoint<Type> interp(fld);
                const interpolationCell<Type> interp(fld);

                return sampledPtr_().sample(interp);
        FatalErrorInFunction
            << "Field " << fieldName << " not found in database"
    return tmp<Field<Type>>::New();
template<class Type, class WeightType>
Type Foam::functionObjects::fieldValues::surfaceFieldValue::
processSameTypeValues
    const Field<Type>& values,
    const vectorField& Sf,
    const Field<WeightType>& weightField
    switch (operation_)
    {
        case opMin:
        {
            result = gMin(values);
            break;
        }
        case opMax:
        {
            result = gMax(values);
            break;
        }
        case opSumMag:
            result = gSum(cmptMag(values));
                tmp<scalarField> weight(weightingFactor(weightField));
                result = gSum(weight*values);
            }
            else
            {
                // Unweighted form
                result = gSum(values);
            }
        case opSumDirectionBalance:
        {
            FatalErrorInFunction
                << "Operation " << operationTypeNames_[operation_]
                << " not available for values of type "
                << pTraits<Type>::typeName
                << exit(FatalError);

            break;
        }
        case opWeightedAverage:
mattijs's avatar
mattijs committed
            {
                const scalarField factor(weightingFactor(weightField));

                result = gSum(factor*values)/(gSum(factor) + ROOTVSMALL);
mattijs's avatar
mattijs committed
            }
            else
            {
                // Unweighted form
                const label n = returnReduce(values.size(), sumOp<label>());
                result = gSum(values)/(scalar(n) + ROOTVSMALL);
            }
        case opAreaAverage:
        case opWeightedAreaAverage:
        case opAbsWeightedAreaAverage:
            if (canWeight(weightField))
            {
                const scalarField factor(weightingFactor(weightField, Sf));
                result = gSum(factor*values)/gSum(factor + ROOTVSMALL);
            }
            else
            {
                // Unweighted form
                const scalarField factor(mag(Sf));
                result = gSum(factor*values)/gSum(factor);
            }
            break;
        }
        case opAreaIntegrate:
        case opWeightedAreaIntegrate:
        case opAbsWeightedAreaIntegrate:
            if (canWeight(weightField))
            {
                tmp<scalarField> factor(weightingFactor(weightField, Sf));
                result = gSum(factor*values);
            }
            else
            {
                // Unweighted form
                tmp<scalarField> factor(mag(Sf));
                result = gSum(factor*values);
            }
            const scalarField magSf(mag(Sf));
            const scalar gSumMagSf = gSum(magSf);

            Type meanValue = gSum(values*magSf)/gSumMagSf;
            for (direction d=0; d < pTraits<Type>::nComponents; ++d)
                tmp<scalarField> vals(values.component(d));
                const scalar mean = component(meanValue, d);
                res =
                    sqrt(gSum(magSf*sqr(vals - mean))/gSumMagSf)
                   /(mean + ROOTVSMALL);
        case opAreaNormalAverage:
        case opAreaNormalIntegrate:

        case opWeightedUniformity:
        case opAbsWeightedUniformity:
        {
            if (canWeight(weightField))
            {
                // Change weighting from vector -> scalar and dispatch again
                return processValues<Type, scalar>
                (
                    values,
                    Sf,
                    weightingFactor(weightField)
                );
            }

            break;
        }
template<class Type, class WeightType>
Type Foam::functionObjects::fieldValues::surfaceFieldValue::processValues
(
    const Field<Type>& values,
    const vectorField& Sf,
    const Field<WeightType>& weightField
) const
{
    return processSameTypeValues(values, Sf, weightField);
}


template<class WeightType>
Foam::tmp<Foam::scalarField>
Foam::functionObjects::fieldValues::surfaceFieldValue::weightingFactor
(
    const Field<WeightType>& weightField
    // The scalar form is specialized.
    // For other types always need mag() to generate a scalar field.
template<class WeightType>
Foam::label Foam::functionObjects::fieldValues::surfaceFieldValue::writeAll
    const Field<WeightType>& weightField,
    const pointField& points,
    const faceList& faces
    for (const word& fieldName : fields_)
            writeValues<scalar>(fieldName, Sf, weightField, points, faces)
         || writeValues<vector>(fieldName, Sf, weightField, points, faces)
         || writeValues<sphericalTensor>
            (
                fieldName, Sf, weightField, points, faces
         || writeValues<symmTensor>(fieldName, Sf, weightField, points, faces)
         || writeValues<tensor>(fieldName, Sf, weightField, points, faces)
            WarningInFunction
                << "Requested field " << fieldName
                << " not found in database and not processed"
                << endl;
    }

    return nProcessed;
}


template<class Type, class WeightType>
bool Foam::functionObjects::fieldValues::surfaceFieldValue::writeValues
(
    const word& fieldName,
    const vectorField& Sf,
    const Field<WeightType>& weightField,
    const pointField& points,
    const faceList& faces
)
{
    const bool ok = validField<Type>(fieldName);

    if (ok)
    {
        Field<Type> values(getFieldValues<Type>(fieldName, true));
        // Write raw values on surface if specified
        if (surfaceWriterPtr_.valid() && surfaceWriterPtr_->enabled())
            Field<Type> allValues(values);
            combineFields(allValues);

            if (Pstream::master())
            {
                    (
                        outputDir()
                      / regionTypeNames_[regionType_] + ("_" + regionName_)
                    ),
                    false  // serial - already merged

                surfaceWriterPtr_->write(fieldName, allValues);

                surfaceWriterPtr_->clear();
            // Apply scale factor
            values *= scaleFactor_;
            Type result = processValues(values, Sf, weightField);
            switch (postOperation_)
            {
                case postOpNone:
                case postOpSqrt:
                {
                    // sqrt: component-wise - doesn't change the type
                    for (direction d=0; d < pTraits<Type>::nComponents; ++d)
                        setComponent(result, d)
                            = sqrt(mag(component(result, d)));
            // Write state/results information
            word prefix, suffix;
            {
                if (postOperation_ != postOpNone)
                {
                    // Adjust result name to include post-operation
                    prefix += postOperationTypeNames_[postOperation_];
                    prefix += '(';
                    suffix += ')';
                }

                prefix += operationTypeNames_[operation_];
                prefix += '(';
                suffix += ')';
            }

            word resultName = prefix + regionName_ + ',' + fieldName + suffix;
            // Write state/results information

            Log << "    " << prefix << regionName_ << suffix
                << " of " << fieldName << " = ";

            // Operation tagged that it always returns scalar?
            const bool alwaysScalar(operation_ & typeScalar);

            if (alwaysScalar)
            {
                const scalar sresult = component(result, 0);

                file()<< tab << sresult;

                Log << sresult << endl;

                this->setResult(resultName, sresult);
            }
            else
            {
                file()<< tab << result;

                Log << result << endl;

                this->setResult(resultName, result);
            }
Foam::functionObjects::fieldValues::surfaceFieldValue::filterField
    const GeometricField<Type, fvPatchField, volMesh>& field
    auto tvalues = tmp<Field<Type>>::New(faceId_.size());
    auto& values = tvalues.ref();
        const label facei = faceId_[i];
        const label patchi = facePatchId_[i];
            values[i] = field.boundaryField()[patchi][facei];
            FatalErrorInFunction
                << type() << " " << name() << ": "
                << regionTypeNames_[regionType_] << "(" << regionName_ << "):"
                << nl
                << "    Unable to process internal faces for volume field "
                << field.name() << nl << abort(FatalError);
        }
    // No need to flip values - all boundary faces point outwards
Foam::functionObjects::fieldValues::surfaceFieldValue::filterField
    const GeometricField<Type, fvsPatchField, surfaceMesh>& field
    auto tvalues = tmp<Field<Type>>::New(faceId_.size());
    auto& values = tvalues.ref();
        const label facei = faceId_[i];
        const label patchi = facePatchId_[i];
            values[i] = field.boundaryField()[patchi][facei];
            values[i] = field[facei];
    if (debug)
    {
        Pout<< "field " << field.name() << " oriented: "
            << field.oriented()() << endl;
    }

    if (field.oriented()())
            if (faceFlip_[i])
            {
                values[i] *= -1;
            }
    }

    return tvalues;
}


// ************************************************************************* //