points0MotionSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2023 OpenFOAM Foundation
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 
26 #include "points0MotionSolver.H"
27 #include "polyDistributionMap.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34 }
35 
36 
38 (
39  const polyMesh& mesh
40 )
41 {
42  const word instance
43  (
45  (
46  ".",
47  "points0",
49  )
50  );
51 
52  if (instance != mesh.time().constant())
53  {
54  // Points0 written to a time folder
55 
56  return pointVectorField
57  (
58  IOobject
59  (
60  "points0",
61  instance,
62  mesh,
65  false
66  ),
68  );
69  }
70  else
71  {
72  // Return copy of original mesh points
73 
75  (
76  IOobject
77  (
78  "points",
79  mesh.time().constant(),
81  mesh,
84  false
85  )
86  );
87 
89  (
90  IOobject
91  (
92  "points0",
93  instance,
94  mesh,
97  false
98  ),
101  );
102 
103  points0.primitiveFieldRef() = points;
104 
105  return points0;
106  }
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111 
113 (
114  const word& name,
115  const polyMesh& mesh,
116  const dictionary& dict,
117  const word& type
118 )
119 :
120  motionSolver(name, mesh, dict, type),
121  points0_(readPoints0(mesh))
122 {
123  if (points0_.size() != mesh.nPoints())
124  {
126  << "Number of points in mesh " << mesh.nPoints()
127  << " differs from number of points " << points0_.size()
128  << " read from file "
130  (
131  "points",
132  mesh.time().constant(),
134  mesh,
137  false
138  ).filePath()
139  << exit(FatalError);
140  }
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {}
154 
155 
157 {
159 }
160 
161 
163 {
164  points0_.primitiveFieldRef() = mesh().points();
165 
166  // The processor boundaries may have changed, so we need to update the
167  // boundary field. There is no data in this field, so we don't need to map
168  // anything. We can just reset it to a freshly created calculated field.
169  points0_.boundaryFieldRef().reset
170  (
172  (
173  points0_.mesh().boundary(),
174  points0_.internalField(),
175  calculatedPointPatchVectorField::typeName
176  )
177  );
178 }
179 
180 
182 (
183  const polyDistributionMap& map
184 )
185 {
186  map.distributePointData(points0_.primitiveFieldRef());
187 
188  // See above
189  points0_.boundaryFieldRef().reset
190  (
192  (
193  points0_.mesh().boundary(),
194  points0_.internalField(),
195  calculatedPointPatchVectorField::typeName
196  )
197  );
198 }
199 
200 
202 {
203  if (points0_.writeOpt() == IOobject::AUTO_WRITE)
204  {
205  points0_.write();
206  }
207 
208  return motionSolver::write();
209 }
210 
211 
212 // ************************************************************************* //
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricBoundaryField class.
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Return the location of "dir" containing the file "name".
Definition: Time.C:658
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
virtual bool write() const
Optionally write motion state information for restart.
Definition: motionSolver.C:128
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:137
const Time & time() const
Return time.
Virtual base class for displacement motion solvers.
points0MotionSolver(const word &name, const polyMesh &, const dictionary &, const word &type)
Construct from mesh and dictionary.
virtual ~points0MotionSolver()
Destructor.
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
static pointVectorField readPoints0(const polyMesh &mesh)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual bool write() const
Write points0 if the mesh topology changed.
pointVectorField points0_
Starting points.
pointField & points0()
Return reference to the reference field.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
void distributePointData(List< T > &lst) const
Distribute list of point data.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:270
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
label nPoints() const
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
fileName filePath() const
Return the path for the file for this Type.
Definition: IOobject.H:563
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
defineTypeNameAndDebug(combustionModel, 0)
PointField< vector > pointVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
error FatalError
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict