cumulativeDisplacement.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2007-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Application
29  cumulativeDisplacement
30 
31 Description
32  Computes and writes the difference between the mesh points in each time
33  instance and the ones in the constant folder. Additionally, the projection
34  of this difference to the normal point vectors of the initial mesh is also
35  written
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "fvCFD.H"
40 #include "emptyFvPatch.H"
41 #include "coupledFvPatch.H"
42 #include "pointMesh.H"
43 #include "pointPatchField.H"
44 #include "pointPatchFieldsFwd.H"
45 #include "syncTools.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
51  timeSelector::addOptions();
52 
53  #include "addRegionOption.H"
54 
55  #include "setRootCase.H"
56  #include "createTime.H"
57  instantList timeDirs = timeSelector::select0(runTime, args);
58  #include "createNamedMesh.H"
59  #include "createFields.H"
60 
61  // polyPatch::pointNormals will give the wrong result for points
62  // belonging to multiple patches or patch-processorPatch intersections.
63  // Keeping a mesh-wide field to allow easy reduction using syncTools.
64  // A bit expensive? Better way?
65  vectorField pointNormals(mesh.nPoints(), vector::zero);
66  for (const fvPatch& patch : mesh.boundary())
67  {
68  // Each local patch point belongs to these local patch faces.
69  // Local numbering
70  const labelListList& patchPointFaces = patch.patch().pointFaces();
71  if (!isA<coupledFvPatch>(patch) && !isA<emptyFvPatch>(patch))
72  {
73  const labelList& meshPoints = patch.patch().meshPoints();
74  const vectorField nf(patch.nf());
75  forAll(meshPoints, ppI)
76  {
77  const labelList& pointFaces = patchPointFaces[ppI];
78  forAll(pointFaces, pfI)
79  {
80  const label& localFaceIndex = pointFaces[pfI];
81  pointNormals[meshPoints[ppI]] += nf[localFaceIndex];
82  }
83  }
84  }
85  }
86  // Sum from all processors
87  syncTools::syncPointList
88  (
89  mesh, pointNormals, plusEqOp<vector>(), vector::zero
90  );
91  pointNormals /= mag(pointNormals) + VSMALL;
92 
93  forAll(timeDirs, timeI)
94  {
95  runTime.setTime(timeDirs[timeI], timeI);
96  Info<< "Time = " << runTime.timeName() << endl;
97  mesh.readUpdate();
98 
99  const pointMesh& pMesh(pointMesh::New(mesh));
100  // Point displacement projected to the
101  // unit point normal of the initial geometry
102  pointScalarField normalDisplacement
103  (
104  IOobject
105  (
106  "normalDisplacement",
107  runTime.timeName(),
108  mesh,
109  IOobject::NO_READ,
110  IOobject::AUTO_WRITE
111  ),
112  pMesh,
114  );
115  // Point displacement as a vector
116  pointVectorField displacement
117  (
118  IOobject
119  (
120  "displacement",
121  runTime.timeName(),
122  mesh,
123  IOobject::NO_READ,
124  IOobject::AUTO_WRITE
125  ),
126  pMesh,
128  );
129 
130  Info<< "Calculating cumulative mesh movement for time "
131  << runTime.timeName() << endl;
132  // Normal displacement
133  const pointField& meshPoints = mesh.points();
134  forAll(mesh.boundary(), pI)
135  {
136  const polyPatch& patch = mesh.boundaryMesh()[pI];
137  const labelList& localPoints = patch.meshPoints();
138  forAll(localPoints, ppI)
139  {
140  label pointI = localPoints[ppI];
141  normalDisplacement[pointI] =
142  (meshPoints[pointI] - points0[pointI])
143  & pointNormals[pointI];
144  }
145  }
146  normalDisplacement.write();
147  // Vectorial displacement
148  displacement.primitiveFieldRef() = meshPoints - points0;
149  displacement.write();
150  }
151 
152  // Print execution time
153  mesh.time().printExecutionTime(Info);
154 
155  Info<< "End\n" << endl;
156 
157  return 0;
158 }
159 
160 
161 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointVectorField
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Definition: pointFieldsFwd.H:61
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
pointPatchField.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:44
pointPatchFieldsFwd.H
syncTools.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
coupledFvPatch.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
addRegionOption.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
createNamedMesh.H
Required Variables.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
emptyFvPatch.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
setRootCase.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
createTime.H
fvCFD.H
args
Foam::argList args(argc, argv)
pointMesh.H
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::pointScalarField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFieldsFwd.H:56