fvMeshSubset.H 10.2 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2017 OpenFOAM Foundation
6
     \\/     M anipulation  | Copyright (C) 2016-2018 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

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"

    - 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.

43
44
45
46
47
48
    - 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)
49
50
51
52
53
54
55
56
57
58

SourceFiles
    fvMeshSubset.C

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

#ifndef fvMeshSubset_H
#define fvMeshSubset_H

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

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

namespace Foam
{

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

class fvMeshSubset
{
    // Private data

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

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

83
84
85
86
87
88
89
90
91
92
93
94
        //- Point mapping array
        labelList pointMap_;

        //- Face mapping array
        labelList faceMap_;

        //- Cell mapping array
        labelList cellMap_;

        //- Patch mapping array
        labelList patchMap_;

95
        //- Optional face mapping array with flip encoded
96
97
        mutable autoPtr<labelList> faceFlipMapPtr_;

98
99
100
101
102
103
104

    // Private Member Functions

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

        //- Mark points (with 0) in labelList
105
106
107
108
109
        static void markPoints
        (
            const labelUList& curPoints,
            labelList& pointMap
        );
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128

        //- 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();

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

136
137
        //- No copy construct
        fvMeshSubset(const fvMeshSubset&) = delete;
138

139
140
        //- No copy assignment
        void operator=(const fvMeshSubset&) = delete;
141

142

143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
public:

    // Constructors

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


    // Member Functions

        // Edit

            //- 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
158
            //  Handles coupled patches by if necessary making coupled patch
159
160
161
162
163
164
165
166
167
            //  face part of patchID (so uncoupled)
            void setLargeCellSubset
            (
                const labelList& region,
                const label currentRegion,
                const label patchID = -1,
                const bool syncCouples = true
            );

168
169
170
            //- setLargeCellSubset but only marking certain cells
            void setLargeCellSubset
            (
171
                const labelUList& globalCellMap,
172
173
174
175
                const label patchID = -1,
                const bool syncPar = true
            );

176
177
178
179
180
181
182
183
184
            //- setLargeCellSubset but with labelHashSet.
            void setLargeCellSubset
            (
                const labelHashSet& globalCellMap,
                const label patchID = -1,
                const bool syncPar = true
            );


185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
            //- 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
                );


211
212
213
214
215
216
217
218
        // Access

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

mattijs's avatar
mattijs committed
219
220
221
            //- Have subMesh?
            bool hasSubMesh() const;

222
223
            //- Return reference to subset mesh
            const fvMesh& subMesh() const;
Mattijs Janssens's avatar
Mattijs Janssens committed
224

225
226
227
228
229
230
231
232
            fvMesh& subMesh();

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

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

233
234
235
            //- Return face map with sign to encode flipped faces
            const labelList& faceFlipMap() const;

236
237
238
239
240
241
242
243
244
245
246
            //- Return cell map
            const labelList& cellMap() const;

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


        // Field mapping

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

            template<class Type>
258
            tmp<GeometricField<Type, fvPatchField, volMesh>>
259
260
261
262
263
            interpolate
            (
                const GeometricField<Type, fvPatchField, volMesh>&
            ) const;

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

            template<class Type>
278
            tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
279
280
            interpolate
            (
281
                const GeometricField<Type, fvsPatchField, surfaceMesh>&
282
            ) const;
Mattijs Janssens's avatar
Mattijs Janssens committed
283
284
285

            //- Map point field
            template<class Type>
286
            static tmp<GeometricField<Type, pointPatchField, pointMesh>>
Mattijs Janssens's avatar
Mattijs Janssens committed
287
288
289
290
291
292
293
294
295
            interpolate
            (
                const GeometricField<Type, pointPatchField, pointMesh>&,
                const pointMesh& sMesh,
                const labelList& patchMap,
                const labelList& pointMap
            );

            template<class Type>
296
            tmp<GeometricField<Type, pointPatchField, pointMesh>>
Mattijs Janssens's avatar
Mattijs Janssens committed
297
298
299
300
            interpolate
            (
                const GeometricField<Type, pointPatchField, pointMesh>&
            ) const;
301

302
            //- Map dimensioned field
303
            template<class Type>
304
            static tmp<DimensionedField<Type, volMesh>>
305
306
307
308
309
310
311
312
            interpolate
            (
                const DimensionedField<Type, volMesh>&,
                const fvMesh& sMesh,
                const labelList& cellMap
            );

            template<class Type>
313
            tmp<DimensionedField<Type, volMesh>>
314
            interpolate(const DimensionedField<Type, volMesh>&) const;
315
316
317
318
319
320
321
322
323
324
};


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

} // End namespace Foam

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

#ifdef NoRepository
325
    #include "fvMeshSubsetInterpolate.C"
326
327
328
329
330
331
332
#endif

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

#endif

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