fvMeshSubset.H 10.6 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
6
     \\/     M anipulation  | Copyright (C) 2016 OpenCFD Ltd.
7
8
9
10
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

11
12
13
14
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
15
16
17
18
19
20
21

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    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
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43

Class
    Foam::fvMeshSubset

Description
    Post-processing mesh subset tool.  Given the original mesh and the
    list of selected cells, it creates the mesh consisting only of the
    desired cells, with the mapping list for points, faces, and cells.

    Puts all exposed internal faces into either
    - a user supplied patch
    - a newly created patch "oldInternalFaces"

    - setCellSubset is for small subsets. Uses Maps to minimize memory.
    - setLargeCellSubset is for largish subsets (>10% of mesh).
      Uses labelLists instead.

    - setLargeCellSubset does coupled patch subsetting as well. If it detects
      a face on a coupled patch 'losing' its neighbour it will move the
      face into the oldInternalFaces patch.

44
45
46
47
48
49
    - if a user supplied patch is used it is up to the destination
      patchField to handle exposed internal faces (mapping from face -1).
      If not provided the default is to assign the internalField. All the
      basic patch field types (e.g. fixedValue) will give a warning and
      preferably derived patch field types should be used that know how to
      handle exposed faces (e.g. use uniformFixedValue instead of fixedValue)
50
51
52
53
54
55
56
57
58
59

SourceFiles
    fvMeshSubset.C

\*---------------------------------------------------------------------------*/

#ifndef fvMeshSubset_H
#define fvMeshSubset_H

#include "fvMesh.H"
Mattijs Janssens's avatar
Mattijs Janssens committed
60
#include "pointMesh.H"
61
#include "GeometricField.H"
62
#include "HashSet.H"
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
#include "surfaceMesh.H"

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{

/*---------------------------------------------------------------------------*\
                        Class fvMeshSubset Declaration
\*---------------------------------------------------------------------------*/

class fvMeshSubset
{

private:

    // Private data

        //- Mesh to subset from
        const fvMesh& baseMesh_;

        //- Subset mesh pointer
Mattijs Janssens's avatar
Mattijs Janssens committed
85
86
        autoPtr<fvMesh> fvMeshSubsetPtr_;

87
88
89
90
91
92
93
94
95
96
97
98
        //- Point mapping array
        labelList pointMap_;

        //- Face mapping array
        labelList faceMap_;

        //- Cell mapping array
        labelList cellMap_;

        //- Patch mapping array
        labelList patchMap_;

99
        //- Optional face mapping array with flip encoded
100
101
        mutable autoPtr<labelList> faceFlipMapPtr_;

102
103
104
105
106
107
108

    // Private Member Functions

        //- Check if subset has been performed
        bool checkCellSubset() const;

        //- Mark points in Map
109
        static void markPoints(const labelList&, Map<label>&);
110
111

        //- Mark points (with 0) in labelList
112
        static void markPoints(const labelList&, labelList&);
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131

        //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
        void doCoupledPatches
        (
            const bool syncPar,
            labelList& nCellsUsingFace
        ) const;

        //- Subset of subset
        static labelList subset
        (
            const label nElems,
            const labelList& selectedElements,
            const labelList& subsetMap
        );

        //- Create zones for submesh
        void subsetZones();

132
133
134
135
136
137
138
        //- Helper: extract cells-to-remove from cells-to-keep
        labelList getCellsToRemove
        (
            const labelList& region,
            const label currentRegion
        ) const;

139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
        //- Disallow default bitwise copy construct
        fvMeshSubset(const fvMeshSubset&);

        //- Disallow default bitwise assignment
        void operator=(const fvMeshSubset&);

public:

    // Constructors

        //- Construct given a mesh to subset
        explicit fvMeshSubset(const fvMesh&);


    // Member Functions

        // Edit

            //- Set the subset. Create "oldInternalFaces" patch for exposed
            //  internal faces (patchID==-1) or use supplied patch.
            //  Does not handle coupled patches correctly if only one side
            //  gets deleted.
            void setCellSubset
            (
                const labelHashSet& globalCellMap,
                const label patchID = -1
            );

            //- Set the subset from all cells with region == currentRegion.
            //  Create "oldInternalFaces" patch for exposed
            //  internal faces (patchID==-1) or use supplied patch.
andy's avatar
andy committed
170
            //  Handles coupled patches by if necessary making coupled patch
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
            //  face part of patchID (so uncoupled)
            void setLargeCellSubset
            (
                const labelList& region,
                const label currentRegion,
                const label patchID = -1,
                const bool syncCouples = true
            );

            //- setLargeCellSubset but with labelHashSet.
            void setLargeCellSubset
            (
                const labelHashSet& globalCellMap,
                const label patchID = -1,
                const bool syncPar = true
            );


189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
            //- Two step subsetting

                //- Get labels of exposed faces.
                //  These are
                //  - internal faces that become boundary faces
                //  - coupled faces that become uncoupled (since one of the
                //    sides gets deleted)
                labelList getExposedFaces
                (
                    const labelList& region,
                    const label currentRegion,
                    const bool syncCouples = true
                ) const;

                //- For every exposed face (from above getExposedFaces)
                //  used supplied (existing!) patch
                void setLargeCellSubset
                (
                    const labelList& region,
                    const label currentRegion,
                    const labelList& exposedFaces,
                    const labelList& patchIDs,
                    const bool syncCouples = true
                );


215
216
217
218
219
220
221
222
        // Access

            //- Original mesh
            const fvMesh& baseMesh() const
            {
                return baseMesh_;
            }

mattijs's avatar
mattijs committed
223
224
225
            //- Have subMesh?
            bool hasSubMesh() const;

226
227
            //- Return reference to subset mesh
            const fvMesh& subMesh() const;
Mattijs Janssens's avatar
Mattijs Janssens committed
228

229
230
231
232
233
234
235
236
            fvMesh& subMesh();

            //- Return point map
            const labelList& pointMap() const;

            //- Return face map
            const labelList& faceMap() const;

237
238
239
            //- Return face map with sign to encode flipped faces
            const labelList& faceFlipMap() const;

240
241
242
243
244
245
246
247
248
249
250
            //- Return cell map
            const labelList& cellMap() const;

            //- Return patch map
            const labelList& patchMap() const;


        // Field mapping

            //- Map volume field
            template<class Type>
251
            static tmp<GeometricField<Type, fvPatchField, volMesh>>
252
253
254
255
256
257
258
259
260
261
            interpolate
            (
                const GeometricField<Type, fvPatchField, volMesh>&,
                const fvMesh& sMesh,
                const labelList& patchMap,
                const labelList& cellMap,
                const labelList& faceMap
            );

            template<class Type>
262
            tmp<GeometricField<Type, fvPatchField, volMesh>>
263
264
265
266
267
            interpolate
            (
                const GeometricField<Type, fvPatchField, volMesh>&
            ) const;

268
269
            //- Map surface field. Optionally negates value if flipping
            //  a face (from exposing an internal face)
270
            template<class Type>
271
            static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
Mattijs Janssens's avatar
Mattijs Janssens committed
272
            interpolate
273
274
275
276
            (
                const GeometricField<Type, fvsPatchField, surfaceMesh>&,
                const fvMesh& sMesh,
                const labelList& patchMap,
277
278
279
                const labelList& cellMap,
                const labelList& faceMap,
                const bool negateIfFlipped = true
280
281
282
            );

            template<class Type>
283
            tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
284
285
            interpolate
            (
286
287
                const GeometricField<Type, fvsPatchField, surfaceMesh>&,
                const bool negateIfFlipped = true
288
            ) const;
Mattijs Janssens's avatar
Mattijs Janssens committed
289
290
291

            //- Map point field
            template<class Type>
292
            static tmp<GeometricField<Type, pointPatchField, pointMesh>>
Mattijs Janssens's avatar
Mattijs Janssens committed
293
294
295
296
297
298
299
300
301
            interpolate
            (
                const GeometricField<Type, pointPatchField, pointMesh>&,
                const pointMesh& sMesh,
                const labelList& patchMap,
                const labelList& pointMap
            );

            template<class Type>
302
            tmp<GeometricField<Type, pointPatchField, pointMesh>>
Mattijs Janssens's avatar
Mattijs Janssens committed
303
304
305
306
            interpolate
            (
                const GeometricField<Type, pointPatchField, pointMesh>&
            ) const;
307

308
            //- Map dimensioned field
309
            template<class Type>
310
            static tmp<DimensionedField<Type, volMesh>>
311
312
313
314
315
316
317
318
            interpolate
            (
                const DimensionedField<Type, volMesh>&,
                const fvMesh& sMesh,
                const labelList& cellMap
            );

            template<class Type>
319
            tmp<DimensionedField<Type, volMesh>>
320
            interpolate(const DimensionedField<Type, volMesh>&) const;
321
322
323
324
325
326
327
328
329
330
};


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

} // End namespace Foam

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#ifdef NoRepository
331
    #include "fvMeshSubsetInterpolate.C"
332
333
334
335
336
337
338
#endif

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#endif

// ************************************************************************* //