surfaceIteratorIso.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) 2019-2020 DLR
9 -------------------------------------------------------------------------------
10 License
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 "surfaceIteratorIso.H"
29 #include "scalarMatrices.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const fvMesh& mesh,
36  scalarField& pointVal,
37  const scalar tol
38 )
39 :
40  mesh_(mesh),
41  ap_(pointVal),
42  cutCell_(mesh_,ap_),
43  surfCellTol_(tol)
44 {}
45 
46 
47 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48 
50 ( const label celli,
51  const scalar alpha1,
52  const scalar tol,
53  const label maxIter
54 )
55 {
56  // Finding cell vertex extrema values
57  const labelList& pLabels = mesh_.cellPoints()[celli];
58  scalarField fvert(pLabels.size());
59  forAll(pLabels, pi)
60  {
61  fvert[pi] = ap_[pLabels[pi]];
62  }
63 
64  const labelList order(Foam::sortedOrder(fvert));
65  scalar f1 = fvert[order.first()];
66  scalar f2 = fvert[order.last()];
67 
68  // Special case where method is given an almost full or empty cell
69  if (alpha1 < tol)
70  {
71  return -1; // is area and centre is zero;
72  }
73  else if (1 - alpha1 < tol)
74  {
75  return 1;
76  }
77 
78  // Finding the two vertices inbetween which the isovalue giving alpha1 lies
79  label L1 = 0;
80  label L2 = fvert.size() - 1;
81  scalar a1 = 1;
82  scalar a2 = 0;
83  scalar L3, f3, a3;
84 
85  while (L2 - L1 > 1)
86  {
87  L3 = round(0.5*(L1 + L2));
88  f3 = fvert[order[L3]];
89  cutCell_.calcSubCell(celli, f3);
90  a3 = cutCell_.VolumeOfFluid();
91  if (a3 > alpha1)
92  {
93  L1 = L3; f1 = f3; a1 = a3;
94  }
95  else if (a3 < alpha1)
96  {
97  L2 = L3; f2 = f3; a2 = a3;
98  }
99  else
100  {
101  return cutCell_.calcSubCell(celli, f3);
102  }
103  }
104 
105  if (mag(f1 - f2) < 10*SMALL)
106  {
107  return cutCell_.calcSubCell(celli, f1);
108  }
109 
110  if (mag(a1 - a2) < tol)
111  {
112  return cutCell_.calcSubCell(celli, 0.5*(f1 + f2));
113  }
114  // Now we know that a(f) = alpha1 is to be found on the f interval
115 
116 
117  // Finding coefficients in 3 deg polynomial alpha(f) from 4 solutions
118  // Finding 2 additional points on 3 deg polynomial
119  f3 = f1 + (f2 - f1)/3;
120  cutCell_.calcSubCell(celli, f3);
121  a3 = cutCell_.VolumeOfFluid();
122 
123  scalar f4 = f1 + 2*(f2 - f1)/3;
124  cutCell_.calcSubCell(celli, f4);
125  scalar a4 = cutCell_.VolumeOfFluid();
126 
127  // Building and solving Vandermonde matrix equation
128  scalarField a(4), f(4), C(4), fOld(4);
129  {
130  a[0] = a1, a[1] = a3, a[2] = a4, a[3] = a2;
131  fOld[0] = f1, fOld[1] = f3, fOld[2] = f4, fOld[3] = f2;
132  f[0] = 0, f[1] = (f3-f1)/(f2-f1), f[2] = (f4-f1)/(f2-f1), f[3] = 1;
134  forAll(f, i)
135  {
136  forAll(f, j)
137  {
138  M[i][j] = pow(f[i], 3 - j);
139  }
140  }
141  // C holds the 4 polynomial coefficients
142  C = a;
143  LUsolve(M, C);
144  }
145 
146  // Finding root with Newton method
147 
148  f3 = f[1]; a3 = a[1];
149  label nIter = 0;
150  scalar res = mag(a3 - alpha1);
151  while (res > tol && nIter < 10*maxIter)
152  {
153  f3 -= (C[0]*pow3(f3) + C[1]*sqr(f3) + C[2]*f3 + C[3] - alpha1)
154  /(3*C[0]*sqr(f3) + 2*C[1]*f3 + C[2]);
155  a3 = C[0]*pow3(f3) + C[1]*sqr(f3) + C[2]*f3 + C[3];
156  res = mag(a3 - alpha1);
157  nIter++;
158  }
159  // Scaling back to original range
160  f3 = f3*(f2 - f1) + f1;
161 
162  //Check result
163  label status = cutCell_.calcSubCell(celli, f3);
164  const scalar VOF = cutCell_.VolumeOfFluid();
165  res = mag(VOF - alpha1);
166 
167  if (res > tol)
168  {
169  }
170  else
171  {
172  return status;
173  }
174 
175  // If tolerance not met use the secant method with f3 as a hopefully very
176  // good initial guess to crank res the last piece down below tol
177 
178  scalar x2 = f3;
179  scalar g2 = VOF - alpha1;
180  scalar x1 = max(1e-3*(f2 - f1), 100*SMALL);
181  x1 = max(x1, f1);
182  x1 = min(x1, f2);
183  cutCell_.calcSubCell(celli, x1);
184  scalar g1 = cutCell_.VolumeOfFluid() - alpha1;
185 
186  nIter = 0;
187  scalar g0(0), x0(0);
188  while (res > tol && nIter < maxIter && g1 != g2)
189  {
190  x0 = (x2*g1 - x1*g2)/(g1 - g2);
191  status = cutCell_.calcSubCell(celli, x0);
192  g0 = cutCell_.VolumeOfFluid() - alpha1;
193  res = mag(g0);
194  x2 = x1; g2 = g1;
195  x1 = x0; g1 = g0;
196  nIter++;
197  }
198 
199  return status;
200 }
201 
202 
203 // ************************************************************************* //
Foam::LUsolve
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Definition: scalarMatricesTemplates.C:215
Foam::surfaceIteratorIso::surfaceIteratorIso
surfaceIteratorIso(const fvMesh &mesh, scalarField &pointVal, const scalar tol)
Construct from fvMesh and a scalarField.
Definition: surfaceIteratorIso.C:34
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
C
volScalarField & C
Definition: readThermalProperties.H:102
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::SquareMatrix< scalar >
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
f
labelList f(nPoints)
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
M
#define M(I)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
scalarMatrices.H
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::surfaceIteratorIso::vofCutCell
label vofCutCell(const label celli, const scalar alpha1, const scalar tol, const label maxIter)
Definition: surfaceIteratorIso.C:50
surfaceIteratorIso.H