removeFaces.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) 2011-2019 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 Application
25  removeFaces
26 
27 Description
28  Utility to remove faces (combines cells on both sides).
29 
30  Takes faceSet of candidates for removal and writes faceSet with faces that
31  will actually be removed. (because e.g. would cause two faces between the
32  same cells). See removeFaces in dynamicMesh library for constraints.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "argList.H"
37 #include "Time.H"
38 #include "polyTopoChange.H"
39 #include "faceSet.H"
40 #include "removeFaces.H"
41 #include "ReadFields.H"
42 #include "volFields.H"
43 #include "surfaceFields.H"
44 #include "pointFields.H"
45 
46 using namespace Foam;
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  #include "addOverwriteOption.H"
53  argList::validArgs.append("faceSet");
55  (
56  "noFields",
57  "do not update fields"
58  );
59 
60 
61  #include "setRootCase.H"
62  #include "createTime.H"
64 
65  const bool overwrite = args.optionFound("overwrite");
66  const bool fields = !args.optionFound("noFields");
67 
68  #include "createMesh.H"
69  const word oldInstance = mesh.pointsInstance();
70 
71  const word setName = args[1];
72 
73  // Read faces
74  faceSet candidateSet(mesh, setName);
75 
76  Pout<< "Read " << candidateSet.size() << " faces to remove" << nl
77  << endl;
78 
79 
80  labelList candidates(candidateSet.toc());
81 
82  // Face removal engine. No checking for not merging boundary faces.
83  removeFaces faceRemover(mesh, 2);
84 
85  // Get compatible set of faces and connected sets of cells.
86  labelList cellRegion;
87  labelList cellRegionMaster;
88  labelList facesToRemove;
89 
90  faceRemover.compatibleRemoves
91  (
92  candidates,
93  cellRegion,
94  cellRegionMaster,
95  facesToRemove
96  );
97 
98  {
99  faceSet compatibleRemoves(mesh, "compatibleRemoves", facesToRemove);
100 
101  Pout<< "Original faces to be removed:" << candidateSet.size() << nl
102  << "New faces to be removed:" << compatibleRemoves.size() << nl
103  << endl;
104 
105  Pout<< "Writing new faces to be removed to faceSet "
106  << compatibleRemoves.instance()
107  /compatibleRemoves.local()
108  /compatibleRemoves.name()
109  << endl;
110 
111  compatibleRemoves.write();
112  }
113 
114 
115  // Read objects in time directory
117 
118  if (fields) Info<< "Reading geometric fields" << nl << endl;
119 
120  #include "readVolFields.H"
121  #include "readSurfaceFields.H"
122  #include "readPointFields.H"
123 
124  Info<< endl;
125 
126 
127  // Topo changes container
128  polyTopoChange meshMod(mesh);
129 
130  // Insert mesh refinement into polyTopoChange.
131  faceRemover.setRefinement
132  (
133  facesToRemove,
134  cellRegion,
135  cellRegionMaster,
136  meshMod
137  );
138 
139  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
140 
141  mesh.updateMesh(morphMap);
142 
143  // Move mesh (since morphing does not do this)
144  if (morphMap().hasMotionPoints())
145  {
146  mesh.movePoints(morphMap().preMotionPoints());
147  }
148 
149  // Update numbering of cells/vertices.
150  faceRemover.updateMesh(morphMap);
151 
152  if (!overwrite)
153  {
154  runTime++;
155  }
156  else
157  {
158  mesh.setInstance(oldInstance);
159  }
160 
161  // Take over refinement levels and write to new time directory.
162  Pout<< "Writing mesh to time " << runTime.timeName() << endl;
163  mesh.write();
164 
165  Pout<< "End\n" << endl;
166 
167  return 0;
168 }
169 
170 
171 // ************************************************************************* //
Foam::surfaceFields.
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:82
Given list of faces to remove insert all the topology changes. Contains helper function to get consis...
Definition: removeFaces.H:62
A list of face labels.
Definition: faceSet.H:48
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
void off()
Switch the function objects off.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
Field reading functions for post-processing utilities.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1035
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:802
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:795
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
static const char nl
Definition: Ostream.H:265
objects
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:716
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:410
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Direct mesh changes based on v1.3 polyTopoChange syntax.
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual Ostream & write(const token &)=0
Write next token to stream.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:117
Foam::argList args(argc, argv)
Namespace for OpenFOAM.