initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / mesh / manipulation / createBaffles / createBaffles.C
blob0ba1af5beaa04787aecb9871b5e0f36a2719b22d
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 Description
26     Makes internal faces into boundary faces. Does not duplicate points, unlike
27     mergeOrSplitBaffles.
29     Note: if any coupled patch face is selected for baffling the opposite
30     member has to be selected for baffling as well. Note that this
31     is the same as repatching. This was added only for convenience so
32     you don't have to filter coupled boundary out of your set.
34 \*---------------------------------------------------------------------------*/
36 #include "syncTools.H"
37 #include "argList.H"
38 #include "Time.H"
39 #include "faceSet.H"
40 #include "polyTopoChange.H"
41 #include "polyModifyFace.H"
42 #include "polyAddFace.H"
43 #include "ReadFields.H"
44 #include "volFields.H"
45 #include "surfaceFields.H"
46 #include "ZoneIDs.H"
48 using namespace Foam;
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 void modifyOrAddFace
54     polyTopoChange& meshMod,
55     const face& f,
56     const label faceI,
57     const label own,
58     const bool flipFaceFlux,
59     const label newPatchI,
60     const label zoneID,
61     const bool zoneFlip,
63     PackedBoolList& modifiedFace
66     if (!modifiedFace[faceI])
67     {
68         // First usage of face. Modify.
69         meshMod.setAction
70         (
71             polyModifyFace
72             (
73                 f,                          // modified face
74                 faceI,                      // label of face
75                 own,                        // owner
76                 -1,                         // neighbour
77                 flipFaceFlux,               // face flip
78                 newPatchI,                  // patch for face
79                 false,                      // remove from zone
80                 zoneID,                     // zone for face
81                 zoneFlip                    // face flip in zone
82             )
83         );
84         modifiedFace[faceI] = 1;
85     }
86     else
87     {
88         // Second or more usage of face. Add.
89         meshMod.setAction
90         (
91             polyAddFace
92             (
93                 f,                          // modified face
94                 own,                        // owner
95                 -1,                         // neighbour
96                 -1,                         // master point
97                 -1,                         // master edge
98                 faceI,                      // master face
99                 flipFaceFlux,               // face flip
100                 newPatchI,                  // patch for face
101                 zoneID,                     // zone for face
102                 zoneFlip                    // face flip in zone
103             )
104         );
105     }
109 label findPatchID(const polyMesh& mesh, const word& name)
111     label patchI = mesh.boundaryMesh().findPatchID(name);
113     if (patchI == -1)
114     {
115         FatalErrorIn("findPatchID(const polyMesh&, const word&)")
116             << "Cannot find patch " << name << endl
117             << "Valid patches are " << mesh.boundaryMesh().names()
118             << exit(FatalError);
119     }
120     return patchI;
124 // Main program:
126 int main(int argc, char *argv[])
128     argList::validArgs.append("faceZone");
129     argList::validArgs.append("patch");
130     argList::validOptions.insert("additionalPatches", "(patch2 .. patchN)");
131     argList::validOptions.insert("internalFacesOnly", "");
132     argList::validOptions.insert("overwrite", "");
134 #   include "setRootCase.H"
135 #   include "createTime.H"
136     runTime.functionObjects().off();
137 #   include "createMesh.H"
138     const word oldInstance = mesh.pointsInstance();
140     const polyBoundaryMesh& patches = mesh.boundaryMesh();
141     const faceZoneMesh& faceZones = mesh.faceZones();
143     // Faces to baffle
144     faceZoneID zoneID(args.additionalArgs()[0], faceZones);
146     Info<< "Converting faces on zone " << zoneID.name()
147         << " into baffles." << nl << endl;
149     const faceZone& fZone = faceZones[zoneID.index()];
151     Info<< "Found " << returnReduce(fZone.size(), sumOp<label>())
152         << " faces on zone " << zoneID.name() << nl << endl;
154     // Make sure patches and zoneFaces are synchronised across couples
155     patches.checkParallelSync(true);
156     fZone.checkParallelSync(true);
158     // Patches to put baffles into
159     DynamicList<label> newPatches(1);
161     word patchName(args.additionalArgs()[1]);
162     newPatches.append(findPatchID(mesh, patchName));
163     Info<< "Using patch " << patchName
164         << " at index " << newPatches[0] << endl;
167     // Additional patches
168     if (args.optionFound("additionalPatches"))
169     {
170         const wordList patchNames
171         (
172             args.optionLookup("additionalPatches")()
173         );
175         newPatches.reserve(patchNames.size() + 1);
176         forAll(patchNames, i)
177         {
178             newPatches.append(findPatchID(mesh, patchNames[i]));
179             Info<< "Using additional patch " << patchNames[i]
180                 << " at index " << newPatches[newPatches.size()-1] << endl;
181         }
182     }
185     bool overwrite = args.optionFound("overwrite");
187     bool internalFacesOnly = args.optionFound("internalFacesOnly");
189     if (internalFacesOnly)
190     {
191         Info<< "Not converting faces on non-coupled patches." << nl << endl;
192     }
195     // Read objects in time directory
196     IOobjectList objects(mesh, runTime.timeName());
198     // Read vol fields.
200     PtrList<volScalarField> vsFlds;
201     ReadFields(mesh, objects, vsFlds);
203     PtrList<volVectorField> vvFlds;
204     ReadFields(mesh, objects, vvFlds);
206     PtrList<volSphericalTensorField> vstFlds;
207     ReadFields(mesh, objects, vstFlds);
209     PtrList<volSymmTensorField> vsymtFlds;
210     ReadFields(mesh, objects, vsymtFlds);
212     PtrList<volTensorField> vtFlds;
213     ReadFields(mesh, objects, vtFlds);
215     // Read surface fields.
217     PtrList<surfaceScalarField> ssFlds;
218     ReadFields(mesh, objects, ssFlds);
220     PtrList<surfaceVectorField> svFlds;
221     ReadFields(mesh, objects, svFlds);
223     PtrList<surfaceSphericalTensorField> sstFlds;
224     ReadFields(mesh, objects, sstFlds);
226     PtrList<surfaceSymmTensorField> ssymtFlds;
227     ReadFields(mesh, objects, ssymtFlds);
229     PtrList<surfaceTensorField> stFlds;
230     ReadFields(mesh, objects, stFlds);
233     // Mesh change container
234     polyTopoChange meshMod(mesh);
237     // Do the actual changes. Note:
238     // - loop in incrementing face order (not necessary if faceZone ordered).
239     //   Preserves any existing ordering on patch faces.
240     // - two passes, do non-flip faces first and flip faces second. This
241     //   guarantees that when e.g. creating a cyclic all faces from one
242     //   side come first and faces from the other side next.
244     // Whether first use of face (modify) or consecutive (add)
245     PackedBoolList modifiedFace(mesh.nFaces());
246     // Never modify coupled faces
247     forAll(patches, patchI)
248     {
249         const polyPatch& pp = patches[patchI];
250         if (pp.coupled())
251         {
252             forAll(pp, i)
253             {
254                 modifiedFace[pp.start()+i] = 1;
255             }
256         }
257     }
258     label nModified = 0;
260     forAll(newPatches, i)
261     {
262         label newPatchI = newPatches[i];
264         // Pass 1. Do selected side of zone
265         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
267         for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
268         {
269             label zoneFaceI = fZone.whichFace(faceI);
271             if (zoneFaceI != -1)
272             {
273                 if (!fZone.flipMap()[zoneFaceI])
274                 {
275                     // Use owner side of face
276                     modifyOrAddFace
277                     (
278                         meshMod,
279                         mesh.faces()[faceI],    // modified face
280                         faceI,                  // label of face
281                         mesh.faceOwner()[faceI],// owner
282                         false,                  // face flip
283                         newPatchI,              // patch for face
284                         zoneID.index(),         // zone for face
285                         false,                  // face flip in zone
286                         modifiedFace            // modify or add status
287                     );
288                 }
289                 else
290                 {
291                     // Use neighbour side of face
292                     modifyOrAddFace
293                     (
294                         meshMod,
295                         mesh.faces()[faceI].reverseFace(),  // modified face
296                         faceI,                      // label of face
297                         mesh.faceNeighbour()[faceI],// owner
298                         true,                       // face flip
299                         newPatchI,                  // patch for face
300                         zoneID.index(),             // zone for face
301                         true,                       // face flip in zone
302                         modifiedFace                // modify or add status
303                     );
304                 }
306                 nModified++;
307             }
308         }
311         // Pass 2. Do other side of zone
312         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
314         for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
315         {
316             label zoneFaceI = fZone.whichFace(faceI);
318             if (zoneFaceI != -1)
319             {
320                 if (!fZone.flipMap()[zoneFaceI])
321                 {
322                     // Use neighbour side of face
323                     modifyOrAddFace
324                     (
325                         meshMod,
326                         mesh.faces()[faceI].reverseFace(),  // modified face
327                         faceI,                              // label of face
328                         mesh.faceNeighbour()[faceI],        // owner
329                         true,                               // face flip
330                         newPatchI,                          // patch for face
331                         zoneID.index(),                     // zone for face
332                         true,                               // face flip in zone
333                         modifiedFace                        // modify or add
334                     );
335                 }
336                 else
337                 {
338                     // Use owner side of face
339                     modifyOrAddFace
340                     (
341                         meshMod,
342                         mesh.faces()[faceI],    // modified face
343                         faceI,                  // label of face
344                         mesh.faceOwner()[faceI],// owner
345                         false,                  // face flip
346                         newPatchI,              // patch for face
347                         zoneID.index(),         // zone for face
348                         false,                  // face flip in zone
349                         modifiedFace            // modify or add status
350                     );
351                 }
352             }
353         }
356         // Modify any boundary faces
357         // ~~~~~~~~~~~~~~~~~~~~~~~~~
359         // Normal boundary:
360         // - move to new patch. Might already be back-to-back baffle
361         // you want to add cyclic to. Do warn though.
362         //
363         // Processor boundary:
364         // - do not move to cyclic
365         // - add normal patches though.
367         // For warning once per patch.
368         labelHashSet patchWarned;
370         forAll(patches, patchI)
371         {
372             const polyPatch& pp = patches[patchI];
374             if (pp.coupled() && patches[newPatchI].coupled())
375             {
376                 // Do not allow coupled faces to be moved to different coupled
377                 // patches.
378             }
379             else if (pp.coupled() || !internalFacesOnly)
380             {
381                 forAll(pp, i)
382                 {
383                     label faceI = pp.start()+i;
385                     label zoneFaceI = fZone.whichFace(faceI);
387                     if (zoneFaceI != -1)
388                     {
389                         if (patchWarned.insert(patchI))
390                         {
391                             WarningIn(args.executable())
392                                 << "Found boundary face (in patch " << pp.name()
393                                 << ") in faceZone " << fZone.name()
394                                 << " to convert to baffle patch "
395                                 << patches[newPatchI].name()
396                                 << endl
397                                 << "    Run with -internalFacesOnly option"
398                                 << " if you don't wish to convert"
399                                 << " boundary faces." << endl;
400                         }
402                         modifyOrAddFace
403                         (
404                             meshMod,
405                             mesh.faces()[faceI],        // modified face
406                             faceI,                      // label of face
407                             mesh.faceOwner()[faceI],    // owner
408                             false,                      // face flip
409                             newPatchI,                  // patch for face
410                             zoneID.index(),             // zone for face
411                             fZone.flipMap()[zoneFaceI], // face flip in zone
412                             modifiedFace                // modify or add status
413                         );
414                         nModified++;
415                     }
416                 }
417             }
418         }
419     }
422     Info<< "Converted " << returnReduce(nModified, sumOp<label>())
423         << " faces into boundary faces on patch " << patchName << nl << endl;
425     if (!overwrite)
426     {
427         runTime++;
428     }
430     // Change the mesh. Change points directly (no inflation).
431     autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
433     // Update fields
434     mesh.updateMesh(map);
436     // Move mesh (since morphing might not do this)
437     if (map().hasMotionPoints())
438     {
439         mesh.movePoints(map().preMotionPoints());
440     }
442     if (overwrite)
443     {
444         mesh.setInstance(oldInstance);
445     }
446     Info<< "Writing mesh to " << runTime.timeName() << endl;
448     mesh.write();
450     Info<< "End\n" << endl;
452     return 0;
456 // ************************************************************************* //