1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
28 Post-processing mesh subset tool. Given the original mesh and the
29 list of selected cells, it creates the mesh consisting only of the
30 desired cells, with the mapping list for points, faces, and cells.
32 Puts all exposed internal faces into either
33 - a user supplied patch
34 - a newly created patch "oldInternalFaces"
36 - setCellSubset is for small subsets. Uses Maps to minimize memory.
37 - setLargeCellSubset is for largish subsets (>10% of mesh).
38 Uses labelLists instead.
40 - setLargeCellSubset does coupled patch subsetting as well. If it detects
41 a face on a coupled patch 'losing' its neighbour it will move the
42 face into the oldInternalFaces patch.
44 - if a user supplied patch is used the mapping becomes a problem.
45 Do the new faces get the value of the internal face they came from?
46 What if e.g. the user supplied patch is a fixedValue 0? So for now
47 they get the face of existing patch face 0.
52 \*---------------------------------------------------------------------------*/
54 #ifndef fvMeshSubset_H
55 #define fvMeshSubset_H
58 #include "pointMesh.H"
59 #include "fvPatchFieldMapper.H"
60 #include "PointPatchFieldMapper.H"
61 #include "GeometricField.H"
63 #include "surfaceMesh.H"
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 /*---------------------------------------------------------------------------*\
71 Class fvMeshSubset Declaration
72 \*---------------------------------------------------------------------------*/
81 //- Patch-field subset interpolation class
82 class patchFieldSubset
84 public fvPatchFieldMapper
86 label sizeBeforeMapping_;
88 labelField directAddressing_;
94 //- Construct from components
95 patchFieldSubset(const label sbm, const labelList& da)
97 sizeBeforeMapping_(sbm),
101 //- Construct given addressing
102 patchFieldSubset(const fvMeshSubset& ms, const label patchNo)
106 ms.baseMesh().boundary()[ms.patchMap()[patchNo]].size()
110 static_cast<const labelField&>
115 ms.subMesh().boundary()[patchNo].size(),
116 ms.subMesh().boundary()[patchNo].patch().start()
119 - ms.baseMesh().boundary()
120 [ms.patchMap()[patchNo]].patch().start()
123 // If patchNo supplied so exposed internal faces or uncoupled
124 // patch faces get into existing patch what to do with
125 // mapping? Truncate to 0 for now.
128 ms.baseMesh().boundary()[ms.patchMap()[patchNo]].size();
130 forAll(directAddressing_, i)
132 label& address = directAddressing_[i];
134 if (address < 0 || address >= sz)
144 virtual ~patchFieldSubset()
152 return directAddressing_.size();
155 label sizeBeforeMapping() const
157 return directAddressing_.size();
165 const unallocLabelList& directAddressing() const
167 return directAddressing_;
172 //- Patch-field subset interpolation class
173 class pointPatchFieldSubset
175 public PointPatchFieldMapper
178 const labelList& directAddressing_;
185 //- Construct given addressing
186 pointPatchFieldSubset(const labelList& directAddressing)
188 directAddressing_(directAddressing)
193 virtual ~pointPatchFieldSubset()
201 return directAddressing_.size();
204 label sizeBeforeMapping() const
206 return directAddressing_.size();
214 const unallocLabelList& directAddressing() const
216 return directAddressing_;
225 //- Reference to original mesh
226 const fvMesh& baseMesh_;
228 //- Subset fvMesh pointer
229 fvMesh* fvMeshSubsetPtr_;
231 //- Subset pointMesh pointer
232 mutable pointMesh* pointMeshSubsetPtr_;
234 //- Point mapping array
237 //- Face mapping array
240 //- Cell mapping array
243 //- Patch mapping array
247 // Private Member Functions
249 //- Disallow default bitwise copy construct
250 fvMeshSubset(const fvMeshSubset&);
252 //- Disallow default bitwise assignment
253 void operator=(const fvMeshSubset&);
256 //- Check if subset has been performed
257 bool checkCellSubset() const;
259 //- Mark points in Map
260 static void markPoints(const labelList&, Map<label>&);
262 //- Mark points (with 0) in labelList
263 static void markPoints(const labelList&, labelList&);
265 //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
266 void doCoupledPatches
269 labelList& nCellsUsingFace
273 static labelList subset
276 const labelList& selectedElements,
277 const labelList& subsetMap
280 //- Create zones for submesh
283 //- Make dictionaries in case they are not present
284 // Creates fvSchemes and fvSolution copying from base
285 void makeFvDictionaries();
292 //- Construct given a mesh
293 explicit fvMeshSubset(const IOobject&, const fvMesh&);
305 //- Set the subset. Create "oldInternalFaces" patch for exposed
306 // internal faces (patchID==-1) or use supplied patch.
307 // Does not handle coupled patches correctly if only one side
311 const labelHashSet& globalCellMap,
312 const label patchID = -1
315 //- Set the subset from all cells with region == currentRegion.
316 // Create "oldInternalFaces" patch for exposed
317 // internal faces (patchID==-1) or use supplied patch.
318 // Handles coupled patches by if nessecary making coupled patch
319 // face part of patchID (so uncoupled)
320 void setLargeCellSubset
322 const labelList& region,
323 const label currentRegion,
324 const label patchID = -1,
325 const bool syncCouples = true
328 //- setLargeCellSubset but with labelHashSet.
329 void setLargeCellSubset
331 const labelHashSet& globalCellMap,
332 const label patchID = -1,
333 const bool syncPar = true
339 //- Return reference to base mesh
340 const fvMesh& baseMesh() const
345 //- Return reference to subset mesh
346 const fvMesh& subMesh() const;
350 //- Return reference to demand-driven subset pointMesh
351 const pointMesh& subPointMesh() const;
353 pointMesh& subPointMesh();
356 const labelList& pointMap() const;
359 const labelList& faceMap() const;
362 const labelList& cellMap() const;
365 const labelList& patchMap() const;
370 //- Map volume field given map addressing
371 // - it does not patch up the exposed patches with surface fields
372 // - it does not handle patch type creation correctly
373 // Please use only in dummy mapping operations
376 static tmp<GeometricField<Type, fvPatchField, volMesh> >
379 const GeometricField<Type, fvPatchField, volMesh>&,
381 const labelList& patchMap,
382 const labelList& cellMap,
383 const labelList& faceMap
388 tmp<GeometricField<Type, fvPatchField, volMesh> >
391 const GeometricField<Type, fvPatchField, volMesh>&
394 //- Map surface field
396 static tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
399 const GeometricField<Type, fvsPatchField, surfaceMesh>&,
401 const labelList& patchMap,
402 const labelList& faceMap
406 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
409 const GeometricField<Type, fvsPatchField, surfaceMesh>&
414 static tmp<GeometricField<Type, pointPatchField, pointMesh> >
417 const GeometricField<Type, pointPatchField, pointMesh>&,
418 const pointMesh& sMesh,
419 const labelList& patchMap,
420 const labelList& pointMap
424 tmp<GeometricField<Type, pointPatchField, pointMesh> >
427 const GeometricField<Type, pointPatchField, pointMesh>&
432 virtual bool writeData(Foam::Ostream&) const
440 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
442 } // End namespace Foam
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
447 # include "fvMeshSubsetInterpolate.C"
450 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
454 // ************************************************************************* //