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 "velocityComponentLaplacianFvMotionSolver.H"
00027 #include <fvMotionSolvers/motionDiffusivity.H>
00028 #include <finiteVolume/fvmLaplacian.H>
00029 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00030 #include <finiteVolume/volPointInterpolation.H>
00031
00032
00033
00034 namespace Foam
00035 {
00036 defineTypeNameAndDebug(velocityComponentLaplacianFvMotionSolver, 0);
00037
00038 addToRunTimeSelectionTable
00039 (
00040 fvMotionSolver,
00041 velocityComponentLaplacianFvMotionSolver,
00042 dictionary
00043 );
00044 }
00045
00046
00047
00048
00049 Foam::velocityComponentLaplacianFvMotionSolver::
00050 velocityComponentLaplacianFvMotionSolver
00051 (
00052 const polyMesh& mesh,
00053 Istream& msData
00054 )
00055 :
00056 fvMotionSolver(mesh),
00057 cmptName_(msData),
00058 cmpt_(0),
00059 pointMotionU_
00060 (
00061 IOobject
00062 (
00063 "pointMotionU" + cmptName_,
00064 fvMesh_.time().timeName(),
00065 fvMesh_,
00066 IOobject::MUST_READ,
00067 IOobject::AUTO_WRITE
00068 ),
00069 pointMesh::New(fvMesh_)
00070 ),
00071 cellMotionU_
00072 (
00073 IOobject
00074 (
00075 "cellMotionU" + cmptName_,
00076 mesh.time().timeName(),
00077 mesh,
00078 IOobject::READ_IF_PRESENT,
00079 IOobject::AUTO_WRITE
00080 ),
00081 fvMesh_,
00082 dimensionedScalar
00083 (
00084 "cellMotionU",
00085 pointMotionU_.dimensions(),
00086 0
00087 ),
00088 cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
00089 ),
00090 diffusivityPtr_
00091 (
00092 motionDiffusivity::New(*this, lookup("diffusivity"))
00093 )
00094 {
00095 if (cmptName_ == "x")
00096 {
00097 cmpt_ = vector::X;
00098 }
00099 else if (cmptName_ == "y")
00100 {
00101 cmpt_ = vector::Y;
00102 }
00103 else if (cmptName_ == "z")
00104 {
00105 cmpt_ = vector::Z;
00106 }
00107 else
00108 {
00109 FatalErrorIn
00110 (
00111 "velocityComponentLaplacianFvMotionSolver::"
00112 "velocityComponentLaplacianFvMotionSolver"
00113 "(const polyMesh& mesh, Istream& msData)"
00114 ) << "Given component name " << cmptName_ << " should be x, y or z"
00115 << exit(FatalError);
00116 }
00117 }
00118
00119
00120
00121
00122 Foam::velocityComponentLaplacianFvMotionSolver::
00123 ~velocityComponentLaplacianFvMotionSolver()
00124 {}
00125
00126
00127
00128
00129 Foam::tmp<Foam::pointField>
00130 Foam::velocityComponentLaplacianFvMotionSolver::curPoints() const
00131 {
00132 volPointInterpolation::New(fvMesh_).interpolate
00133 (
00134 cellMotionU_,
00135 pointMotionU_
00136 );
00137
00138 tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
00139
00140 tcurPoints().replace
00141 (
00142 cmpt_,
00143 tcurPoints().component(cmpt_)
00144 + fvMesh_.time().deltaT().value()*pointMotionU_.internalField()
00145 );
00146
00147 twoDCorrectPoints(tcurPoints());
00148
00149 return tcurPoints;
00150 }
00151
00152
00153 void Foam::velocityComponentLaplacianFvMotionSolver::solve()
00154 {
00155
00156
00157 movePoints(fvMesh_.points());
00158
00159 diffusivityPtr_->correct();
00160 pointMotionU_.boundaryField().updateCoeffs();
00161
00162 Foam::solve
00163 (
00164 fvm::laplacian
00165 (
00166 diffusivityPtr_->operator()(),
00167 cellMotionU_,
00168 "laplacian(diffusivity,cellMotionU)"
00169 )
00170 );
00171 }
00172
00173
00174 void Foam::velocityComponentLaplacianFvMotionSolver::updateMesh
00175 (
00176 const mapPolyMesh& mpm
00177 )
00178 {
00179 fvMotionSolver::updateMesh(mpm);
00180
00181
00182
00183 diffusivityPtr_.reset(NULL);
00184 diffusivityPtr_ = motionDiffusivity::New(*this, lookup("diffusivity"));
00185 }
00186
00187
00188