Loading...
Searching...
No Matches
dynamicRefineFvMeshTemplates.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2018-2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "surfaceFields.H"
29
30// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31
32template<class T>
34(
35 const labelList& faceMap,
37)
38{
40
41 //- Make flat field for ease of looping
42 Field<T> tsFld(this->nFaces(), Zero);
43 SubField<T>(tsFld, this->nInternalFaces()) = sFld.internalField();
44
45 const typename GeoField::Boundary& bFld = sFld.boundaryField();
46 forAll(bFld, patchi)
47 {
48 label facei = this->boundaryMesh()[patchi].start();
49 for (const T& val : bFld[patchi])
50 {
51 tsFld[facei++] = val;
52 }
53 }
54
55 const labelUList& owner = this->faceOwner();
56 const labelUList& neighbour = this->faceNeighbour();
57 const cellList& cells = this->cells();
58
59 for (label facei = 0; facei < nInternalFaces(); facei++)
60 {
61 label oldFacei = faceMap[facei];
62
63 // Map surface field on newly generated faces by obtaining the
64 // hull of the outside faces
65 if (oldFacei == -1)
66 {
67 // Loop over all owner/neighbour cell faces
68 // and find already mapped ones (master-faces):
69 T tmpValue(pTraits<T>::zero);
70 label counter = 0;
71
72 const cell& ownFaces = cells[owner[facei]];
73 for (auto ownFacei : ownFaces)
74 {
75 if (faceMap[ownFacei] != -1)
76 {
77 tmpValue += tsFld[ownFacei];
78 counter++;
79 }
80 }
81
82 const cell& neiFaces = cells[neighbour[facei]];
83 for (auto neiFacei : neiFaces)
84 {
85 if (faceMap[neiFacei] != -1)
86 {
87 tmpValue += tsFld[neiFacei];
88 counter++;
89 }
90 }
91
92 if (counter > 0)
93 {
94 sFld[facei] = tmpValue/counter;
95 }
96 }
97 }
98}
99
100
101template<class T>
103(
104 const labelList& faceMap
105)
106{
107 typedef GeometricField<T, fvsPatchField, surfaceMesh> GeoField;
108
109 for (GeoField& sFld : this->objectRegistry::sorted<GeoField>())
110 {
111 //if (mapSurfaceFields_.found(sFld.name()))
112 {
114 << "dynamicRefineFvMesh::mapNewInternalFaces():"
115 << " Mapping new internal faces by interpolation on "
116 << sFld.name()<< endl;
117
118 if (sFld.is_oriented())
119 {
120 WarningInFunction << "Ignoring mapping oriented field "
121 << sFld.name() << " since of type " << sFld.type() << endl;
122 }
123 else
124 {
125 mapNewInternalFaces(faceMap, sFld);
126 }
128 }
129}
130
131template<class T>
133(
134 const surfaceVectorField& Sf,
135 const surfaceScalarField& magSf,
136 const labelList& faceMap
137)
138{
140
141 typedef GeometricField
142 <
146 > NormalGeoField;
147
148
149 for (GeoField& sFld : this->objectRegistry::sorted<GeoField>())
150 {
151 //if (mapSurfaceFields_.found(sFld.name()))
152 {
154 << "dynamicRefineFvMesh::mapNewInternalFaces():"
155 << " Mapping new internal faces by interpolation on "
156 << sFld.name() << endl;
157
158 if (sFld.is_oriented())
159 {
161 << "dynamicRefineFvMesh::mapNewInternalFaces(): "
162 << "Converting oriented field " << sFld.name()
163 << " to intensive field and mapping" << endl;
164
165 // Assume any oriented field is face area weighted (i.e. a flux)
166 // Convert to intensive (& oriented) before mapping. Untested.
167
168 // Convert to intensive and non oriented
169 NormalGeoField fFld(sFld*Sf/Foam::sqr(magSf));
170
171 // Interpolate
172 mapNewInternalFaces(faceMap, fFld);
173
174 // Convert back to extensive and oriented
175 sFld = (fFld & Sf);
176 }
177 else
178 {
179 mapNewInternalFaces(faceMap, sFld);
180 }
181 }
182 }
183}
184
185
186// ************************************************************************* //
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition SubField.H:58
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
void mapNewInternalFaces(const labelList &faceMap, GeometricField< T, fvsPatchField, surfaceMesh > &)
Map single non-flux surface<Type>Field.
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition fvMesh.H:572
const labelUList & neighbour() const
Internal face neighbour.
Definition fvMesh.H:580
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
UPtrList< Type > sorted()
Return sorted list of objects with a class satisfying isA<Type> or isType<Type> (with Strict).
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
label nInternalFaces() const noexcept
Number of internal faces.
label nFaces() const noexcept
Number of mesh faces.
Mesh data needed to do the Finite Volume discretisation.
Definition surfaceMesh.H:47
const cellShapeList & cells
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.