FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

pointPatchInterpolate.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "pointPatchInterpolation.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <OpenFOAM/pointFields.H>
00029 #include <finiteVolume/emptyFvPatch.H>
00030 #include <OpenFOAM/valuePointPatchField.H>
00031 #include <OpenFOAM/coupledPointPatchField.H>
00032 #include <OpenFOAM/coupledFacePointPatch.H>
00033 #include <OpenFOAM/transform.H>
00034 
00035 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00036 
00037 namespace Foam
00038 {
00039 
00040 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00041 
00042 template<class Type>
00043 void pointPatchInterpolation::interpolate
00044 (
00045     const GeometricField<Type, fvPatchField, volMesh>& vf,
00046     GeometricField<Type, pointPatchField, pointMesh>& pf,
00047     bool overrideFixedValue
00048 ) const
00049 {
00050     if (debug)
00051     {
00052         Info<< "pointPatchInterpolation::interpolate("
00053             << "const GeometricField<Type, fvPatchField, volMesh>&, "
00054             << "GeometricField<Type, pointPatchField, pointMesh>&) : "
00055             << "interpolating field from cells to points"
00056             << endl;
00057     }
00058 
00059     // Interpolate patch values: over-ride the internal values for the points
00060     // on the patch with the interpolated point values from the faces of the
00061     // patch
00062 
00063     const fvBoundaryMesh& bm = fvMesh_.boundary();
00064     const pointBoundaryMesh& pbm = pointMesh::New(fvMesh_).boundary();
00065 
00066     forAll(bm, patchi)
00067     {
00068         if (!isA<emptyFvPatch>(bm[patchi]) && !bm[patchi].coupled())
00069         {
00070             pointPatchField<Type>& ppf = pf.boundaryField()[patchi];
00071 
00072             // Only map the values corresponding to the points associated with
00073             // faces, not "lone" points due to decomposition
00074             ppf.setInInternalField
00075             (
00076                 pf.internalField(),
00077                 patchInterpolators_[patchi]
00078                .faceToPointInterpolate(vf.boundaryField()[patchi])(),
00079                 bm[patchi].patch().meshPoints()
00080             );
00081 
00082             if
00083             (
00084                 overrideFixedValue
00085              && isA<valuePointPatchField<Type> >(ppf)
00086             )
00087             {
00088                 refCast<valuePointPatchField<Type> >(ppf) = ppf;
00089             }
00090         }
00091         else if (bm[patchi].coupled())
00092         {
00093             // Initialise the "lone" points on the coupled patch to zero,
00094             // these values are obtained from the couple-transfer
00095 
00096             const labelList& loneMeshPoints =
00097                 refCast<const coupledFacePointPatch>(pbm[patchi])
00098                .loneMeshPoints();
00099 
00100             forAll(loneMeshPoints, i)
00101             {
00102                 pf[loneMeshPoints[i]] = pTraits<Type>::zero;
00103             }
00104         }
00105 
00106     }
00107 
00108 
00109     // Correct patch-patch boundary points by interpolation "around" corners
00110     const labelListList& PointFaces = fvMesh_.pointFaces();
00111 
00112     forAll(patchPatchPoints_, pointi)
00113     {
00114         const label curPoint = patchPatchPoints_[pointi];
00115         const labelList& curFaces = PointFaces[curPoint];
00116 
00117         label fI = 0;
00118 
00119         // Reset the boundary value before accumulation
00120         pf[curPoint] = pTraits<Type>::zero;
00121 
00122         // Go through all the faces
00123         forAll(curFaces, facei)
00124         {
00125             if (!fvMesh_.isInternalFace(curFaces[facei]))
00126             {
00127                 label patchi =
00128                     fvMesh_.boundaryMesh().whichPatch(curFaces[facei]);
00129 
00130                 if (!isA<emptyFvPatch>(bm[patchi]) && !bm[patchi].coupled())
00131                 {
00132                     label faceInPatchi =
00133                         bm[patchi].patch().whichFace(curFaces[facei]);
00134 
00135                     pf[curPoint] +=
00136                         patchPatchPointWeights_[pointi][fI]
00137                        *vf.boundaryField()[patchi][faceInPatchi];
00138 
00139                     fI++;
00140                 }
00141             }
00142         }
00143     }
00144 
00145     // Update coupled boundaries
00146     forAll(pf.boundaryField(), patchi)
00147     {
00148         if (pf.boundaryField()[patchi].coupled())
00149         {
00150             refCast<coupledPointPatchField<Type> >(pf.boundaryField()[patchi])
00151                 .initSwapAdd(pf.internalField());
00152         }
00153     }
00154 
00155     forAll(pf.boundaryField(), patchi)
00156     {
00157         if (pf.boundaryField()[patchi].coupled())
00158         {
00159             refCast<coupledPointPatchField<Type> >(pf.boundaryField()[patchi])
00160                 .swapAdd(pf.internalField());
00161         }
00162     }
00163 
00164 
00165     // Override constrained pointPatchField types with the constraint value.
00166     // This relys on only constrained pointPatchField implementing the evaluate
00167     // function
00168     pf.correctBoundaryConditions();
00169 
00170 
00171     // Apply multiple constraints on edge/corner points
00172     applyCornerConstraints(pf);
00173 
00174 
00175     if (debug)
00176     {
00177         Info<< "pointPatchInterpolation::interpolate("
00178             << "const GeometricField<Type, fvPatchField, volMesh>&, "
00179             << "GeometricField<Type, pointPatchField, pointMesh>&) : "
00180             << "finished interpolating field from cells to points"
00181             << endl;
00182     }
00183 }
00184 
00185 
00186 template<class Type>
00187 void pointPatchInterpolation::applyCornerConstraints
00188 (
00189     GeometricField<Type, pointPatchField, pointMesh>& pf
00190 ) const
00191 {
00192     forAll(patchPatchPointConstraintPoints_, pointi)
00193     {
00194         pf[patchPatchPointConstraintPoints_[pointi]] = transform
00195         (
00196             patchPatchPointConstraintTensors_[pointi],
00197             pf[patchPatchPointConstraintPoints_[pointi]]
00198         );
00199     }
00200 }
00201 
00202 
00203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00204 
00205 } // End namespace Foam
00206 
00207 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines