Go to the documentation of this file.00001 #include <OpenFOAM/scalar.H>
00002 #include <OpenFOAM/vector.H>
00003 #include "curveTools.H"
00004
00005
00006
00007 namespace Foam
00008 {
00009
00010
00011
00012 scalar distance(const vector& p1, const vector& p2)
00013 {
00014 return mag(p2 - p1);
00015 }
00016
00017
00018 bool stepForwardsToNextPoint
00019 (
00020 const vector& o,
00021 vector& n,
00022 label& i,
00023 label& ip1,
00024 scalar l,
00025 const curve& Curve
00026 )
00027 {
00028 label ip1n = ip1-1;
00029 while (++ip1n < Curve.size() && distance(o, Curve[ip1n]) < l);
00030 label in = ip1n - 1;
00031
00032 bool eoc = true;
00033
00034 if (ip1n < Curve.size() && in >= 0)
00035 {
00036 eoc = interpolate(Curve[in], Curve[ip1n], o, n, l);
00037
00038 i = in;
00039 ip1 = ip1n;
00040 }
00041
00042 return eoc;
00043 }
00044
00045
00046 bool stepBackwardsToNextPoint
00047 (
00048 const vector& o,
00049 vector& n,
00050 label& i,
00051 label& ip1,
00052 scalar l,
00053 const curve& Curve
00054 )
00055 {
00056 label ip1n = ip1+1;
00057 while (--ip1n >= 0 && distance(o, Curve[ip1n]) < l);
00058 label in = ip1n + 1;
00059
00060 bool eoc = true;
00061
00062 if (ip1n >= 0 && in < Curve.size())
00063 {
00064 eoc = interpolate(Curve[in], Curve[ip1n], o, n, l);
00065
00066 i = in;
00067 ip1 = ip1n;
00068 }
00069
00070 return eoc;
00071 }
00072
00073
00074 bool interpolate
00075 (
00076 const vector& p1,
00077 const vector& p2,
00078 const vector& o,
00079 vector& n,
00080 scalar l
00081 )
00082 {
00083 vector D = p1 - p2;
00084 scalar a = magSqr(D);
00085 scalar b = 2.0*(D&(p2 - o));
00086 scalar c = magSqr(p2) + (o&(o - 2.0*p2)) - l*l;
00087
00088 scalar b2m4ac = b*b - 4.0*a*c;
00089
00090 if (b2m4ac >= 0.0)
00091 {
00092 scalar srb2m4ac = sqrt(b2m4ac);
00093
00094 scalar lamda = (-b - srb2m4ac)/(2.0*a);
00095
00096 if (lamda > 1.0+curveSmall || lamda < -curveSmall)
00097 {
00098 lamda = (-b + srb2m4ac)/(2.0*a);
00099 }
00100
00101 if (lamda < 1.0+curveSmall && lamda > -curveSmall)
00102 {
00103 n = p2 + lamda*(p1 - p2);
00104
00105 return false;
00106 }
00107 else
00108 {
00109 return true;
00110 }
00111 }
00112 else
00113 {
00114 return true;
00115 }
00116 }
00117
00118
00119
00120 bool XstepForwardsToNextPoint
00121 (
00122 const vector& o,
00123 vector& n,
00124 label& i,
00125 label& ip1,
00126 scalar l,
00127 const curve& Curve
00128 )
00129 {
00130 label ip1n = ip1-1;
00131 while (++ip1n < Curve.size() && mag(o.x() - Curve[ip1n].x()) < l);
00132 label in = ip1n - 1;
00133
00134 bool eoc = true;
00135
00136 if (ip1n < Curve.size() && in >= 0)
00137 {
00138 eoc = Xinterpolate(Curve[in], Curve[ip1n], o, n, l);
00139
00140 i = in;
00141 ip1 = ip1n;
00142 }
00143
00144 return eoc;
00145 }
00146
00147
00148
00149 bool Xinterpolate
00150 (
00151 const vector& p1,
00152 const vector& p2,
00153 const vector& o,
00154 vector& n,
00155 scalar l
00156 )
00157 {
00158 n.x() = o.x() + l;
00159
00160 if (p2.x() < o.x() + l - curveSmall && p2.x() > o.x() - l + curveSmall)
00161 {
00162 return true;
00163 }
00164
00165 if (p2.x() < o.x() + l)
00166 {
00167 n.x() = o.x() - l;
00168 }
00169
00170 vector D = p2 - p1;
00171 scalar lamda = (n.x() - p1.x())/D.x();
00172 n.y() = p1.y() + lamda*D.y();
00173 return false;
00174
00175 }
00176
00177
00178
00179
00180 }
00181
00182