/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd. ------------------------------------------------------------------------------- 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 . Class Foam::volPointInterpolation Description Interpolate from cell centres to points (vertices) using inverse distance weighting SourceFiles volPointInterpolation.C volPointInterpolate.C \*---------------------------------------------------------------------------*/ #ifndef volPointInterpolation_H #define volPointInterpolation_H #include "MeshObject.H" #include "scalarList.H" #include "volFields.H" #include "pointFields.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { class fvMesh; class pointMesh; /*---------------------------------------------------------------------------*\ Class volPointInterpolation Declaration \*---------------------------------------------------------------------------*/ class volPointInterpolation : public MeshObject { // Private data //- Interpolation scheme weighting factor array. scalarListList pointWeights_; // Boundary handling //- Boundary addressing autoPtr boundaryPtr_; //- Per boundary face whether is on non-coupled patch boolList boundaryIsPatchFace_; //- Per mesh(!) point whether is on non-coupled patch (on any // processor) boolList isPatchPoint_; //- Per boundary point the weights per pointFaces. scalarListList boundaryPointWeights_; // Private Member Functions //- Construct addressing over all boundary faces void calcBoundaryAddressing(); //- Make weights for internal and coupled-only boundarypoints void makeInternalWeights(scalarField& sumWeights); //- Make weights for points on uncoupled patches void makeBoundaryWeights(scalarField& sumWeights); //- Construct all point weighting factors void makeWeights(); //- Helper: push master point data to collocated points template void pushUntransformedData(List&) const; //- Get boundary field in same order as boundary faces. Field is // zero on all coupled and empty patches template tmp > flatBoundaryField ( const GeometricField& vf ) const; //- Add separated contributions template void addSeparated ( GeometricField& ) const; //- Disallow default bitwise copy construct volPointInterpolation(const volPointInterpolation&); //- Disallow default bitwise assignment void operator=(const volPointInterpolation&); public: // Declare name of the class and its debug switch ClassName("volPointInterpolation"); // Constructors //- Constructor given fvMesh and pointMesh. explicit volPointInterpolation(const fvMesh&); //- Destructor ~volPointInterpolation(); // Member functions // Edit //- Update mesh topology using the morph engine void updateMesh(const mapPolyMesh&); //- Correct weighting factors for moving mesh. bool movePoints(); // Interpolation functions //- Interpolate volField using inverse distance weighting // returning pointField template tmp > interpolate ( const GeometricField& ) const; //- Interpolate tmp using inverse distance weighting // returning pointField template tmp > interpolate ( const tmp >& ) const; //- Interpolate volField using inverse distance weighting // returning pointField with the same patchField types. Assigns // to any fixedValue boundary conditions to make them consistent // with internal field template tmp > interpolate ( const GeometricField&, const wordList& patchFieldTypes ) const; //- Interpolate tmp using inverse distance weighting // returning pointField with the same patchField types. Assigns // to any fixedValue boundary conditions to make them consistent // with internal field template tmp > interpolate ( const tmp >&, const wordList& patchFieldTypes ) const; //- Interpolate dimensionedField using inverse distance weighting // returning pointField template tmp > interpolate ( const DimensionedField& ) const; //- Interpolate tmp using inverse distance // weighting returning pointField template tmp > interpolate ( const tmp >& ) const; // Low level //- Interpolate internal field from volField to pointField // using inverse distance weighting template void interpolateInternalField ( const GeometricField&, GeometricField& ) const; //- Interpolate boundary field without applying constraints/boundary // conditions template void interpolateBoundaryField ( const GeometricField& vf, GeometricField& pf ) const; //- Interpolate boundary with constraints/boundary conditions template void interpolateBoundaryField ( const GeometricField& vf, GeometricField& pf, const bool overrideFixedValue ) const; //- Interpolate from volField to pointField // using inverse distance weighting template void interpolate ( const GeometricField&, GeometricField& ) const; //- Interpolate volField using inverse distance weighting // returning pointField with name. Optionally caches template tmp > interpolate ( const GeometricField&, const word& name, const bool cache ) const; //- Interpolate dimensioned internal field from cells to points // using inverse distance weighting template void interpolateDimensionedInternalField ( const DimensionedField& vf, DimensionedField& pf ) const; //- Interpolate dimensionedField using inverse distance weighting // returning pointField with name. Optionally caches template tmp > interpolate ( const DimensionedField&, const word& name, const bool cache ) const; // Interpolation for displacement (applies 2D correction) //- Interpolate from volField to pointField // using inverse distance weighting void interpolateDisplacement ( const volVectorField&, pointVectorField& ) const; }; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace Foam // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #ifdef NoRepository # include "volPointInterpolate.C" #endif // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #endif // ************************************************************************* //