removeFaces.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 | Copyright (C) 2011-2015 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 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 
50 int main(int argc, char *argv[])
51 {
52  #include "addOverwriteOption.H"
53  argList::validArgs.append("faceSet");
54 
55  #include "setRootCase.H"
56  #include "createTime.H"
57  runTime.functionObjects().off();
58  #include "createMesh.H"
59  const word oldInstance = mesh.pointsInstance();
60 
61  const word setName = args[1];
62  const bool overwrite = args.optionFound("overwrite");
63 
64  // Read faces
65  faceSet candidateSet(mesh, setName);
66 
67  Pout<< "Read " << candidateSet.size() << " faces to remove" << nl
68  << endl;
69 
70 
71  labelList candidates(candidateSet.toc());
72 
73  // Face removal engine. No checking for not merging boundary faces.
74  removeFaces faceRemover(mesh, 2);
75 
76  // Get compatible set of faces and connected sets of cells.
77  labelList cellRegion;
78  labelList cellRegionMaster;
79  labelList facesToRemove;
80 
81  faceRemover.compatibleRemoves
82  (
83  candidates,
84  cellRegion,
85  cellRegionMaster,
86  facesToRemove
87  );
88 
89  {
90  faceSet compatibleRemoves(mesh, "compatibleRemoves", facesToRemove);
91 
92  Pout<< "Original faces to be removed:" << candidateSet.size() << nl
93  << "New faces to be removed:" << compatibleRemoves.size() << nl
94  << endl;
95 
96  Pout<< "Writing new faces to be removed to faceSet "
97  << compatibleRemoves.instance()
98  /compatibleRemoves.local()
99  /compatibleRemoves.name()
100  << endl;
101 
102  compatibleRemoves.write();
103  }
104 
105 
106  // Read objects in time directory
107  IOobjectList objects(mesh, runTime.timeName());
108 
109  // Read vol fields.
111  ReadFields(mesh, objects, vsFlds);
112 
114  ReadFields(mesh, objects, vvFlds);
115 
117  ReadFields(mesh, objects, vstFlds);
118 
119  PtrList<volSymmTensorField> vsymtFlds;
120  ReadFields(mesh, objects, vsymtFlds);
121 
123  ReadFields(mesh, objects, vtFlds);
124 
125  // Read surface fields.
127  ReadFields(mesh, objects, ssFlds);
128 
130  ReadFields(mesh, objects, svFlds);
131 
133  ReadFields(mesh, objects, sstFlds);
134 
136  ReadFields(mesh, objects, ssymtFlds);
137 
139  ReadFields(mesh, objects, stFlds);
140 
141 
142  // Topo changes container
143  polyTopoChange meshMod(mesh);
144 
145  // Insert mesh refinement into polyTopoChange.
146  faceRemover.setRefinement
147  (
148  facesToRemove,
149  cellRegion,
150  cellRegionMaster,
151  meshMod
152  );
153 
154  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
155 
156  mesh.updateMesh(morphMap);
157 
158  // Move mesh (since morphing does not do this)
159  if (morphMap().hasMotionPoints())
160  {
161  mesh.movePoints(morphMap().preMotionPoints());
162  }
163 
164  // Update numbering of cells/vertices.
165  faceRemover.updateMesh(morphMap);
166 
167  if (!overwrite)
168  {
169  runTime++;
170  }
171  else
172  {
173  mesh.setInstance(oldInstance);
174  }
175 
176  // Take over refinement levels and write to new time directory.
177  Pout<< "Writing mesh to time " << runTime.timeName() << endl;
178  mesh.write();
179 
180  Pout<< "End\n" << endl;
181 
182  return 0;
183 }
184 
185 
186 // ************************************************************************* //
Foam::surfaceFields.
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
Field reading functions for post-processing utilities.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:800
dynamicFvMesh & mesh
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
A class for handling words, derived from string.
Definition: word.H:59
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
static const char nl
Definition: Ostream.H:262
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:721
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
Direct mesh changes based on v1.3 polyTopoChange syntax.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::argList args(argc, argv)
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
Namespace for OpenFOAM.