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 "displacementComponentLaplacianFvMotionSolver.H"
00027 #include <fvMotionSolvers/motionDiffusivity.H>
00028 #include <finiteVolume/fvmLaplacian.H>
00029 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00030 #include <OpenFOAM/mapPolyMesh.H>
00031 #include <finiteVolume/volPointInterpolation.H>
00032
00033
00034
00035 namespace Foam
00036 {
00037 defineTypeNameAndDebug(displacementComponentLaplacianFvMotionSolver, 0);
00038
00039 addToRunTimeSelectionTable
00040 (
00041 fvMotionSolver,
00042 displacementComponentLaplacianFvMotionSolver,
00043 dictionary
00044 );
00045 }
00046
00047
00048
00049
00050 Foam::direction Foam::displacementComponentLaplacianFvMotionSolver::cmpt
00051 (
00052 const word& cmptName
00053 ) const
00054 {
00055 if (cmptName == "x")
00056 {
00057 return vector::X;
00058 }
00059 else if (cmptName == "y")
00060 {
00061 return vector::Y;
00062 }
00063 else if (cmptName == "z")
00064 {
00065 return vector::Z;
00066 }
00067 else
00068 {
00069 FatalErrorIn
00070 (
00071 "displacementComponentLaplacianFvMotionSolver::"
00072 "displacementComponentLaplacianFvMotionSolver"
00073 "(const polyMesh& mesh, Istream& msData)"
00074 ) << "Given component name " << cmptName << " should be x, y or z"
00075 << exit(FatalError);
00076
00077 return 0;
00078 }
00079 }
00080
00081
00082
00083
00084 Foam::displacementComponentLaplacianFvMotionSolver::
00085 displacementComponentLaplacianFvMotionSolver
00086 (
00087 const polyMesh& mesh,
00088 Istream& msData
00089 )
00090 :
00091 fvMotionSolver(mesh),
00092 cmptName_(msData),
00093 cmpt_(cmpt(cmptName_)),
00094 points0_
00095 (
00096 pointIOField
00097 (
00098 IOobject
00099 (
00100 "points",
00101 time().constant(),
00102 polyMesh::meshSubDir,
00103 mesh,
00104 IOobject::MUST_READ,
00105 IOobject::NO_WRITE,
00106 false
00107 )
00108 ).component(cmpt_)
00109 ),
00110 pointDisplacement_
00111 (
00112 IOobject
00113 (
00114 "pointDisplacement" + cmptName_,
00115 fvMesh_.time().timeName(),
00116 fvMesh_,
00117 IOobject::MUST_READ,
00118 IOobject::AUTO_WRITE
00119 ),
00120 pointMesh::New(fvMesh_)
00121 ),
00122 cellDisplacement_
00123 (
00124 IOobject
00125 (
00126 "cellDisplacement" + cmptName_,
00127 mesh.time().timeName(),
00128 mesh,
00129 IOobject::READ_IF_PRESENT,
00130 IOobject::AUTO_WRITE
00131 ),
00132 fvMesh_,
00133 dimensionedScalar
00134 (
00135 "cellDisplacement",
00136 pointDisplacement_.dimensions(),
00137 0
00138 ),
00139 cellMotionBoundaryTypes<scalar>(pointDisplacement_.boundaryField())
00140 ),
00141 pointLocation_(NULL),
00142 diffusivityPtr_
00143 (
00144 motionDiffusivity::New(*this, lookup("diffusivity"))
00145 ),
00146 frozenPointsZone_
00147 (
00148 found("frozenPointsZone")
00149 ? fvMesh_.pointZones().findZoneID(lookup("frozenPointsZone"))
00150 : -1
00151 )
00152 {
00153 IOobject io
00154 (
00155 "pointLocation",
00156 fvMesh_.time().timeName(),
00157 fvMesh_,
00158 IOobject::MUST_READ,
00159 IOobject::AUTO_WRITE
00160 );
00161
00162 if (debug)
00163 {
00164 Info<< "displacementComponentLaplacianFvMotionSolver:" << nl
00165 << " diffusivity : " << diffusivityPtr_().type() << nl
00166 << " frozenPoints zone : " << frozenPointsZone_ << endl;
00167 }
00168
00169 if (io.headerOk())
00170 {
00171 pointLocation_.reset
00172 (
00173 new pointVectorField
00174 (
00175 io,
00176 pointMesh::New(fvMesh_)
00177 )
00178 );
00179
00180 if (debug)
00181 {
00182 Info<< "displacementComponentLaplacianFvMotionSolver :"
00183 << " Read pointVectorField "
00184 << io.name() << " to be used for boundary conditions on points."
00185 << nl
00186 << "Boundary conditions:"
00187 << pointLocation_().boundaryField().types() << endl;
00188 }
00189 }
00190 }
00191
00192
00193
00194
00195 Foam::displacementComponentLaplacianFvMotionSolver::
00196 ~displacementComponentLaplacianFvMotionSolver()
00197 {}
00198
00199
00200
00201
00202 Foam::tmp<Foam::pointField>
00203 Foam::displacementComponentLaplacianFvMotionSolver::curPoints() const
00204 {
00205 volPointInterpolation::New(fvMesh_).interpolate
00206 (
00207 cellDisplacement_,
00208 pointDisplacement_
00209 );
00210
00211 if (pointLocation_.valid())
00212 {
00213 if (debug)
00214 {
00215 Info<< "displacementComponentLaplacianFvMotionSolver : applying "
00216 << " boundary conditions on " << pointLocation_().name()
00217 << " to new point location."
00218 << endl;
00219 }
00220
00221
00222
00223 pointLocation_().internalField() = fvMesh_.points();
00224
00225 pointLocation_().internalField().replace
00226 (
00227 cmpt_,
00228 points0_ + pointDisplacement_.internalField()
00229 );
00230
00231 pointLocation_().correctBoundaryConditions();
00232
00233
00234 if (frozenPointsZone_ != -1)
00235 {
00236 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
00237
00238 forAll(pz, i)
00239 {
00240 label pointI = pz[i];
00241
00242 pointLocation_()[pointI][cmpt_] = points0_[pointI];
00243 }
00244 }
00245
00246 twoDCorrectPoints(pointLocation_().internalField());
00247
00248 return tmp<pointField>(pointLocation_().internalField());
00249 }
00250 else
00251 {
00252 tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
00253
00254 tcurPoints().replace
00255 (
00256 cmpt_,
00257 points0_ + pointDisplacement_.internalField()
00258 );
00259
00260
00261 if (frozenPointsZone_ != -1)
00262 {
00263 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
00264
00265 forAll(pz, i)
00266 {
00267 label pointI = pz[i];
00268
00269 tcurPoints()[pointI][cmpt_] = points0_[pointI];
00270 }
00271 }
00272
00273 twoDCorrectPoints(tcurPoints());
00274
00275 return tcurPoints;
00276 }
00277 }
00278
00279
00280 void Foam::displacementComponentLaplacianFvMotionSolver::solve()
00281 {
00282
00283
00284 movePoints(fvMesh_.points());
00285
00286 diffusivityPtr_->correct();
00287 pointDisplacement_.boundaryField().updateCoeffs();
00288
00289 Foam::solve
00290 (
00291 fvm::laplacian
00292 (
00293 diffusivityPtr_->operator()(),
00294 cellDisplacement_,
00295 "laplacian(diffusivity,cellDisplacement)"
00296 )
00297 );
00298 }
00299
00300
00301 void Foam::displacementComponentLaplacianFvMotionSolver::updateMesh
00302 (
00303 const mapPolyMesh& mpm
00304 )
00305 {
00306 fvMotionSolver::updateMesh(mpm);
00307
00308
00309
00310
00311
00312
00313 const scalarField points
00314 (
00315 mpm.hasMotionPoints()
00316 ? mpm.preMotionPoints().component(cmpt_)
00317 : fvMesh_.points().component(cmpt_)
00318 );
00319
00320
00321 const scalar scale =
00322 (gMax(points0_)-gMin(points0_))
00323 /(gMax(points)-gMin(points));
00324
00325 scalarField newPoints0(mpm.pointMap().size());
00326
00327 forAll(newPoints0, pointI)
00328 {
00329 label oldPointI = mpm.pointMap()[pointI];
00330
00331 if (oldPointI >= 0)
00332 {
00333 label masterPointI = mpm.reversePointMap()[oldPointI];
00334
00335 if (masterPointI == pointI)
00336 {
00337 newPoints0[pointI] = points0_[oldPointI];
00338 }
00339 else
00340 {
00341
00342 newPoints0[pointI] =
00343 points0_[oldPointI]
00344 + scale*(points[pointI]-points[masterPointI]);
00345 }
00346 }
00347 else
00348 {
00349 FatalErrorIn
00350 (
00351 "displacementLaplacianFvMotionSolver::updateMesh"
00352 "(const mapPolyMesh& mpm)"
00353 ) << "Cannot work out coordinates of introduced vertices."
00354 << " New vertex " << pointI << " at coordinate "
00355 << points[pointI] << exit(FatalError);
00356 }
00357 }
00358 points0_.transfer(newPoints0);
00359
00360
00361
00362 diffusivityPtr_.reset(NULL);
00363 diffusivityPtr_ = motionDiffusivity::New(*this, lookup("diffusivity"));
00364 }
00365
00366
00367