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{
36namespace ITHACAsurfacetools
37{
38
39 template<typename T>
40 List<label> surfaceIndexInt(T& field, const label patchInt, const label patchExt)
41 {
42 List<label>* result = new List<label>;
43
44 for (size_t i = 0; i < field.mesh().boundaryMesh()[patchInt].size(); i++) {
45 result->append(field.mesh().boundaryMesh()[patchInt].faceCells()[i]);
46 }
47
48 return *result;
49 }
50
51 template List<label> surfaceIndexInt(volScalarField& field, const label patchInt, const label patchExt);
52 template List<label> surfaceIndexInt(volVectorField& field, const label patchInt, const label patchExt);
53 template List<label> surfaceIndexInt(volTensorField& field, const label patchInt, const label patchExt);
54
55 template<typename T>
56 List<label> surfaceIndexExt(T& field, const label patchInt, const label patchExt)
57 {
58 return surfaceIndexInt(field, patchExt, patchInt);
59 }
60
61 template List<label> surfaceIndexExt(volScalarField& field, const label patchInt, const label patchExt);
62 template List<label> surfaceIndexExt(volVectorField& field, const label patchInt, const label patchExt);
63 template List<label> surfaceIndexExt(volTensorField& field, const label patchInt, const label patchExt);
64
65 template<typename T, typename V>
66 void surfaceValuesInt(T& field, const label patchInt, const label patchExt, List<V>& result)
67 {
68 List<label> indexes = surfaceIndexInt(field, patchInt, patchExt);
69
70 for (size_t i = 0; i < indexes.size(); i++) {
71 result.append(field[indexes[i]]);
72 }
73 }
74
75 template void surfaceValuesInt(volScalarField& field, const label patchInt, const label patchExt, List<scalar>& result);
76 template void surfaceValuesInt(volVectorField& field, const label patchInt, const label patchExt, List<Foam::Vector<scalar>>& result);
77 template void surfaceValuesInt(volTensorField& field, const label patchInt, const label patchExt, List<Foam::Tensor<scalar>>& result);
78
79 template<typename T, typename V>
80 void surfaceValuesExt(T& field, const label patchInt, const label patchExt, List<V>& result)
81 {
82 surfaceValuesInt(field, patchExt, patchInt, result);
83 }
84
85 template void surfaceValuesExt(volScalarField& field, const label patchInt, const label patchExt, List<scalar>& result);
86 template void surfaceValuesExt(volVectorField& field, const label patchInt, const label patchExt, List<Foam::Vector<scalar>>& result);
87 template void surfaceValuesExt(volTensorField& field, const label patchInt, const label patchExt, List<Foam::Tensor<scalar>>& result);
88
89 template<typename T>
90 Foam::Vector<scalar> surfaceFindMirrorPoint(T& field, const label patchInt, const label patchExt, const label cellID)
91 {
92 Foam::Vector<scalar> result = *(new Foam::Vector<scalar>);
93 List<label> indexesInt = surfaceIndexInt(field, patchInt, patchExt);
94
95 result = 2.0 * field.mesh().boundaryMesh()[patchInt].faceCentres()[cellID] - field.mesh().C()[indexesInt[cellID]];
96
97 return result;
98 }
99
100 template Foam::Vector<scalar> surfaceFindMirrorPoint(volScalarField& field, const label patchInt, const label patchExt, const label cellID);
101 template Foam::Vector<scalar> surfaceFindMirrorPoint(volVectorField& field, const label patchInt, const label patchExt, const label cellID);
102 template Foam::Vector<scalar> surfaceFindMirrorPoint(volTensorField& field, const label patchInt, const label patchExt, const label cellID);
103
104 template<typename T>
105 label surfaceFindClosest(T& field, const label patchInt, const label patchExt, Foam::Vector<scalar> point)
106 {
107 label result = 0;
108 scalar dist = 0;
109 List<label> indexesExt = surfaceIndexExt(field, patchInt, patchExt);
110
111 dist = mag(point - field.mesh().C()[indexesExt[0]]);
112
113 for(int i = 1; i < indexesExt.size(); i++)
114 {
115 scalar temp = mag(point - field.mesh().C()[indexesExt[i]]);
116 if(temp <= dist)
117 {
118 dist = temp;
119 result = i;
120 }
121 }
122 result = indexesExt[result];
123 return result;
124 }
125
126 template label surfaceFindClosest(volScalarField& field, const label patchInt, const label patchExt, Foam::Vector<scalar> point);
127 template label surfaceFindClosest(volVectorField& field, const label patchInt, const label patchExt, Foam::Vector<scalar> point);
128 template label surfaceFindClosest(volTensorField& field, const label patchInt, const label patchExt, Foam::Vector<scalar> point);
129
130 template<typename T, typename V>
131 void surfaceAverage(T& field, const label patchInt, const label patchExt, List<V>& result)
132 {
133 List<label> indexesInt = surfaceIndexInt(field, patchInt, patchExt);
134 result.resize(0);
135
136 for(int i = 0; i < indexesInt.size(); i++)
137 {
138 Foam::Vector<scalar> mirror = surfaceFindMirrorPoint(field, patchInt, patchExt, i);
139 label closest = surfaceFindClosest(field, patchInt, patchExt, mirror);
140 result.append(0.5 * field[closest] + 0.5 * field[indexesInt[i]]);
141 }
142 }
143
144 template void surfaceAverage(volScalarField& field, const label patchInt, const label patchExt, List<scalar>& result);
145 template void surfaceAverage(volVectorField& field, const label patchInt, const label patchExt, List<Foam::Vector<scalar>>& result);
146 template void surfaceAverage(volTensorField& field, const label patchInt, const label patchExt, List<Foam::Tensor<scalar>>& result);
147
148 template<typename T, typename V>
149 void surfaceJump(T& field, const label patchInt, const label patchExt, List<V>& result)
150 {
151 List<label> indexesInt = surfaceIndexInt(field, patchInt, patchExt);
152 result.resize(0);
153
154 for(int i = 0; i < indexesInt.size(); i++)
155 {
156 Foam::Vector<scalar> mirror = surfaceFindMirrorPoint(field, patchInt, patchExt, i);
157 label closest = surfaceFindClosest(field, patchInt, patchExt, mirror);
158 result.append(field[closest] - field[indexesInt[i]]);
159 }
160 }
161
162 template void surfaceJump(volScalarField& field, const label patchInt, const label patchExt, List<scalar>& result);
163 template void surfaceJump(volVectorField& field, const label patchInt, const label patchExt, List<Foam::Vector<scalar>>& result);
164 template void surfaceJump(volTensorField& field, const label patchInt, const label patchExt, List<Foam::Tensor<scalar>>& result);
165} // End namespace ITHACAsurfacetools
166} // End namespace ITHACAutilities
167
168// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169
170// ************************************************************************* //
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