Minimum changes to compile everything with gcc-4.3.0

This commit is contained in:
henry 2008-05-26 19:27:23 +01:00
parent f6107f4033
commit 0f687ccd76
22 changed files with 158 additions and 135 deletions

View File

@ -58,7 +58,7 @@ int main(int argc, char *argv[])
<< "reactions" << cr.reactions() << ';' << endl;
OFstream thermoFile(FOAMThermodynamicsFileName);
thermoFile<< cr.specieThermo() << endl;
thermoFile<< cr.speciesThermo() << endl;
Info << "End\n" << endl;

View File

@ -86,7 +86,14 @@ set WM_COMPILER_INST=OpenFOAM
switch ("$WM_COMPILER_INST")
case OpenFOAM:
setenv WM_COMPILER_DIR $FOAM_INST_DIR/$WM_ARCH/gcc-4.2.2$WM_COMPILER_ARCH
switch ("$WM_COMPILER")
case Gcc43:
setenv WM_COMPILER_DIR $FOAM_INST_DIR/$WM_ARCH/gcc-4.3.0$WM_COMPILER_ARCH
breaksw
case Gcc:
setenv WM_COMPILER_DIR $FOAM_INST_DIR/$WM_ARCH/gcc-4.2.2$WM_COMPILER_ARCH
breaksw
endsw
# Check that the compiler directory can be found
if ( ! -d "$WM_COMPILER_DIR" ) then

View File

@ -318,7 +318,8 @@ bool HashTable<T, Key, Hash>::erase(const iterator& cit)
delete it.elmtPtr_;
// Search back for previous non-zero table entry
while (--it.hashIndex_ >= 0 && !table_[it.hashIndex_]);
while (--it.hashIndex_ >= 0 && !table_[it.hashIndex_])
{}
if (it.hashIndex_ >= 0)
{

View File

@ -220,7 +220,8 @@ HashTable<T, Key, Hash>::iterator::operator++()
(
++hashIndex_ < curHashTable_.tableSize_
&& !(elmtPtr_ = curHashTable_.table_[hashIndex_])
);
)
{}
if (hashIndex_ == curHashTable_.tableSize_)
{
@ -259,7 +260,8 @@ HashTable<T, Key, Hash>::begin()
{
label i = 0;
while (table_ && !table_[i] && ++i < tableSize_);
while (table_ && !table_[i] && ++i < tableSize_)
{}
if (i == tableSize_)
{
@ -396,7 +398,8 @@ HashTable<T, Key, Hash>::const_iterator::operator++()
(
++hashIndex_ < curHashTable_.tableSize_
&& !(elmtPtr_ = curHashTable_.table_[hashIndex_])
);
)
{}
}
return *this;
@ -430,7 +433,8 @@ HashTable<T, Key, Hash>::begin() const
{
label i = 0;
while (table_ && !table_[i] && ++i < tableSize_);
while (table_ && !table_[i] && ++i < tableSize_)
{}
if (i == tableSize_)
{

View File

@ -34,6 +34,7 @@ License
#include "JobInfo.H"
#include "labelList.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
@ -389,7 +390,8 @@ Foam::argList::argList
rootPath_/globalCase_/"processor"
+ name(++nProcDirs)
)
);
)
{}
if (nProcDirs != Pstream::nProcs())
{

View File

@ -29,6 +29,7 @@ Description
#include <new>
#include <iostream>
#include <cstdlib>
namespace Foam
{

View File

@ -64,7 +64,8 @@ Type interpolateXY
label n = xOld.size();
label lo = 0;
for (lo=0; lo<n && xOld[lo]>x; ++lo);
for (lo=0; lo<n && xOld[lo]>x; ++lo)
{}
label low = lo;
if (low < n)
@ -79,7 +80,8 @@ Type interpolateXY
}
label hi = 0;
for (hi=0; hi<n && xOld[hi]<x; ++hi);
for (hi=0; hi<n && xOld[hi]<x; ++hi)
{}
label high = hi;
if (high < n)

View File

@ -272,7 +272,7 @@ inline vector triangle<Point, PointRef>::circumCentre() const
scalar c = c1 + c2 + c3;
return
return
(
((c2 + c3)*a_ + (c3 + c1)*b_ + (c1 + c2)*c_)/(2*c)
);
@ -397,7 +397,7 @@ inline pointHit triangle<Point, PointRef>::ray
bool eligible =
alg == intersection::FULL_RAY
|| alg == intersection::HALF_RAY && dist > -planarPointTol
|| (alg == intersection::HALF_RAY && dist > -planarPointTol)
|| (
alg == intersection::VISIBLE
&& ((q1 & normal()) < -VSMALL)
@ -434,7 +434,7 @@ template<class Point, class PointRef>
inline pointHit triangle<Point, PointRef>::intersection
(
const point& orig,
const vector& dir,
const vector& dir,
const intersection::algorithm alg
) const
{
@ -617,7 +617,7 @@ inline bool triangle<Point, PointRef>::classify
scalar beta = 0;
bool hit = false;
if (Foam::mag(u1) < SMALL)
{
beta = u0/u2;

View File

@ -73,7 +73,8 @@ long long readLongLong(Istream& is)
static const label zeroOffset = int('0');
// Get next non-whitespace character
while (is.read(c) && isspace(c));
while (is.read(c) && isspace(c))
{}
do
{

View File

@ -49,6 +49,8 @@ SourceFiles
#include "char.H"
#include <string>
#include <cstring>
#include <cstdlib>
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

View File

@ -103,7 +103,7 @@ inline bool Foam::string::stripInvalid(string& s)
++nValid;
}
}
s.resize(nValid);
return true;

View File

@ -180,7 +180,7 @@ public:
return meshInfo_;
}
const label size() const
label size() const
{
return IDLList<ParticleType>::size();
};

View File

@ -1055,7 +1055,7 @@ pointIndexHit indexedOctree<Type>::findLine
direction startBit = treeBb.posBits(start);
direction endBit = treeBb.posBits(end);
if (startBit&endBit != 0)
if ((startBit & endBit) != 0)
{
// Both start and end outside domain and in same block.
return pointIndexHit(false, vector::zero, -1);

View File

@ -4,6 +4,7 @@
# include <fstream>
# include <cmath>
# include <ctime>
# include <cstring>
using namespace std;
@ -24,7 +25,7 @@ double d_epsilon ( void )
// D_EPSILON is a number R which is a power of 2 with the property that,
// to the precision of the computer's arithmetic,
// 1 < 1 + R
// but
// but
// 1 = ( 1 + R / 2 )
//
// Modified:
@ -79,7 +80,7 @@ double d_max ( double x, double y )
if ( y < x )
{
return x;
}
}
else
{
return y;
@ -113,7 +114,7 @@ double d_min ( double x, double y )
if ( y < x )
{
return y;
}
}
else
{
return x;
@ -623,7 +624,7 @@ void d2vec_sort_quick_a ( int n, double a[] )
}
//******************************************************************************
int diaedg ( double x0, double y0, double x1, double y1, double x2, double y2,
int diaedg ( double x0, double y0, double x1, double y1, double x2, double y2,
double x3, double y3 )
//******************************************************************************
@ -660,7 +661,7 @@ int diaedg ( double x0, double y0, double x1, double y1, double x2, double y2,
//
// Parameters:
//
// Input, double X0, Y0, X1, Y1, X2, Y2, X3, Y3, the coordinates of the
// Input, double X0, Y0, X1, Y1, X2, Y2, X3, Y3, the coordinates of the
// vertices of a quadrilateral, given in counter clockwise order.
//
// Output, int DIAEDG, chooses a diagonal:
@ -696,12 +697,12 @@ int diaedg ( double x0, double y0, double x1, double y1, double x2, double y2,
dx32 = x3 - x2;
dy32 = y3 - y2;
tola = tol * d_max ( fabs ( dx10 ),
d_max ( fabs ( dy10 ),
tola = tol * d_max ( fabs ( dx10 ),
d_max ( fabs ( dy10 ),
d_max ( fabs ( dx30 ), fabs ( dy30 ) ) ) );
tolb = tol * d_max ( fabs ( dx12 ),
d_max ( fabs ( dy12 ),
tolb = tol * d_max ( fabs ( dx12 ),
d_max ( fabs ( dy12 ),
d_max ( fabs ( dx32 ), fabs ( dy32 ) ) ) );
ca = dx10 * dx30 + dy10 * dy30;
@ -718,7 +719,7 @@ int diaedg ( double x0, double y0, double x1, double y1, double x2, double y2,
else
{
tola = d_max ( tola, tolb );
s = ( dx10 * dy30 - dx30 * dy10 ) * cb
s = ( dx10 * dy30 - dx30 * dy10 ) * cb
+ ( dx32 * dy12 - dx12 * dy32 ) * ca;
if ( tola < s )
@ -771,7 +772,7 @@ void dmat_transpose_print ( int m, int n, double a[], const char *title )
}
//******************************************************************************
void dmat_transpose_print_some ( int m, int n, double a[], int ilo, int jlo,
void dmat_transpose_print_some ( int m, int n, double a[], int ilo, int jlo,
int ihi, int jhi, const char *title )
//******************************************************************************
@ -909,7 +910,7 @@ void dmat_uniform ( int m, int n, double b, double c, int *seed, double r[] )
//
// Input/output, int *SEED, the "seed" value. Normally, this
// value should not be 0, otherwise the output value of SEED
// will still be 0, and D_UNIFORM will be 0. On output, SEED has
// will still be 0, and D_UNIFORM will be 0. On output, SEED has
// been updated.
//
// Output, double DMAT_UNIFORM[M*N], a matrix of pseudorandom values.
@ -943,7 +944,7 @@ void dmat_uniform ( int m, int n, double b, double c, int *seed, double r[] )
}
//******************************************************************************
int dtris2 ( int point_num, double point_xy[], int *tri_num,
int dtris2 ( int point_num, double point_xy[], int *tri_num,
int tri_vert[], int tri_nabe[] )
//******************************************************************************
@ -1046,11 +1047,11 @@ int dtris2 ( int point_num, double point_xy[], int *tri_num,
for ( j = 0; j <= 1; j++ )
{
cmax = d_max ( fabs ( point_xy[2*(m-1)+j] ),
cmax = d_max ( fabs ( point_xy[2*(m-1)+j] ),
fabs ( point_xy[2*(m1-1)+j] ) );
if ( tol * ( cmax + 1.0 )
< fabs ( point_xy[2*(m-1)+j] - point_xy[2*(m1-1)+j] ) )
if ( tol * ( cmax + 1.0 )
< fabs ( point_xy[2*(m-1)+j] - point_xy[2*(m1-1)+j] ) )
{
k = j;
break;
@ -1095,7 +1096,7 @@ int dtris2 ( int point_num, double point_xy[], int *tri_num,
m = j;
lr = lrline ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1],
point_xy[2*(m1-1)+0], point_xy[2*(m1-1)+1],
point_xy[2*(m1-1)+0], point_xy[2*(m1-1)+1],
point_xy[2*(m2-1)+0], point_xy[2*(m2-1)+1], 0.0 );
if ( lr != 0 )
@ -1182,8 +1183,8 @@ int dtris2 ( int point_num, double point_xy[], int *tri_num,
m2 = tri_vert[3*(ltri-1)+0];
}
lr = lrline ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1],
point_xy[2*(m1-1)+0], point_xy[2*(m1-1)+1],
lr = lrline ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1],
point_xy[2*(m1-1)+0], point_xy[2*(m1-1)+1],
point_xy[2*(m2-1)+0], point_xy[2*(m2-1)+1], 0.0 );
if ( 0 < lr )
@ -1199,7 +1200,7 @@ int dtris2 ( int point_num, double point_xy[], int *tri_num,
redg = (l % 3) + 1;
}
vbedg ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1], point_num,
vbedg ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1], point_num,
point_xy, *tri_num, tri_vert, tri_nabe, &ltri, &ledg, &rtri, &redg );
n = *tri_num + 1;
@ -1255,7 +1256,7 @@ int dtris2 ( int point_num, double point_xy[], int *tri_num,
ltri = n;
ledg = 2;
error = swapec ( m, &top, &ltri, &ledg, point_num, point_xy, *tri_num,
error = swapec ( m, &top, &ltri, &ledg, point_num, point_xy, *tri_num,
tri_vert, tri_nabe, stack );
if ( error != 0 )
@ -1473,9 +1474,9 @@ void dvec_print ( int n, double a[], const char *title )
}
cout << "\n";
for ( i = 0; i <= n-1; i++ )
for ( i = 0; i <= n-1; i++ )
{
cout << setw(6) << i + 1 << " "
cout << setw(6) << i + 1 << " "
<< setw(14) << a[i] << "\n";
}
@ -1600,8 +1601,8 @@ int i_modp ( int i, int j )
//
// Formula:
//
// If
// NREM = I_MODP ( I, J )
// If
// NREM = I_MODP ( I, J )
// NMULT = ( I - NREM ) / J
// then
// I = J * NMULT + NREM
@ -1620,7 +1621,7 @@ int i_modp ( int i, int j )
// Examples:
//
// I J MOD I_MODP I_MODP Factorization
//
//
// 107 50 7 7 107 = 2 * 50 + 7
// 107 -50 7 7 107 = -2 * -50 + 7
// -107 50 -7 43 -107 = -3 * 50 + 43
@ -1640,7 +1641,7 @@ int i_modp ( int i, int j )
//
// Input, int J, the number that divides I.
//
// Output, int I_MODP, the nonnegative remainder when I is
// Output, int I_MODP, the nonnegative remainder when I is
// divided by J.
//
{
@ -1692,7 +1693,7 @@ int i_sign ( int i )
//
// Output, int I_SIGN, the sign of I.
{
if ( i < 0 )
if ( i < 0 )
{
return (-1);
}
@ -1809,7 +1810,7 @@ void imat_transpose_print ( int m, int n, int a[], const char *title )
}
//******************************************************************************
void imat_transpose_print_some ( int m, int n, int a[], int ilo, int jlo,
void imat_transpose_print_some ( int m, int n, int a[], int ilo, int jlo,
int ihi, int jhi, const char *title )
//******************************************************************************
@ -1962,7 +1963,7 @@ void ivec_heap_d ( int n, int a[] )
// Only nodes (N/2)-1 down to 0 can be "parent" nodes.
//
for ( i = (n/2)-1; 0 <= i; i-- )
{
{
//
// Copy the value out of the parent node.
// Position IFREE is now "open".
@ -2180,7 +2181,7 @@ void ivec_sorted_unique ( int n, int a[], int *nuniq )
for ( i = 1; i < n; i++ )
{
if ( a[i] != a[*nuniq] )
if ( a[i] != a[*nuniq] )
{
*nuniq = *nuniq + 1;
a[*nuniq] = a[i];
@ -2192,7 +2193,7 @@ void ivec_sorted_unique ( int n, int a[], int *nuniq )
}
//******************************************************************************
int lrline ( double xu, double yu, double xv1, double yv1, double xv2,
int lrline ( double xu, double yu, double xv1, double yv1, double xv2,
double yv2, double dv )
//******************************************************************************
@ -2253,9 +2254,9 @@ int lrline ( double xu, double yu, double xv1, double yv1, double xv2,
dxu = xu - xv1;
dyu = yu - yv1;
tolabs = tol * d_max ( fabs ( dx ),
d_max ( fabs ( dy ),
d_max ( fabs ( dxu ),
tolabs = tol * d_max ( fabs ( dx ),
d_max ( fabs ( dy ),
d_max ( fabs ( dxu ),
d_max ( fabs ( dyu ), fabs ( dv ) ) ) ) );
t = dy * dxu - dx * dyu + dv * sqrt ( dx * dx + dy * dy );
@ -2426,7 +2427,7 @@ int *points_delaunay_naive_2d ( int n, double p[], int *ntri )
//********************************************************************
//
// Purpose:
// Purpose:
//
// POINTS_DELAUNAY_NAIVE_2D computes the Delaunay triangulation in 2D.
//
@ -2456,7 +2457,7 @@ int *points_delaunay_naive_2d ( int n, double p[], int *ntri )
//
// Joseph O'Rourke,
// Computational Geometry,
// Cambridge University Press,
// Cambridge University Press,
// Second Edition, 1998, page 187.
//
// Parameters:
@ -2467,7 +2468,7 @@ int *points_delaunay_naive_2d ( int n, double p[], int *ntri )
//
// Output, int *NTRI, the number of triangles.
//
// Output, int POINTS_DELAUNAY_NAIVE_2D[3*NTRI], the indices of the
// Output, int POINTS_DELAUNAY_NAIVE_2D[3*NTRI], the indices of the
// nodes making each triangle.
//
{
@ -2506,34 +2507,34 @@ int *points_delaunay_naive_2d ( int n, double p[], int *ntri )
//
// For each triple (I,J,K):
//
for ( i = 0; i < n - 2; i++ )
for ( i = 0; i < n - 2; i++ )
{
for ( j = i+1; j < n; j++ )
for ( j = i+1; j < n; j++ )
{
for ( k = i+1; k < n; k++ )
for ( k = i+1; k < n; k++ )
{
if ( j != k )
if ( j != k )
{
xn = ( p[1+j*2] - p[1+i*2] ) * ( z[k] - z[i] )
xn = ( p[1+j*2] - p[1+i*2] ) * ( z[k] - z[i] )
- ( p[1+k*2] - p[1+i*2] ) * ( z[j] - z[i] );
yn = ( p[0+k*2] - p[0+i*2] ) * ( z[j] - z[i] )
yn = ( p[0+k*2] - p[0+i*2] ) * ( z[j] - z[i] )
- ( p[0+j*2] - p[0+i*2] ) * ( z[k] - z[i] );
zn = ( p[0+j*2] - p[0+i*2] ) * ( p[1+k*2] - p[1+i*2] )
zn = ( p[0+j*2] - p[0+i*2] ) * ( p[1+k*2] - p[1+i*2] )
- ( p[0+k*2] - p[0+i*2] ) * ( p[1+j*2] - p[1+i*2] );
flag = ( zn < 0 );
if ( flag )
if ( flag )
{
for ( m = 0; m < n; m++ )
for ( m = 0; m < n; m++ )
{
flag = flag && ( ( p[0+m*2] - p[0+i*2] ) * xn
+ ( p[1+m*2] - p[1+i*2] ) * yn
flag = flag && ( ( p[0+m*2] - p[0+i*2] ) * xn
+ ( p[1+m*2] - p[1+i*2] ) * yn
+ ( z[m] - z[i] ) * zn <= 0 );
}
}
if ( flag )
if ( flag )
{
if ( pass == 2 )
{
@ -2545,8 +2546,8 @@ int *points_delaunay_naive_2d ( int n, double p[], int *ntri )
}
}
}
}
}
}
}
}
@ -2587,7 +2588,7 @@ int s_len_trim ( const char *s )
n = strlen ( s );
t = const_cast<char*>(s) + n - 1;
while ( 0 < n )
while ( 0 < n )
{
if ( *t != ' ' )
{
@ -2601,8 +2602,8 @@ int s_len_trim ( const char *s )
}
//******************************************************************************
int swapec ( int i, int *top, int *btri, int *bedg, int point_num,
double point_xy[], int tri_num, int tri_vert[], int tri_nabe[],
int swapec ( int i, int *top, int *btri, int *bedg, int point_num,
double point_xy[], int tri_num, int tri_vert[], int tri_nabe[],
int stack[] )
//******************************************************************************
@ -2698,7 +2699,7 @@ int swapec ( int i, int *top, int *btri, int *bedg, int point_num,
for ( ; ; )
{
if ( *top <= 0 )
if ( *top <= 0 )
{
break;
}
@ -2741,7 +2742,7 @@ int swapec ( int i, int *top, int *btri, int *bedg, int point_num,
c = tri_vert[3*(u-1)+1];
}
swap = diaedg ( x, y,
swap = diaedg ( x, y,
point_xy[2*(a-1)+0], point_xy[2*(a-1)+1],
point_xy[2*(c-1)+0], point_xy[2*(c-1)+1],
point_xy[2*(b-1)+0], point_xy[2*(b-1)+1] );
@ -3021,16 +3022,16 @@ double *triangle_circumcenter_2d ( double t[] )
center = new double[DIM_NUM];
asq = ( t[0+1*2] - t[0+0*2] ) * ( t[0+1*2] - t[0+0*2] )
asq = ( t[0+1*2] - t[0+0*2] ) * ( t[0+1*2] - t[0+0*2] )
+ ( t[1+1*2] - t[1+0*2] ) * ( t[1+1*2] - t[1+0*2] );
csq = ( t[0+2*2] - t[0+0*2] ) * ( t[0+2*2] - t[0+0*2] )
csq = ( t[0+2*2] - t[0+0*2] ) * ( t[0+2*2] - t[0+0*2] )
+ ( t[1+2*2] - t[1+0*2] ) * ( t[1+2*2] - t[1+0*2] );
top1 = ( t[1+1*2] - t[1+0*2] ) * csq - ( t[1+2*2] - t[1+0*2] ) * asq;
top2 = ( t[0+1*2] - t[0+0*2] ) * csq - ( t[0+2*2] - t[0+0*2] ) * asq;
bot = ( t[1+1*2] - t[1+0*2] ) * ( t[0+2*2] - t[0+0*2] )
bot = ( t[1+1*2] - t[1+0*2] ) * ( t[0+2*2] - t[0+0*2] )
- ( t[1+2*2] - t[1+0*2] ) * ( t[0+1*2] - t[0+0*2] );
center[0] = t[0+0*2] + 0.5 * top1 / bot;
@ -3042,7 +3043,7 @@ double *triangle_circumcenter_2d ( double t[] )
}
//******************************************************************************
bool triangulation_plot_eps ( const char *file_out_name, int g_num, double g_xy[],
bool triangulation_plot_eps ( const char *file_out_name, int g_num, double g_xy[],
int tri_num, int nod_tri[] )
//******************************************************************************
@ -3140,7 +3141,7 @@ bool triangulation_plot_eps ( const char *file_out_name, int g_num, double g_xy[
file_out << "%%Title: " << file_out_name << "\n";
file_out << "%%CreationDate: " << date_time << "\n";
file_out << "%%Pages: 1\n";
file_out << "%%Bounding Box: " << x_ps_min << " " << y_ps_min << " "
file_out << "%%Bounding Box: " << x_ps_min << " " << y_ps_min << " "
<< x_ps_max << " " << y_ps_max << "\n";
file_out << "%%Document-Fonts: Times-Roman\n";
file_out << "%%LanguageLevel: 1\n";
@ -3199,17 +3200,17 @@ bool triangulation_plot_eps ( const char *file_out_name, int g_num, double g_xy[
for ( g = 0; g < g_num; g++ )
{
x_ps = int(
( ( x_max - g_xy[0+g*2] ) * double( x_ps_min )
+ ( g_xy[0+g*2] - x_min ) * double( x_ps_max ) )
x_ps = int(
( ( x_max - g_xy[0+g*2] ) * double( x_ps_min )
+ ( g_xy[0+g*2] - x_min ) * double( x_ps_max ) )
/ ( x_max - x_min ) );
y_ps = int(
( ( y_max - g_xy[1+g*2] ) * double( y_ps_min )
+ ( g_xy[1+g*2] - y_min ) * double( y_ps_max ) )
y_ps = int(
( ( y_max - g_xy[1+g*2] ) * double( y_ps_min )
+ ( g_xy[1+g*2] - y_min ) * double( y_ps_max ) )
/ ( y_max - y_min ) );
file_out << "newpath " << x_ps << " "
file_out << "newpath " << x_ps << " "
<< y_ps << " 5 0 360 arc closepath fill\n";
}
@ -3231,14 +3232,14 @@ bool triangulation_plot_eps ( const char *file_out_name, int g_num, double g_xy[
k = nod_tri[3*(t-1)+e-1];
x_ps = int(
( ( x_max - g_xy[0+(k-1)*2] ) * double( x_ps_min )
+ ( g_xy[0+(k-1)*2] - x_min ) * double( x_ps_max ) )
x_ps = int(
( ( x_max - g_xy[0+(k-1)*2] ) * double( x_ps_min )
+ ( g_xy[0+(k-1)*2] - x_min ) * double( x_ps_max ) )
/ ( x_max - x_min ) );
y_ps = int(
( ( y_max - g_xy[1+(k-1)*2] ) * double( y_ps_min )
+ ( g_xy[1+(k-1)*2] - y_min ) * double( y_ps_max ) )
y_ps = int(
( ( y_max - g_xy[1+(k-1)*2] ) * double( y_ps_min )
+ ( g_xy[1+(k-1)*2] - y_min ) * double( y_ps_max ) )
/ ( y_max - y_min ) );
if ( j == 1 )
@ -3269,7 +3270,7 @@ bool triangulation_plot_eps ( const char *file_out_name, int g_num, double g_xy[
}
//******************************************************************************
void triangulation_print ( int point_num, double xc[], int tri_num,
void triangulation_print ( int point_num, double xc[], int tri_num,
int tri_vert[], int tri_nabe[] )
//******************************************************************************
@ -3356,7 +3357,7 @@ void triangulation_print ( int point_num, double xc[], int tri_num,
imat_transpose_print ( 3, tri_num, tri_nabe, " Triangle neighbors" );
//
// Determine VERTEX_NUM, the number of vertices. This is not
// Determine VERTEX_NUM, the number of vertices. This is not
// the same as the number of points!
//
vertex_list = new int[3*tri_num];
@ -3439,7 +3440,7 @@ void triangulation_print ( int point_num, double xc[], int tri_num,
}
//******************************************************************************
void vbedg ( double x, double y, int point_num, double point_xy[], int tri_num,
void vbedg ( double x, double y, int point_num, double point_xy[], int tri_num,
int tri_vert[], int tri_nabe[], int *ltri, int *ledg, int *rtri, int *redg )
//******************************************************************************

View File

@ -118,7 +118,7 @@ public:
// Member Functions
virtual const speciesTable& species() const = 0;
virtual const HashPtrTable<reactionThermo>& specieThermo() const = 0;
virtual const HashPtrTable<reactionThermo>& speciesThermo() const = 0;
virtual const SLPtrList<reaction>& reactions() const = 0;
};

View File

@ -607,15 +607,15 @@ bool finishReaction = false;
HashPtrTable<reactionThermo>::iterator specieThermoIter
(
specieThermo_.find(currentSpecieName)
speciesThermo_.find(currentSpecieName)
);
if (specieThermoIter != specieThermo_.end())
if (specieThermoIter != speciesThermo_.end())
{
specieThermo_.erase(specieThermoIter);
speciesThermo_.erase(specieThermoIter);
}
specieThermo_.insert
speciesThermo_.insert
(
currentSpecieName,
new reactionThermo
@ -742,7 +742,7 @@ bool finishReaction = false;
<< "Plasma momentum-transfer in reaction on line "
<< lineNo_ << "not yet supported"
<< exit(FatalError);
BEGIN(readReactionKeyword);
break;
}
@ -894,7 +894,7 @@ bool finishReaction = false;
FatalErrorIn("chemkinReader::lex()")
<< "HIGH keyword given for a chemically"
" activated bimolecular reaction which does not"
" contain a pressure dependent specie"
" contain a pressure dependent specie"
<< " on line " << lineNo_
<< exit(FatalError);
}
@ -935,7 +935,7 @@ bool finishReaction = false;
<< exit(FatalError);
}
if
if
(
fofType == unknownFallOffFunctionType
|| fofType == Lindemann
@ -969,7 +969,7 @@ bool finishReaction = false;
<< exit(FatalError);
}
if
if
(
fofType == unknownFallOffFunctionType
|| fofType == Lindemann
@ -1046,7 +1046,7 @@ bool finishReaction = false;
}
rrType = LandauTeller;
reactionCoeffsName =
reactionCoeffsName =
word(reactionTypeNames[rType])
+ reactionRateTypeNames[rrType];
BEGIN(readReactionCoeffs);
@ -1316,7 +1316,7 @@ bool finishReaction = false;
<readPDependentSpecie>{pDependentSpecie} {
word rhsPDependentSpecieName = pDependentSpecieName;
pDependentSpecieName =
pDependentSpecieName =
foamName(foamSpecieString(YYText()));
pDependentSpecieName =
pDependentSpecieName(0, pDependentSpecieName.size() - 1);
@ -1412,7 +1412,7 @@ bool finishReaction = false;
else
{
FatalErrorIn("chemkinReader::lex()")
<< "unknown specie " << currentSpecieName
<< "unknown specie " << currentSpecieName
<< " given in reaction-order specification"
<< " on line " << lineNo_ << nl
<< "Valid species are : " << nl
@ -1449,7 +1449,7 @@ bool finishReaction = false;
}
FatalErrorIn("chemkinReader::lex()")
<< "Specie " << currentSpecieName
<< "Specie " << currentSpecieName
<< " on line " << lineNo_
<< " not present in " << side << " of reaction " << nl << lrhs
<< exit(FatalError);

View File

@ -52,7 +52,7 @@ namespace Foam
/* * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * */
const char* Foam::chemkinReader::reactionTypeNames[4] =
const char* Foam::chemkinReader::reactionTypeNames[4] =
{
"irreversible",
"reversible",
@ -154,7 +154,7 @@ void Foam::chemkinReader::checkCoeffs
<< " rate expression on line "
<< lineNo_-1 << ", should be "
<< nCoeffs << " but " << reactionCoeffs.size() << " supplied." << nl
<< "Coefficients are "
<< "Coefficients are "
<< reactionCoeffs << nl
<< exit(FatalError);
}
@ -182,7 +182,7 @@ void Foam::chemkinReader::addReactionType
speciesTable_,
lhs.shrink(),
rhs.shrink(),
specieThermo_
speciesThermo_
),
rr
)
@ -201,7 +201,7 @@ void Foam::chemkinReader::addReactionType
speciesTable_,
lhs.shrink(),
rhs.shrink(),
specieThermo_
speciesThermo_
),
rr
)
@ -291,7 +291,7 @@ void Foam::chemkinReader::addPressureDependentReaction
<< "Wrong number of coefficients for Troe rate expression"
" on line " << lineNo_-1 << ", should be 3 or 4 but "
<< TroeCoeffs.size() << " supplied." << nl
<< "Coefficients are "
<< "Coefficients are "
<< TroeCoeffs << nl
<< exit(FatalError);
}
@ -347,7 +347,7 @@ void Foam::chemkinReader::addPressureDependentReaction
<< "Wrong number of coefficients for SRI rate expression"
" on line " << lineNo_-1 << ", should be 3 or 5 but "
<< SRICoeffs.size() << " supplied." << nl
<< "Coefficients are "
<< "Coefficients are "
<< SRICoeffs << nl
<< exit(FatalError);
}
@ -397,7 +397,7 @@ void Foam::chemkinReader::addPressureDependentReaction
if (fofType < 4)
{
FatalErrorIn("chemkinReader::addPressureDependentReaction")
<< "Fall-off function type "
<< "Fall-off function type "
<< fallOffFunctionNames[fofType]
<< " on line " << lineNo_-1
<< " not implemented"
@ -459,7 +459,7 @@ void Foam::chemkinReader::addReaction
// Calculate the unit conversion factor for the A coefficient
// for the change from mol/cm^3 to kmol/m^3 concentraction units
const scalar concFactor = 0.001;
const scalar concFactor = 0.001;
scalar sumExp = 0.0;
forAll (lhs, i)
{
@ -500,7 +500,7 @@ void Foam::chemkinReader::addReaction
speciesTable_,
lhs.shrink(),
rhs.shrink(),
specieThermo_
speciesThermo_
),
ArrheniusReactionRate
(
@ -553,7 +553,7 @@ void Foam::chemkinReader::addReaction
speciesTable_,
lhs.shrink(),
rhs.shrink(),
specieThermo_
speciesThermo_
),
thirdBodyArrheniusReactionRate
(
@ -658,7 +658,7 @@ void Foam::chemkinReader::addReaction
speciesTable_,
lhs.shrink(),
rhs.shrink(),
specieThermo_
speciesThermo_
),
LandauTellerReactionRate
(
@ -814,7 +814,8 @@ void Foam::chemkinReader::read
yy_buffer_state* bufferPtr(yy_create_buffer(&thermoStream, yyBufSize));
yy_switch_to_buffer(bufferPtr);
while(lex() != 0);
while(lex() != 0)
{}
yy_delete_buffer(bufferPtr);
@ -838,7 +839,8 @@ void Foam::chemkinReader::read
initReactionKeywordTable();
while(lex() != 0);
while(lex() != 0)
{}
yy_delete_buffer(bufferPtr);
}
@ -874,7 +876,7 @@ Foam::chemkinReader::chemkinReader(const dictionary& thermoDict)
);
fileName thermoFileName = fileName::null;
if (thermoDict.found("CHEMKINThermoFile"))
{
thermoFileName =

View File

@ -198,7 +198,7 @@ private:
HashTable<phase> speciePhase_;
//- Table of the thermodynamic data given in the CHEMKIN file
HashPtrTable<reactionThermo> specieThermo_;
HashPtrTable<reactionThermo> speciesThermo_;
//- Table of species composition
HashTable<List<specieElement> > specieComposition_;
@ -363,9 +363,9 @@ public:
}
//- Table of the thermodynamic data given in the CHEMKIN file
const HashPtrTable<reactionThermo>& specieThermo() const
const HashPtrTable<reactionThermo>& speciesThermo() const
{
return specieThermo_;
return speciesThermo_;
}
//- Table of species composition

View File

@ -45,12 +45,12 @@ Foam::foamChemistryReader::foamChemistryReader
const fileName& thermoFileName
)
:
specieThermo_(IFstream(thermoFileName)()),
speciesThermo_(IFstream(thermoFileName)()),
speciesTable_(dictionary(IFstream(reactionsFileName)()).lookup("species")),
reactions_
(
dictionary(IFstream(reactionsFileName)()).lookup("reactions"),
reaction::iNew(speciesTable_, specieThermo_)
reaction::iNew(speciesTable_, speciesThermo_)
)
{}
@ -58,7 +58,7 @@ Foam::foamChemistryReader::foamChemistryReader
// Construct from components
Foam::foamChemistryReader::foamChemistryReader(const dictionary& thermoDict)
:
specieThermo_
speciesThermo_
(
IFstream
(
@ -84,7 +84,7 @@ Foam::foamChemistryReader::foamChemistryReader(const dictionary& thermoDict)
fileName(thermoDict.lookup("foamChemistryFile")).expand()
)()
).lookup("reactions"),
reaction::iNew(speciesTable_, specieThermo_)
reaction::iNew(speciesTable_, speciesThermo_)
)
{}

View File

@ -60,7 +60,7 @@ class foamChemistryReader
public chemistryReader
{
//- Table of the thermodynamic data given in the foamChemistry file
HashPtrTable<reactionThermo> specieThermo_;
HashPtrTable<reactionThermo> speciesThermo_;
//- Table of species
speciesTable speciesTable_;
@ -113,9 +113,9 @@ public:
}
//- Table of the thermodynamic data given in the foamChemistry file
const HashPtrTable<reactionThermo>& specieThermo() const
const HashPtrTable<reactionThermo>& speciesThermo() const
{
return specieThermo_;
return speciesThermo_;
}
//- List of the reactions

View File

@ -45,7 +45,7 @@ reactingMixture::reactingMixture
(
thermoDict,
autoPtr<chemistryReader>::operator()().species(),
autoPtr<chemistryReader>::operator()().specieThermo(),
autoPtr<chemistryReader>::operator()().speciesThermo(),
mesh
),
PtrList<chemistryReader::reaction>