-
Notifications
You must be signed in to change notification settings - Fork 11
/
Copy pathmanifold.cc
186 lines (168 loc) · 7.58 KB
/
manifold.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
#include <Sacado.hpp>
#include<cmath>
#include "manifold.h"
using namespace dealii;
//----------------------------------------------------------------------------------------------------
//This is a curve boundary manifold object used in ringleb's flow.
//it's used on the left and right side of the domain, if necessary, can also be used on interior cells
//----------------------------------------------------------------------------------------------------
//given point in the real space, return point in the reference domain
Point<2> RinglebSide::pull_back(const Point<2>& x) const
{
double Vmin=0.2, Vmax=1.8; //just used as initial guess for iteration process to numerically solve V
unsigned int imax;
double b, rho, L, V, psi,Vpsi, Vp, tol, theta;
imax =500;
tol = 1.0e-15;
//Initial guess, Vp=1
Vp = 0.5*( Vmin + Vmax );
//Start iterations. Eq.(7.13.76). given point x, return V
for(unsigned int i=1; i<imax; ++i){
b = std::sqrt(1-0.2*Vp*Vp);
rho = std::pow(b,5);
L = 1./b+1./(3*b*b*b)+1./(5*rho)-0.5*std::log((1+b)/(1-b));
V = std::sqrt(1/2./rho/std::sqrt((x[0]-L/2)*(x[0]-L/2)+x[1]*x[1])); // Eq.(7.11.76)
if(std::fabs(V-Vp)<tol)
break;
if(i == imax){
std::cout<<"did not converge...sorry."<<std::endl;
exit(0);
}
Vp = V;
}
psi = std::sqrt(0.5/V/V-rho*(x[0]-L/2));
Vpsi = std::sqrt(0.5-V*V*(x[0]-L/2)*rho);
if ( std::fabs( Vpsi - 1 ) < 1.0e-15 )
theta = 0.5*3.141592653589793238; // Do not rely on `asin' here.
else
theta = std::asin( psi*V );
/*std::cout<<"given real_space point, (x,y) = "<<x[0]<<","<<x[1]
<<std::endl<<"return reference_space point, (psi, V)= "<<psi<<","<<V
<<std::endl;*/ //used in test.
return Point<2>(psi, theta);
}
//given point in the reference domain, return point in the real space
Point<2> RinglebSide::push_forward(const Point<2>& psi_theta) const
{
const double psi = psi_theta[0];
const double theta = psi_theta[1];
const double V = std::sin(theta)/psi;
const double b = std::sqrt(1-0.2*V*V);
const double rho = std::pow(b,5);
const double L = 1./b+1./(3*b*b*b)+1./(5*rho)-0.5*std::log((1+b)/(1-b));
if(std::fabs(psi*V-1) < 1.0e-15){
return Point<2> ((1./rho)*(0.5/V/V- psi*psi) + 0.5*L, //x, Eq.(7.11.66)
0); //y, Eq.(7.11.67)), to avoid sqrt(negative) below.
}
else{
return Point<2> ((1./rho)*(0.5/V/V- psi*psi) + 0.5*L, //x, Eq.(7.11.66)
psi/(rho*V)*std::sqrt(1.-psi*psi*V*V)); //y, Eq.(7.11.67))
}
}
//----------------------------------------------------------------------------------------------------
// Another Manifold
// this is used at the top boundary of rinleb's flow domain(local coordinate now is psi_V, not psi_theta)
//----------------------------------------------------------------------------------------------------
//given point in the real space, return point in the reference domain
Point<2> RinglebTop::pull_back(const Point<2>& x) const
{
double Vmin=0.2, Vmax=1.8; //these are just used in initial guess for iteration process to numerically solve V
unsigned int imax;
double b, rho, L, V, psi, Vp, tol;
imax =500;
tol = 1.0e-15;
//Initial guess, Vp=1
Vp = 0.5*( Vmin + Vmax );
//Start iterations. Eq.(7.13.76). given point x, return V
for(unsigned int i=1; i<imax; ++i){
b = std::sqrt(1-0.2*Vp*Vp);
rho = std::pow(b,5);
L = 1./b+1./(3*b*b*b)+1./(5*rho)-0.5*std::log((1+b)/(1-b));
V = std::sqrt(1/2./rho/std::sqrt((x[0]-L/2)*(x[0]-L/2)+x[1]*x[1])); // Eq.(7.11.76)
if(std::fabs(V-Vp)<tol)
break;
if(i == imax){
std::cout<<"did not converge...sorry."<<std::endl;
exit(0);
}
Vp = V;
}
psi = std::sqrt(0.5/V/V-rho*(x[0]-L/2));
/*std::cout<<"given real_space point, (x,y) = "<<x[0]<<","<<x[1]
<<std::endl<<"return reference_space point, (psi, V)= "<<psi<<","<<V
<<std::endl;*/ //used in test
return Point<2>(psi, V);
}
//given point in the reference domain, return point in the real space
Point<2> RinglebTop::push_forward(const Point<2>& psi_V) const
{
const double psi = psi_V[0];
const double V = psi_V[1];
const double b = std::sqrt(1-0.2*V*V);
const double rho = std::pow(b,5);
const double L = 1./b+1./(3*b*b*b)+1./(5*rho)-0.5*std::log((1+b)/(1-b));
if(std::fabs(psi*V-1) < 1.0e-15){
return Point<2> ((1./rho)*(0.5/V/V- psi*psi) + 0.5*L, //x, Eq.(7.11.66)
0); //y, Eq.(7.11.67)), to avoid sqrt(negative) below.
}
else{
/* std::cout<<"given reference_space point, (psi, v) = "<<psi_V[0]<<","<<psi_V[1]
<<std::endl<<"return real_space point, (x,y)= "<<(1./rho)*(0.5/V/V- psi*psi) + 0.5*L
<<","<<psi/(rho*V)*std::sqrt(1.-psi*psi*V*V)
<<std::endl;*/
return Point<2> ((1./rho)*(0.5/V/V- psi*psi) + 0.5*L, //x, Eq.(7.11.66)
psi/(rho*V)*std::sqrt(1.-psi*psi*V*V)); //y, Eq.(7.11.67))
}
}
//----------------------------------------------------------------------------------------------------
// NACA0012 airfoil Manifold. Just project the candidate point onto the true boundary in the normal
// direction to the boundary edge. see the related documentation
//----------------------------------------------------------------------------------------------------
//used in newton method for projection to manifold. m=dy/dx, P is the candidate point,
//x is P's x_coordinate_value. this function will modify f and d
void NACA0012::func(const double& m, const Point<2>& P, const double& x, double& f, double& d) const
{
//set and mark independent variable
Sacado::Fad::DFad<double> xd = x;
xd.diff(0,1); //set xd to be dof 0, in a 1-dof system.
Sacado::Fad::DFad<double> y = naca0012(xd);
if(P(1)<0) y=-y; //get the corresponding y_coordinate_value
Sacado::Fad::DFad<double> F = m * (y-P[1]) + (xd - P[0]);
f = F.val();
d = F.fastAccessDx(0); //dF/d(xd)
}
//note that candidate point is calculated as the convex combination of the given points
//by Manifold::get_intermediate_point()
Point<2> NACA0012::project_to_manifold (const std::vector<Point<2>>& surrounding_points,
const Point<2>& candidate) const
{
const double GEOMTOL = 1.0e-13;
Assert((surrounding_points[0][1] > -GEOMTOL && surrounding_points[1][1] >-GEOMTOL) ||
(surrounding_points[0][1] < GEOMTOL && surrounding_points[1][1] < GEOMTOL),
ExcMessage("End points not on same side of airfoil"));
const double dx = surrounding_points[1][0] - surrounding_points[0][0];
const double dy = surrounding_points[1][1] - surrounding_points[0][1];
Assert (std::fabs(dx) > GEOMTOL, ExcMessage("dx is too small"));
const double m = dy/dx;
//initial guess for newton. x_coordinate should be between 0 and 1.
double x = candidate[0];
x = std::min(x, 1.0);
x = std::max(x, 0.0);
double f, d;
func(m, candidate, x, f, d); //change f and d
unsigned int i = 0, ITMAX = 10;
while(std::fabs(f) > GEOMTOL)
{
double s = -f/d;
while(x+s < 0 || x+s>1 ) s *= 0.5;
x = x + s; //modify x until converge
func(m, candidate, x, f, d);
++i;
AssertThrow(i < ITMAX, ExcMessage("Newton did not converge"));
}
//now x is computed
double y = naca0012(x);
if(candidate[1] < 0) y = -y;
Point<2> p(x, y);
return p;
}