Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "volPointInterpolation.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <OpenFOAM/pointFields.H>
00030 #include <OpenFOAM/demandDrivenData.H>
00031 #include <OpenFOAM/coupledPointPatchFields.H>
00032 #include <OpenFOAM/pointConstraint.H>
00033
00034
00035
00036 namespace Foam
00037 {
00038
00039
00040
00041 defineTypeNameAndDebug(volPointInterpolation, 0);
00042
00043
00044
00045
00046 void volPointInterpolation::makeWeights()
00047 {
00048 if (debug)
00049 {
00050 Info<< "volPointInterpolation::makeWeights() : "
00051 << "constructing weighting factors"
00052 << endl;
00053 }
00054
00055 const pointField& points = mesh().points();
00056 const labelListList& pointCells = mesh().pointCells();
00057 const vectorField& cellCentres = mesh().cellCentres();
00058
00059
00060 pointWeights_.clear();
00061 pointWeights_.setSize(points.size());
00062
00063 forAll(pointWeights_, pointi)
00064 {
00065 pointWeights_[pointi].setSize(pointCells[pointi].size());
00066 }
00067
00068 pointScalarField sumWeights
00069 (
00070 IOobject
00071 (
00072 "volPointSumWeights",
00073 mesh().polyMesh::instance(),
00074 mesh()
00075 ),
00076 pointMesh::New(mesh()),
00077 dimensionedScalar("zero", dimless, 0)
00078 );
00079
00080
00081
00082 forAll(points, pointi)
00083 {
00084 scalarList& pw = pointWeights_[pointi];
00085 const labelList& pcp = pointCells[pointi];
00086
00087 forAll(pcp, pointCelli)
00088 {
00089 pw[pointCelli] =
00090 1.0/mag(points[pointi] - cellCentres[pcp[pointCelli]]);
00091
00092 sumWeights[pointi] += pw[pointCelli];
00093 }
00094 }
00095
00096 forAll(sumWeights.boundaryField(), patchi)
00097 {
00098 if (sumWeights.boundaryField()[patchi].coupled())
00099 {
00100 refCast<coupledPointPatchScalarField>
00101 (sumWeights.boundaryField()[patchi]).initSwapAdd
00102 (
00103 sumWeights.internalField()
00104 );
00105 }
00106 }
00107
00108 forAll(sumWeights.boundaryField(), patchi)
00109 {
00110 if (sumWeights.boundaryField()[patchi].coupled())
00111 {
00112 refCast<coupledPointPatchScalarField>
00113 (sumWeights.boundaryField()[patchi]).swapAdd
00114 (
00115 sumWeights.internalField()
00116 );
00117 }
00118 }
00119
00120 forAll(points, pointi)
00121 {
00122 scalarList& pw = pointWeights_[pointi];
00123
00124 forAll(pw, pointCelli)
00125 {
00126 pw[pointCelli] /= sumWeights[pointi];
00127 }
00128 }
00129
00130 if (debug)
00131 {
00132 Info<< "volPointInterpolation::makeWeights() : "
00133 << "finished constructing weighting factors"
00134 << endl;
00135 }
00136 }
00137
00138
00139
00140
00141 volPointInterpolation::volPointInterpolation(const fvMesh& vm)
00142 :
00143 MeshObject<fvMesh, volPointInterpolation>(vm),
00144 boundaryInterpolator_(vm)
00145 {
00146 updateMesh();
00147 }
00148
00149
00150
00151
00152 volPointInterpolation::~volPointInterpolation()
00153 {}
00154
00155
00156
00157
00158 void volPointInterpolation::updateMesh()
00159 {
00160 makeWeights();
00161 boundaryInterpolator_.updateMesh();
00162 }
00163
00164
00165 bool volPointInterpolation::movePoints()
00166 {
00167 makeWeights();
00168 boundaryInterpolator_.movePoints();
00169
00170 return true;
00171 }
00172
00173
00174
00175
00176 }
00177
00178