ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / fvMotionSolver / fvMotionSolvers / displacement / displacementFvMotionSolver / displacementFvMotionSolver.C
blob690c4e1a5565b9e8c37d9a2af08c179e2fb6b30f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "displacementFvMotionSolver.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "mapPolyMesh.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 namespace Foam
34     defineTypeNameAndDebug(displacementFvMotionSolver, 0);
38 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
40 Foam::displacementFvMotionSolver::displacementFvMotionSolver
42     const polyMesh& mesh,
43     Istream&
46     fvMotionSolver(mesh),
47     points0_
48     (
49         pointIOField
50         (
51             IOobject
52             (
53                 "points",
54                 mesh.time().constant(),
55                 polyMesh::meshSubDir,
56                 mesh,
57                 IOobject::MUST_READ,
58                 IOobject::NO_WRITE,
59                 false
60             )
61         )
62     )
64     if (points0_.size() != mesh.nPoints())
65     {
66         FatalErrorIn
67         (
68             "displacementFvMotionSolver::displacementFvMotionSolver\n"
69             "(\n"
70             "    const polyMesh&,\n"
71             "    Istream&\n"
72             ")"
73         )   << "Number of points in mesh " << mesh.nPoints()
74             << " differs from number of points " << points0_.size()
75             << " read from file "
76             <<
77                 IOobject
78                 (
79                     "points",
80                     mesh.time().constant(),
81                     polyMesh::meshSubDir,
82                     mesh,
83                     IOobject::MUST_READ,
84                     IOobject::NO_WRITE,
85                     false
86                 ).filePath()
87             << exit(FatalError);
88     }
92 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
94 Foam::displacementFvMotionSolver::~displacementFvMotionSolver()
98 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
100 void Foam::displacementFvMotionSolver::updateMesh(const mapPolyMesh& mpm)
102     fvMotionSolver::updateMesh(mpm);
104     // Map points0_. Bit special since we somehow have to come up with
105     // a sensible points0 position for introduced points.
106     // Find out scaling between points0 and current points
108     // Get the new points either from the map or the mesh
109     const pointField& points =
110     (
111         mpm.hasMotionPoints()
112       ? mpm.preMotionPoints()
113       : fvMesh_.points()
114     );
116     // Note: boundBox does reduce
117     const vector span0 = boundBox(points0_).span();
118     const vector span  = boundBox(points).span();
120     vector scaleFactors(cmptDivide(span0, span));
122     pointField newPoints0(mpm.pointMap().size());
124     forAll(newPoints0, pointI)
125     {
126         label oldPointI = mpm.pointMap()[pointI];
128         if (oldPointI >= 0)
129         {
130             label masterPointI = mpm.reversePointMap()[oldPointI];
132             if (masterPointI == pointI)
133             {
134                 newPoints0[pointI] = points0_[oldPointI];
135             }
136             else
137             {
138                 // New point. Assume motion is scaling.
139                 newPoints0[pointI] = points0_[oldPointI] + cmptMultiply
140                 (
141                     scaleFactors,
142                     points[pointI]-points[masterPointI]
143                 );
144             }
145         }
146         else
147         {
148             FatalErrorIn
149             (
150                 "displacementLaplacianFvMotionSolver::updateMesh"
151                 "(const mapPolyMesh& mpm)"
152             )   << "Cannot work out coordinates of introduced vertices."
153                 << " New vertex " << pointI << " at coordinate "
154                 << points[pointI] << exit(FatalError);
155         }
156     }
157     points0_.transfer(newPoints0);
161 // ************************************************************************* //