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 "displacementLaplacianFvMotionSolver.H"
00027 #include <fvMotionSolvers/motionDiffusivity.H>
00028 #include <finiteVolume/fvmLaplacian.H>
00029 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00030 #include <OpenFOAM/OFstream.H>
00031 #include <meshTools/meshTools.H>
00032 #include <OpenFOAM/mapPolyMesh.H>
00033 #include <finiteVolume/volPointInterpolation.H>
00034
00035
00036
00037 namespace Foam
00038 {
00039 defineTypeNameAndDebug(displacementLaplacianFvMotionSolver, 0);
00040
00041 addToRunTimeSelectionTable
00042 (
00043 fvMotionSolver,
00044 displacementLaplacianFvMotionSolver,
00045 dictionary
00046 );
00047 }
00048
00049
00050
00051
00052 Foam::displacementLaplacianFvMotionSolver::displacementLaplacianFvMotionSolver
00053 (
00054 const polyMesh& mesh,
00055 Istream& is
00056 )
00057 :
00058 displacementFvMotionSolver(mesh, is),
00059 pointDisplacement_
00060 (
00061 IOobject
00062 (
00063 "pointDisplacement",
00064 fvMesh_.time().timeName(),
00065 fvMesh_,
00066 IOobject::MUST_READ,
00067 IOobject::AUTO_WRITE
00068 ),
00069 pointMesh::New(fvMesh_)
00070 ),
00071 cellDisplacement_
00072 (
00073 IOobject
00074 (
00075 "cellDisplacement",
00076 mesh.time().timeName(),
00077 mesh,
00078 IOobject::READ_IF_PRESENT,
00079 IOobject::AUTO_WRITE
00080 ),
00081 fvMesh_,
00082 dimensionedVector
00083 (
00084 "cellDisplacement",
00085 pointDisplacement_.dimensions(),
00086 vector::zero
00087 ),
00088 cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
00089 ),
00090 pointLocation_(NULL),
00091 diffusivityPtr_
00092 (
00093 motionDiffusivity::New(*this, lookup("diffusivity"))
00094 ),
00095 frozenPointsZone_
00096 (
00097 found("frozenPointsZone")
00098 ? fvMesh_.pointZones().findZoneID(lookup("frozenPointsZone"))
00099 : -1
00100 )
00101 {
00102 IOobject io
00103 (
00104 "pointLocation",
00105 fvMesh_.time().timeName(),
00106 fvMesh_,
00107 IOobject::MUST_READ,
00108 IOobject::AUTO_WRITE
00109 );
00110
00111 if (debug)
00112 {
00113 Info<< "displacementLaplacianFvMotionSolver:" << nl
00114 << " diffusivity : " << diffusivityPtr_().type() << nl
00115 << " frozenPoints zone : " << frozenPointsZone_ << endl;
00116 }
00117
00118
00119 if (io.headerOk())
00120 {
00121 pointLocation_.reset
00122 (
00123 new pointVectorField
00124 (
00125 io,
00126 pointMesh::New(fvMesh_)
00127 )
00128 );
00129
00130 if (debug)
00131 {
00132 Info<< "displacementLaplacianFvMotionSolver :"
00133 << " Read pointVectorField "
00134 << io.name() << " to be used for boundary conditions on points."
00135 << nl
00136 << "Boundary conditions:"
00137 << pointLocation_().boundaryField().types() << endl;
00138 }
00139 }
00140 }
00141
00142
00143
00144
00145 Foam::displacementLaplacianFvMotionSolver::
00146 ~displacementLaplacianFvMotionSolver()
00147 {}
00148
00149
00150
00151
00152 Foam::tmp<Foam::pointField>
00153 Foam::displacementLaplacianFvMotionSolver::curPoints() const
00154 {
00155 volPointInterpolation::New(fvMesh_).interpolate
00156 (
00157 cellDisplacement_,
00158 pointDisplacement_
00159 );
00160
00161 if (pointLocation_.valid())
00162 {
00163 if (debug)
00164 {
00165 Info<< "displacementLaplacianFvMotionSolver : applying "
00166 << " boundary conditions on " << pointLocation_().name()
00167 << " to new point location."
00168 << endl;
00169 }
00170
00171 pointLocation_().internalField() =
00172 points0()
00173 + pointDisplacement_.internalField();
00174
00175 pointLocation_().correctBoundaryConditions();
00176
00177
00178 if (frozenPointsZone_ != -1)
00179 {
00180 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
00181
00182 forAll(pz, i)
00183 {
00184 pointLocation_()[pz[i]] = points0()[pz[i]];
00185 }
00186 }
00187
00188 twoDCorrectPoints(pointLocation_().internalField());
00189
00190 return tmp<pointField>(pointLocation_().internalField());
00191 }
00192 else
00193 {
00194 tmp<pointField> tcurPoints
00195 (
00196 points0() + pointDisplacement_.internalField()
00197 );
00198
00199
00200 if (frozenPointsZone_ != -1)
00201 {
00202 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
00203
00204 forAll(pz, i)
00205 {
00206 tcurPoints()[pz[i]] = points0()[pz[i]];
00207 }
00208 }
00209
00210 twoDCorrectPoints(tcurPoints());
00211
00212 return tcurPoints;
00213 }
00214 }
00215
00216
00217 void Foam::displacementLaplacianFvMotionSolver::solve()
00218 {
00219
00220
00221 movePoints(fvMesh_.points());
00222
00223 diffusivityPtr_->correct();
00224 pointDisplacement_.boundaryField().updateCoeffs();
00225
00226 Foam::solve
00227 (
00228 fvm::laplacian
00229 (
00230 diffusivityPtr_->operator()(),
00231 cellDisplacement_,
00232 "laplacian(diffusivity,cellDisplacement)"
00233 )
00234 );
00235 }
00236
00237
00238 void Foam::displacementLaplacianFvMotionSolver::updateMesh
00239 (
00240 const mapPolyMesh& mpm
00241 )
00242 {
00243 displacementFvMotionSolver::updateMesh(mpm);
00244
00245
00246
00247 diffusivityPtr_.reset(NULL);
00248 diffusivityPtr_ = motionDiffusivity::New(*this, lookup("diffusivity"));
00249 }
00250
00251
00252