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 "velocityLaplacianFvMotionSolver.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(velocityLaplacianFvMotionSolver, 0);
00037
00038 addToRunTimeSelectionTable
00039 (
00040 fvMotionSolver,
00041 velocityLaplacianFvMotionSolver,
00042 dictionary
00043 );
00044 }
00045
00046
00047
00048
00049 Foam::velocityLaplacianFvMotionSolver::velocityLaplacianFvMotionSolver
00050 (
00051 const polyMesh& mesh,
00052 Istream&
00053 )
00054 :
00055 fvMotionSolver(mesh),
00056 pointMotionU_
00057 (
00058 IOobject
00059 (
00060 "pointMotionU",
00061 fvMesh_.time().timeName(),
00062 fvMesh_,
00063 IOobject::MUST_READ,
00064 IOobject::AUTO_WRITE
00065 ),
00066 pointMesh::New(fvMesh_)
00067 ),
00068 cellMotionU_
00069 (
00070 IOobject
00071 (
00072 "cellMotionU",
00073 mesh.time().timeName(),
00074 mesh,
00075 IOobject::READ_IF_PRESENT,
00076 IOobject::AUTO_WRITE
00077 ),
00078 fvMesh_,
00079 dimensionedVector
00080 (
00081 "cellMotionU",
00082 pointMotionU_.dimensions(),
00083 vector::zero
00084 ),
00085 cellMotionBoundaryTypes<vector>(pointMotionU_.boundaryField())
00086 ),
00087 diffusivityPtr_
00088 (
00089 motionDiffusivity::New(*this, lookup("diffusivity"))
00090 )
00091 {}
00092
00093
00094
00095
00096 Foam::velocityLaplacianFvMotionSolver::~velocityLaplacianFvMotionSolver()
00097 {}
00098
00099
00100
00101
00102 Foam::tmp<Foam::pointField>
00103 Foam::velocityLaplacianFvMotionSolver::curPoints() const
00104 {
00105 volPointInterpolation::New(fvMesh_).interpolate
00106 (
00107 cellMotionU_,
00108 pointMotionU_
00109 );
00110
00111 tmp<pointField> tcurPoints
00112 (
00113 fvMesh_.points()
00114 + fvMesh_.time().deltaT().value()*pointMotionU_.internalField()
00115 );
00116
00117 twoDCorrectPoints(tcurPoints());
00118
00119 return tcurPoints;
00120 }
00121
00122
00123 void Foam::velocityLaplacianFvMotionSolver::solve()
00124 {
00125
00126
00127 movePoints(fvMesh_.points());
00128
00129 diffusivityPtr_->correct();
00130 pointMotionU_.boundaryField().updateCoeffs();
00131
00132 Foam::solve
00133 (
00134 fvm::laplacian
00135 (
00136 diffusivityPtr_->operator()(),
00137 cellMotionU_,
00138 "laplacian(diffusivity,cellMotionU)"
00139 )
00140 );
00141 }
00142
00143
00144 void Foam::velocityLaplacianFvMotionSolver::updateMesh
00145 (
00146 const mapPolyMesh& mpm
00147 )
00148 {
00149 fvMotionSolver::updateMesh(mpm);
00150
00151
00152
00153 diffusivityPtr_.reset(NULL);
00154 diffusivityPtr_ = motionDiffusivity::New(*this, lookup("diffusivity"));
00155 }
00156
00157
00158