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 "pointPatchInterpolation.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <OpenFOAM/pointFields.H>
00030 #include <finiteVolume/emptyFvPatch.H>
00031 #include <OpenFOAM/demandDrivenData.H>
00032 #include <OpenFOAM/coupledPointPatchFields.H>
00033 #include <OpenFOAM/pointConstraint.H>
00034
00035
00036
00037 namespace Foam
00038 {
00039
00040
00041
00042 defineTypeNameAndDebug(pointPatchInterpolation, 0);
00043
00044
00045
00046
00047 void pointPatchInterpolation::makePatchPatchAddressing()
00048 {
00049 if (debug)
00050 {
00051 Info<< "pointPatchInterpolation::makePatchPatchAddressing() : "
00052 << "constructing boundary addressing"
00053 << endl;
00054 }
00055
00056 const fvBoundaryMesh& bm = fvMesh_.boundary();
00057 const pointBoundaryMesh& pbm = pointMesh::New(fvMesh_).boundary();
00058
00059
00060
00061 label nPatchPatchPoints = 0;
00062
00063 forAll(bm, patchi)
00064 {
00065 if(!isA<emptyFvPatch>(bm[patchi]) && !bm[patchi].coupled())
00066 {
00067 nPatchPatchPoints += bm[patchi].patch().boundaryPoints().size();
00068 }
00069 }
00070
00071
00072
00073 Map<label> patchPatchPointSet(2*nPatchPatchPoints);
00074
00075 patchPatchPoints_.setSize(nPatchPatchPoints);
00076
00077 List<pointConstraint> patchPatchPointConstraints(nPatchPatchPoints);
00078
00079 label pppi = 0;
00080
00081 forAll(bm, patchi)
00082 {
00083 if(!isA<emptyFvPatch>(bm[patchi]) && !bm[patchi].coupled())
00084 {
00085 const labelList& bp = bm[patchi].patch().boundaryPoints();
00086 const labelList& meshPoints = bm[patchi].patch().meshPoints();
00087
00088 forAll(bp, pointi)
00089 {
00090 label ppp = meshPoints[bp[pointi]];
00091
00092 Map<label>::iterator iter = patchPatchPointSet.find(ppp);
00093
00094 if (iter == patchPatchPointSet.end())
00095 {
00096 patchPatchPointSet.insert(ppp, pppi);
00097 patchPatchPoints_[pppi] = ppp;
00098
00099 pbm[patchi].applyConstraint
00100 (
00101 bp[pointi],
00102 patchPatchPointConstraints[pppi]
00103 );
00104 pppi++;
00105 }
00106 else
00107 {
00108 pbm[patchi].applyConstraint
00109 (
00110 bp[pointi],
00111 patchPatchPointConstraints[iter()]
00112 );
00113 }
00114 }
00115 }
00116 }
00117
00118 nPatchPatchPoints = pppi;
00119 patchPatchPoints_.setSize(nPatchPatchPoints);
00120 patchPatchPointConstraints.setSize(nPatchPatchPoints);
00121
00122 patchPatchPointConstraintPoints_.setSize(nPatchPatchPoints);
00123 patchPatchPointConstraintTensors_.setSize(nPatchPatchPoints);
00124
00125 label nConstraints = 0;
00126
00127 forAll(patchPatchPointConstraints, i)
00128 {
00129 if (patchPatchPointConstraints[i].first() != 0)
00130 {
00131 patchPatchPointConstraintPoints_[nConstraints] =
00132 patchPatchPoints_[i];
00133
00134 patchPatchPointConstraintTensors_[nConstraints] =
00135 patchPatchPointConstraints[i].constraintTransformation();
00136
00137 nConstraints++;
00138 }
00139 }
00140
00141 patchPatchPointConstraintPoints_.setSize(nConstraints);
00142 patchPatchPointConstraintTensors_.setSize(nConstraints);
00143
00144
00145 patchInterpolators_.clear();
00146 patchInterpolators_.setSize(bm.size());
00147
00148 forAll(bm, patchi)
00149 {
00150 patchInterpolators_.set
00151 (
00152 patchi,
00153 new primitivePatchInterpolation(bm[patchi].patch())
00154 );
00155 }
00156
00157 if (debug)
00158 {
00159 Info<< "pointPatchInterpolation::makePatchPatchAddressing() : "
00160 << "finished constructing boundary addressing"
00161 << endl;
00162 }
00163 }
00164
00165
00166 void pointPatchInterpolation::makePatchPatchWeights()
00167 {
00168 if (debug)
00169 {
00170 Info<< "pointPatchInterpolation::makePatchPatchWeights() : "
00171 << "constructing boundary weighting factors"
00172 << endl;
00173 }
00174
00175 patchPatchPointWeights_.clear();
00176 patchPatchPointWeights_.setSize(patchPatchPoints_.size());
00177
00178 const labelListList& pf = fvMesh_.pointFaces();
00179 const volVectorField& centres = fvMesh_.C();
00180 const fvBoundaryMesh& bm = fvMesh_.boundary();
00181
00182 pointScalarField sumWeights
00183 (
00184 IOobject
00185 (
00186 "sumWeights",
00187 fvMesh_.polyMesh::instance(),
00188 fvMesh_
00189 ),
00190 pointMesh::New(fvMesh_),
00191 dimensionedScalar("zero", dimless, 0)
00192 );
00193
00194 forAll(patchPatchPoints_, pointi)
00195 {
00196 const label curPoint = patchPatchPoints_[pointi];
00197 const labelList& curFaces = pf[curPoint];
00198
00199 patchPatchPointWeights_[pointi].setSize(curFaces.size());
00200 scalarList& pw = patchPatchPointWeights_[pointi];
00201
00202 label nFacesAroundPoint = 0;
00203
00204 const vector& pointLoc = fvMesh_.points()[curPoint];
00205
00206 forAll(curFaces, facei)
00207 {
00208 if (!fvMesh_.isInternalFace(curFaces[facei]))
00209 {
00210 label patchi =
00211 fvMesh_.boundaryMesh().whichPatch(curFaces[facei]);
00212
00213 if (!isA<emptyFvPatch>(bm[patchi]) && !bm[patchi].coupled())
00214 {
00215 vector d =
00216 pointLoc
00217 - centres.boundaryField()[patchi]
00218 [bm[patchi].patch().whichFace(curFaces[facei])];
00219
00220 pw[nFacesAroundPoint] = 1.0/(mag(d)+VSMALL);
00221
00222 nFacesAroundPoint++;
00223 }
00224 }
00225 }
00226
00227
00228 pw.setSize(nFacesAroundPoint);
00229
00230
00231 sumWeights[curPoint] += sum(pw);
00232 }
00233
00234
00235
00236
00237 forAll(sumWeights.boundaryField(), patchi)
00238 {
00239 if (sumWeights.boundaryField()[patchi].coupled())
00240 {
00241 refCast<coupledPointPatchScalarField>
00242 (sumWeights.boundaryField()[patchi]).initSwapAdd
00243 (
00244 sumWeights.internalField()
00245 );
00246 }
00247 }
00248
00249 forAll(sumWeights.boundaryField(), patchi)
00250 {
00251 if (sumWeights.boundaryField()[patchi].coupled())
00252 {
00253 refCast<coupledPointPatchScalarField>
00254 (sumWeights.boundaryField()[patchi]).swapAdd
00255 (
00256 sumWeights.internalField()
00257 );
00258 }
00259 }
00260
00261
00262
00263 forAll(patchPatchPoints_, pointi)
00264 {
00265 scalarList& pw = patchPatchPointWeights_[pointi];
00266 scalar sumw = sumWeights[patchPatchPoints_[pointi]];
00267
00268 forAll(pw, facei)
00269 {
00270 pw[facei] /= sumw;
00271 }
00272 }
00273
00274
00275 if (debug)
00276 {
00277 Info<< "pointPatchInterpolation::makePatchPatchWeights() : "
00278 << "finished constructing boundary weighting factors"
00279 << endl;
00280 }
00281 }
00282
00283
00284
00285
00286 pointPatchInterpolation::pointPatchInterpolation(const fvMesh& vm)
00287 :
00288 fvMesh_(vm)
00289 {
00290 updateMesh();
00291 }
00292
00293
00294
00295
00296 pointPatchInterpolation::~pointPatchInterpolation()
00297 {}
00298
00299
00300
00301
00302 void pointPatchInterpolation::updateMesh()
00303 {
00304 makePatchPatchAddressing();
00305 makePatchPatchWeights();
00306 }
00307
00308
00309 bool pointPatchInterpolation::movePoints()
00310 {
00311 forAll(patchInterpolators_, patchi)
00312 {
00313 patchInterpolators_[patchi].movePoints();
00314 }
00315
00316 makePatchPatchWeights();
00317
00318 return true;
00319 }
00320
00321
00322
00323
00324 template<>
00325 void pointPatchInterpolation::applyCornerConstraints<scalar>
00326 (
00327 GeometricField<scalar, pointPatchField, pointMesh>& pf
00328 ) const
00329 {}
00330
00331
00332
00333
00334 }
00335
00336