BUG: If no processor*/constant a parallel run would search up and read the undecompos...
[OpenFOAM-1.6.x.git] / src / fvMotionSolver / fvMotionSolvers / displacement / displacementFvMotionSolver / displacementFvMotionSolver.C
blob61e2f9ad3a03acb661eb6d841512087d54a63bd0
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 "displacementFvMotionSolver.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "mapPolyMesh.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 namespace Foam
35     defineTypeNameAndDebug(displacementFvMotionSolver, 0);
39 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
41 Foam::displacementFvMotionSolver::displacementFvMotionSolver
43     const polyMesh& mesh,
44     Istream&
47     fvMotionSolver(mesh),
48     points0_
49     (
50         pointIOField
51         (
52             IOobject
53             (
54                 "points",
55                 mesh.time().constant(),
56                 polyMesh::meshSubDir,
57                 mesh,
58                 IOobject::MUST_READ,
59                 IOobject::NO_WRITE,
60                 false
61             )
62         )
63     )
65     if (points0_.size() != mesh.nPoints())
66     {
67         FatalErrorIn
68         (
69             "displacementFvMotionSolver::displacementFvMotionSolver\n"
70             "(\n"
71             "    const polyMesh&,\n"
72             "    Istream&\n"
73             ")"
74         )   << "Number of points in mesh " << mesh.nPoints()
75             << " differs from number of points " << points0_.size()
76             << " read from file "
77             <<
78                 IOobject
79                 (
80                     "points",
81                     mesh.time().constant(),
82                     polyMesh::meshSubDir,
83                     mesh,
84                     IOobject::MUST_READ,
85                     IOobject::NO_WRITE,
86                     false
87                 ).filePath()
88             << exit(FatalError);
89     }
93 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
95 Foam::displacementFvMotionSolver::~displacementFvMotionSolver()
99 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
101 void Foam::displacementFvMotionSolver::updateMesh(const mapPolyMesh& mpm)
103     fvMotionSolver::updateMesh(mpm);
105     // Map points0_. Bit special since we somehow have to come up with
106     // a sensible points0 position for introduced points.
107     // Find out scaling between points0 and current points
109     // Get the new points either from the map or the mesh
110     const pointField& points =
111     (
112         mpm.hasMotionPoints()
113       ? mpm.preMotionPoints()
114       : fvMesh_.points()
115     );
117     // Note: boundBox does reduce
118     const vector span0 = boundBox(points0_).span();
119     const vector span  = boundBox(points).span();
121     vector scaleFactors(cmptDivide(span0, span));
123     pointField newPoints0(mpm.pointMap().size());
125     forAll(newPoints0, pointI)
126     {
127         label oldPointI = mpm.pointMap()[pointI];
129         if (oldPointI >= 0)
130         {
131             label masterPointI = mpm.reversePointMap()[oldPointI];
133             if (masterPointI == pointI)
134             {
135                 newPoints0[pointI] = points0_[oldPointI];
136             }
137             else
138             {
139                 // New point. Assume motion is scaling.
140                 newPoints0[pointI] = points0_[oldPointI] + cmptMultiply
141                 (
142                     scaleFactors,
143                     points[pointI]-points[masterPointI]
144                 );
145             }
146         }
147         else
148         {
149             FatalErrorIn
150             (
151                 "displacementLaplacianFvMotionSolver::updateMesh"
152                 "(const mapPolyMesh& mpm)"
153             )   << "Cannot work out coordinates of introduced vertices."
154                 << " New vertex " << pointI << " at coordinate "
155                 << points[pointI] << exit(FatalError);
156         }
157     }
158     points0_.transfer(newPoints0);
162 // ************************************************************************* //