Loading...
Searching...
No Matches
ITHACAsurfacetools.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26
27 You should have received a copy of the GNU Lesser General Public License
28 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
29\*---------------------------------------------------------------------------*/
30
31#include "ITHACAsurfacetools.H"
32
33
34namespace ITHACAutilities
35{
37{
38
39template<typename T>
40List<label> surfaceIndexInt(T& field, const label patchInt,
41 const label patchExt)
42{
43 List<label>* result = new List<label>;
44
45 for (size_t i = 0; i < field.mesh().boundaryMesh()[patchInt].size(); i++)
46 {
47 result->append(field.mesh().boundaryMesh()[patchInt].faceCells()[i]);
48 }
49
50 return *result;
51}
52
53template List<label> surfaceIndexInt(volScalarField& field,
54 const label patchInt, const label patchExt);
55template List<label> surfaceIndexInt(volVectorField& field,
56 const label patchInt, const label patchExt);
57template List<label> surfaceIndexInt(volTensorField& field,
58 const label patchInt, const label patchExt);
59
60template<typename T>
61List<label> surfaceIndexExt(T& field, const label patchInt,
62 const label patchExt)
63{
64 return surfaceIndexInt(field, patchExt, patchInt);
65}
66
67template List<label> surfaceIndexExt(volScalarField& field,
68 const label patchInt, const label patchExt);
69template List<label> surfaceIndexExt(volVectorField& field,
70 const label patchInt, const label patchExt);
71template List<label> surfaceIndexExt(volTensorField& field,
72 const label patchInt, const label patchExt);
73
74template<typename T, typename V>
75void surfaceValuesInt(T& field, const label patchInt, const label patchExt,
76 List<V>& result)
77{
78 List<label> indexes = surfaceIndexInt(field, patchInt, patchExt);
79
80 for (size_t i = 0; i < indexes.size(); i++)
81 {
82 result.append(field[indexes[i]]);
83 }
84}
85
86template void surfaceValuesInt(volScalarField& field, const label patchInt,
87 const label patchExt, List<scalar>& result);
88template void surfaceValuesInt(volVectorField& field, const label patchInt,
89 const label patchExt, List<Foam::Vector<scalar >> & result);
90template void surfaceValuesInt(volTensorField& field, const label patchInt,
91 const label patchExt, List<Foam::Tensor<scalar >>& result);
92
93template<typename T, typename V>
94void surfaceValuesExt(T& field, const label patchInt, const label patchExt,
95 List<V>& result)
96{
97 surfaceValuesInt(field, patchExt, patchInt, result);
98}
99
100template void surfaceValuesExt(volScalarField& field, const label patchInt,
101 const label patchExt, List<scalar>& result);
102template void surfaceValuesExt(volVectorField& field, const label patchInt,
103 const label patchExt, List<Foam::Vector<scalar >> & result);
104template void surfaceValuesExt(volTensorField& field, const label patchInt,
105 const label patchExt, List<Foam::Tensor<scalar >>& result);
106
107template<typename T>
108Foam::Vector<scalar> surfaceFindMirrorPoint(T& field, const label patchInt,
109 const label patchExt, const label cellID)
110{
111 Foam::Vector<scalar> result = * (new Foam::Vector<scalar>);
112 List<label> indexesInt = surfaceIndexInt(field, patchInt, patchExt);
113 result = 2.0 * field.mesh().boundaryMesh()[patchInt].faceCentres()[cellID] -
114 field.mesh().C()[indexesInt[cellID]];
115 return result;
116}
117
118template Foam::Vector<scalar> surfaceFindMirrorPoint(volScalarField& field,
119 const label patchInt, const label patchExt, const label cellID);
120template Foam::Vector<scalar> surfaceFindMirrorPoint(volVectorField& field,
121 const label patchInt, const label patchExt, const label cellID);
122template Foam::Vector<scalar> surfaceFindMirrorPoint(volTensorField& field,
123 const label patchInt, const label patchExt, const label cellID);
124
125template<typename T>
126label surfaceFindClosest(T& field, const label patchInt, const label patchExt,
127 Foam::Vector<scalar> point)
128{
129 label result = 0;
130 scalar dist = 0;
131 List<label> indexesExt = surfaceIndexExt(field, patchInt, patchExt);
132 dist = mag(point - field.mesh().C()[indexesExt[0]]);
133
134 for (int i = 1; i < indexesExt.size(); i++)
135 {
136 scalar temp = mag(point - field.mesh().C()[indexesExt[i]]);
137
138 if (temp <= dist)
139 {
140 dist = temp;
141 result = i;
142 }
143 }
144
145 result = indexesExt[result];
146 return result;
147}
148
149template label surfaceFindClosest(volScalarField& field, const label patchInt,
150 const label patchExt, Foam::Vector<scalar> point);
151template label surfaceFindClosest(volVectorField& field, const label patchInt,
152 const label patchExt, Foam::Vector<scalar> point);
153template label surfaceFindClosest(volTensorField& field, const label patchInt,
154 const label patchExt, Foam::Vector<scalar> point);
155
156template<typename T, typename V>
157void surfaceAverage(T& field, const label patchInt, const label patchExt,
158 List<V>& result)
159{
160 List<label> indexesInt = surfaceIndexInt(field, patchInt, patchExt);
161 result.resize(0);
162
163 for (int i = 0; i < indexesInt.size(); i++)
164 {
165 Foam::Vector<scalar> mirror = surfaceFindMirrorPoint(field, patchInt, patchExt,
166 i);
167 label closest = surfaceFindClosest(field, patchInt, patchExt, mirror);
168 result.append(0.5 * field[closest] + 0.5 * field[indexesInt[i]]);
169 }
170}
171
172template void surfaceAverage(volScalarField& field, const label patchInt,
173 const label patchExt, List<scalar>& result);
174template void surfaceAverage(volVectorField& field, const label patchInt,
175 const label patchExt, List<Foam::Vector<scalar >> & result);
176template void surfaceAverage(volTensorField& field, const label patchInt,
177 const label patchExt, List<Foam::Tensor<scalar >>& result);
178
179template<typename T, typename V>
180void surfaceJump(T& field, const label patchInt, const label patchExt,
181 List<V>& result)
182{
183 List<label> indexesInt = surfaceIndexInt(field, patchInt, patchExt);
184 result.resize(0);
185
186 for (int i = 0; i < indexesInt.size(); i++)
187 {
188 Foam::Vector<scalar> mirror = surfaceFindMirrorPoint(field, patchInt, patchExt,
189 i);
190 label closest = surfaceFindClosest(field, patchInt, patchExt, mirror);
191 result.append(field[closest] - field[indexesInt[i]]);
192 }
193}
194
195template void surfaceJump(volScalarField& field, const label patchInt,
196 const label patchExt, List<scalar>& result);
197template void surfaceJump(volVectorField& field, const label patchInt,
198 const label patchExt, List<Foam::Vector<scalar >> & result);
199template void surfaceJump(volTensorField& field, const label patchInt,
200 const label patchExt, List<Foam::Tensor<scalar >>& result);
201} // End namespace ITHACAsurfacetools
202} // End namespace ITHACAutilities
203
204// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205
206// ************************************************************************* //
volScalarField & T
Definition createT.H:46
void surfaceValuesInt(T &field, const label patchInt, const label patchExt, List< V > &result)
List< label > surfaceIndexInt(T &field, const label patchInt, const label patchExt)
Foam::Vector< scalar > surfaceFindMirrorPoint(T &field, const label patchInt, const label patchExt, const label cellID)
label surfaceFindClosest(T &field, const label patchInt, const label patchExt, Foam::Vector< scalar > point)
void surfaceAverage(T &field, const label patchInt, const label patchExt, List< V > &result)
void surfaceValuesExt(T &field, const label patchInt, const label patchExt, List< V > &result)
void surfaceJump(T &field, const label patchInt, const label patchExt, List< V > &result)
List< label > surfaceIndexExt(T &field, const label patchInt, const label patchExt)
Namespace to implement some useful assign operation of OF fields.
label i
Definition pEqn.H:46