Loading...
Searching...
No Matches
zoneBlendedTemplates.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) 2024 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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29
30template<class Type>
31void Foam::zoneBlended<Type>::setSchemes
32(
33 const fvMesh& mesh,
34 const dictionary& dict
35)
36{
37 zoneNames_.resize(dict.size());
38 schemePtrs_.resize(dict.size());
39
40 schemePtrs_[0] =
41 surfaceInterpolationScheme<Type>::New(mesh, dict.lookup("default"));
42
43 zoneNames_[0] = "default";
44
45 corrected_ = scheme(0).corrected();
46
47 label schemei = 1;
48
49 for (const auto& e : dict)
50 {
51 if (e.isDict())
52 {
54 << "Entries must be given in zoneName-scheme pairs"
55 << abort(FatalIOError);
56 }
57
58 const word& key = e.keyword();
59
60 if (key == "default")
61 {
62 // Background scheme - already handled at index 0
63 }
64 else
65 {
66 zoneNames_[schemei] = key;
67
68 schemePtrs_[schemei] =
69 surfaceInterpolationScheme<Type>::New(mesh, e.stream());
70
71 corrected_ = corrected_ || scheme(schemei).corrected();
72
73 ++schemei;
74 }
75 }
76}
77
78
79template<class Type>
80void Foam::zoneBlended<Type>::setSchemes
81(
82 const fvMesh& mesh,
83 const surfaceScalarField& faceFlux,
84 const dictionary& dict
85)
86{
87 zoneNames_.resize(dict.size());
88 schemePtrs_.resize(dict.size());
89
90 schemePtrs_[0] =
91 surfaceInterpolationScheme<Type>::New
92 (
93 mesh,
94 faceFlux,
95 dict.lookup("default")
96 );
97
98 zoneNames_[0] = "default";
99
100 corrected_ = scheme(0).corrected();
101
102 label schemei = 1;
103
104 for (const auto& e : dict)
105 {
106 if (e.isDict())
107 {
109 << "Entries must be given in faceZoneName-scheme pairs"
110 << abort(FatalIOError);
111 }
112
113 const word& key = e.keyword();
114
115 if (key == "default")
116 {
117 // Background scheme - already handled at index 0
118 }
119 else
120 {
121 zoneNames_[schemei] = key;
122
123 schemePtrs_[schemei] =
124 surfaceInterpolationScheme<Type>::New
125 (
126 mesh,
127 faceFlux,
128 e.stream()
129 );
130
131 corrected_ = corrected_ || scheme(schemei).corrected();
132
133 ++schemei;
134 }
135 }
136}
137
138
139template<class Type>
140template<class FieldType>
141void Foam::zoneBlended<Type>::setFaceZoneValues
142(
143 FieldType& dest,
144 const FieldType& src,
145 const faceZone& fz
146) const
147{
148 const auto& mesh = dest.mesh();
149 const auto& pbm = mesh.boundaryMesh();
150 const auto& srcBf = src.boundaryField();
151 auto& destBf = dest.boundaryFieldRef();
152
153 for (const label facei : fz)
154 {
155 if (mesh.isInternalFace(facei))
156 {
157 dest[facei] = src[facei];
158 }
159 else
160 {
161 const labelPair pf = pbm.whichPatchFace(facei);
162 auto& pdest = destBf[pf.first()];
163 if (pdest.size())
164 {
165 pdest[pf.second()] = srcBf[pf.first()][pf.second()];
166 }
167 }
168 }
169}
170
171
172template<class Type>
173template<class FieldType>
174void Foam::zoneBlended<Type>::zeroFaceZoneValues
175(
176 FieldType& dest,
177 const faceZone& fz
178) const
179{
180 const auto& mesh = dest.mesh();
181 const auto& pbm = mesh.boundaryMesh();
182 auto& destBf = dest.boundaryFieldRef();
183
184 for (const label facei : fz)
185 {
186 if (mesh.isInternalFace(facei))
187 {
188 dest[facei] = pTraits<Type>::zero;
189 }
190 else
191 {
192 const labelPair pf = pbm.whichPatchFace(facei);
193 auto& pdest = destBf[pf.first()];
194 if (pdest.size())
195 {
196 pdest[pf.second()] = pTraits<Type>::zero;
197 }
198 }
199 }
200}
201
202
203// ************************************************************************* //
const polyBoundaryMesh & pbm
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
errorManip< error > abort(error &err)
Definition errorManip.H:139
dictionary dict
volScalarField & e