FreeFOAM The Cross-Platform CFD Toolkit
displacementLaplacianFvMotionSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
30 #include <OpenFOAM/OFstream.H>
31 #include <meshTools/meshTools.H>
32 #include <OpenFOAM/mapPolyMesh.H>
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(displacementLaplacianFvMotionSolver, 0);
40 
42  (
43  fvMotionSolver,
44  displacementLaplacianFvMotionSolver,
45  dictionary
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 Foam::displacementLaplacianFvMotionSolver::displacementLaplacianFvMotionSolver
53 (
54  const polyMesh& mesh,
55  Istream& is
56 )
57 :
59  pointDisplacement_
60  (
61  IOobject
62  (
63  "pointDisplacement",
64  fvMesh_.time().timeName(),
65  fvMesh_,
66  IOobject::MUST_READ,
67  IOobject::AUTO_WRITE
68  ),
69  pointMesh::New(fvMesh_)
70  ),
71  cellDisplacement_
72  (
73  IOobject
74  (
75  "cellDisplacement",
76  mesh.time().timeName(),
77  mesh,
78  IOobject::READ_IF_PRESENT,
79  IOobject::AUTO_WRITE
80  ),
81  fvMesh_,
83  (
84  "cellDisplacement",
85  pointDisplacement_.dimensions(),
86  vector::zero
87  ),
88  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
89  ),
90  pointLocation_(NULL),
91  diffusivityPtr_
92  (
93  motionDiffusivity::New(*this, lookup("diffusivity"))
94  ),
95  frozenPointsZone_
96  (
97  found("frozenPointsZone")
98  ? fvMesh_.pointZones().findZoneID(lookup("frozenPointsZone"))
99  : -1
100  )
101 {
102  IOobject io
103  (
104  "pointLocation",
105  fvMesh_.time().timeName(),
106  fvMesh_,
107  IOobject::MUST_READ,
108  IOobject::AUTO_WRITE
109  );
110 
111  if (debug)
112  {
113  Info<< "displacementLaplacianFvMotionSolver:" << nl
114  << " diffusivity : " << diffusivityPtr_().type() << nl
115  << " frozenPoints zone : " << frozenPointsZone_ << endl;
116  }
117 
118 
119  if (io.headerOk())
120  {
121  pointLocation_.reset
122  (
123  new pointVectorField
124  (
125  io,
126  pointMesh::New(fvMesh_)
127  )
128  );
129 
130  if (debug)
131  {
132  Info<< "displacementLaplacianFvMotionSolver :"
133  << " Read pointVectorField "
134  << io.name() << " to be used for boundary conditions on points."
135  << nl
136  << "Boundary conditions:"
137  << pointLocation_().boundaryField().types() << endl;
138  }
139  }
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
154 {
156  (
157  cellDisplacement_,
158  pointDisplacement_
159  );
160 
161  if (pointLocation_.valid())
162  {
163  if (debug)
164  {
165  Info<< "displacementLaplacianFvMotionSolver : applying "
166  << " boundary conditions on " << pointLocation_().name()
167  << " to new point location."
168  << endl;
169  }
170 
171  pointLocation_().internalField() =
172  points0()
173  + pointDisplacement_.internalField();
174 
175  pointLocation_().correctBoundaryConditions();
176 
177  // Implement frozen points
178  if (frozenPointsZone_ != -1)
179  {
180  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
181 
182  forAll(pz, i)
183  {
184  pointLocation_()[pz[i]] = points0()[pz[i]];
185  }
186  }
187 
188  twoDCorrectPoints(pointLocation_().internalField());
189 
190  return tmp<pointField>(pointLocation_().internalField());
191  }
192  else
193  {
194  tmp<pointField> tcurPoints
195  (
196  points0() + pointDisplacement_.internalField()
197  );
198 
199  // Implement frozen points
200  if (frozenPointsZone_ != -1)
201  {
202  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
203 
204  forAll(pz, i)
205  {
206  tcurPoints()[pz[i]] = points0()[pz[i]];
207  }
208  }
209 
210  twoDCorrectPoints(tcurPoints());
211 
212  return tcurPoints;
213  }
214 }
215 
216 
218 {
219  // The points have moved so before interpolation update
220  // the fvMotionSolver accordingly
221  movePoints(fvMesh_.points());
222 
223  diffusivityPtr_->correct();
224  pointDisplacement_.boundaryField().updateCoeffs();
225 
227  (
229  (
230  diffusivityPtr_->operator()(),
231  cellDisplacement_,
232  "laplacian(diffusivity,cellDisplacement)"
233  )
234  );
235 }
236 
237 
239 (
240  const mapPolyMesh& mpm
241 )
242 {
244 
245  // Update diffusivity. Note two stage to make sure old one is de-registered
246  // before creating/registering new one.
247  diffusivityPtr_.reset(NULL);
248  diffusivityPtr_ = motionDiffusivity::New(*this, lookup("diffusivity"));
249 }
250 
251 
252 // ************************ vim: set sw=4 sts=4 et: ************************ //