initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / fvMotionSolver / fvMotionSolvers / displacement / laplacian / displacementLaplacianFvMotionSolver.C
blobe57a17cce9baa0d341c724ce801832b1ea0717e2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "displacementLaplacianFvMotionSolver.H"
28 #include "motionDiffusivity.H"
29 #include "fvmLaplacian.H"
30 #include "addToRunTimeSelectionTable.H"
31 #include "OFstream.H"
32 #include "meshTools.H"
33 #include "mapPolyMesh.H"
34 #include "volPointInterpolation.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 namespace Foam
40     defineTypeNameAndDebug(displacementLaplacianFvMotionSolver, 0);
42     addToRunTimeSelectionTable
43     (
44         fvMotionSolver,
45         displacementLaplacianFvMotionSolver,
46         dictionary
47     );
51 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
53 Foam::displacementLaplacianFvMotionSolver::displacementLaplacianFvMotionSolver
55     const polyMesh& mesh,
56     Istream& is
59     displacementFvMotionSolver(mesh, is),
60     pointDisplacement_
61     (
62         IOobject
63         (
64             "pointDisplacement",
65             fvMesh_.time().timeName(),
66             fvMesh_,
67             IOobject::MUST_READ,
68             IOobject::AUTO_WRITE
69         ),
70         pointMesh::New(fvMesh_)
71     ),
72     cellDisplacement_
73     (
74         IOobject
75         (
76             "cellDisplacement",
77             mesh.time().timeName(),
78             mesh,
79             IOobject::READ_IF_PRESENT,
80             IOobject::AUTO_WRITE
81         ),
82         fvMesh_,
83         dimensionedVector
84         (
85             "cellDisplacement",
86             pointDisplacement_.dimensions(),
87             vector::zero
88         ),
89         cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
90     ),
91     pointLocation_(NULL),
92     diffusivityPtr_
93     (
94         motionDiffusivity::New(*this, lookup("diffusivity"))
95     ),
96     frozenPointsZone_
97     (
98         found("frozenPointsZone")
99       ? fvMesh_.pointZones().findZoneID(lookup("frozenPointsZone"))
100       : -1
101     )
103     IOobject io
104     (
105         "pointLocation",
106         fvMesh_.time().timeName(),
107         fvMesh_,
108         IOobject::MUST_READ,
109         IOobject::AUTO_WRITE
110     );
112     if (debug)
113     {
114         Info<< "displacementLaplacianFvMotionSolver:" << nl
115             << "    diffusivity       : " << diffusivityPtr_().type() << nl
116             << "    frozenPoints zone : " << frozenPointsZone_ << endl;
117     }
120     if (io.headerOk())
121     {
122         pointLocation_.reset
123         (
124             new pointVectorField
125             (
126                 io,
127                 pointMesh::New(fvMesh_)
128             )
129         );
131         if (debug)
132         {
133             Info<< "displacementLaplacianFvMotionSolver :"
134                 << " Read pointVectorField "
135                 << io.name() << " to be used for boundary conditions on points."
136                 << nl
137                 << "Boundary conditions:"
138                 << pointLocation_().boundaryField().types() << endl;
139         }
140     }
144 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
146 Foam::displacementLaplacianFvMotionSolver::
147 ~displacementLaplacianFvMotionSolver()
151 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
153 Foam::tmp<Foam::pointField>
154 Foam::displacementLaplacianFvMotionSolver::curPoints() const
156     volPointInterpolation::New(fvMesh_).interpolate
157     (
158         cellDisplacement_,
159         pointDisplacement_
160     );
162     if (pointLocation_.valid())
163     {
164         if (debug)
165         {
166             Info<< "displacementLaplacianFvMotionSolver : applying "
167                 << " boundary conditions on " << pointLocation_().name()
168                 << " to new point location."
169                 << endl;
170         }
172         pointLocation_().internalField() =
173             points0()
174           + pointDisplacement_.internalField();
176         pointLocation_().correctBoundaryConditions();
178         // Implement frozen points
179         if (frozenPointsZone_ != -1)
180         {
181             const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
183             forAll(pz, i)
184             {
185                 pointLocation_()[pz[i]] = points0()[pz[i]];
186             }
187         }
189         twoDCorrectPoints(pointLocation_().internalField());
191         return tmp<pointField>(pointLocation_().internalField());
192     }
193     else
194     {
195         tmp<pointField> tcurPoints
196         (
197             points0() + pointDisplacement_.internalField()
198         );
200         // Implement frozen points
201         if (frozenPointsZone_ != -1)
202         {
203             const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
205             forAll(pz, i)
206             {
207                 tcurPoints()[pz[i]] = points0()[pz[i]];
208             }
209         }
211         twoDCorrectPoints(tcurPoints());
213         return tcurPoints;
214     }
218 void Foam::displacementLaplacianFvMotionSolver::solve()
220     // The points have moved so before interpolation update
221     // the fvMotionSolver accordingly
222     movePoints(fvMesh_.points());
224     diffusivityPtr_->correct();
225     pointDisplacement_.boundaryField().updateCoeffs();
227     Foam::solve
228     (
229         fvm::laplacian
230         (
231             diffusivityPtr_->operator()(),
232             cellDisplacement_,
233             "laplacian(diffusivity,cellDisplacement)"
234         )
235     );
239 void Foam::displacementLaplacianFvMotionSolver::updateMesh
241     const mapPolyMesh& mpm
244     displacementFvMotionSolver::updateMesh(mpm);
246     // Update diffusivity. Note two stage to make sure old one is de-registered
247     // before creating/registering new one.
248     diffusivityPtr_.reset(NULL);
249     diffusivityPtr_ = motionDiffusivity::New(*this, lookup("diffusivity"));
253 // ************************************************************************* //