From e42fb2445026a01f22c9919fe75f63f85ba78b00 Mon Sep 17 00:00:00 2001 From: Tomas Oberhuber <tomas.oberhuber@fjfi.cvut.cz> Date: Fri, 3 May 2013 14:28:06 +0200 Subject: [PATCH] Implementing The Lax Fridrichs scheme for deformed grids. Changing indexing in grids. --- examples/navier-stokes/Makefile | 11 +-- .../navier-stokes/navierStokesSolver_impl.h | 24 +++--- src/implementation/mesh/tnlGrid2D_impl.h | 46 ++++++++--- src/implementation/mesh/tnlGrid3D_impl.h | 16 +--- .../mesh/tnlIdenticalGridGeometry_impl.h | 46 +++++------ .../schemes/euler/fvm/tnlLaxFridrichs_impl.h | 81 ++++++++++++++++++- src/mesh/tnlGrid.h | 23 +++--- src/mesh/tnlIdenticalGridGeometry.h | 45 +++++------ tests/unit-tests/core/tnlArrayTest.cpp | 4 +- 9 files changed, 195 insertions(+), 101 deletions(-) diff --git a/examples/navier-stokes/Makefile b/examples/navier-stokes/Makefile index 834d274fd7..d25b1b9779 100644 --- a/examples/navier-stokes/Makefile +++ b/examples/navier-stokes/Makefile @@ -9,11 +9,12 @@ INSTALL_DIR = ${HOME}/local CXX = g++ CUDA_CXX = nvcc OMP_FLAGS = -DHAVE_OPENMP -fopenmp -#CXX_FLAGS = -std=gnu++0x -I$(TNL_INCLUDE_DIR) -O0 -g -DDEBUG $(OMP_FLAGS) -#CXX_FLAGS = -std=gnu++0x -I$(TNL_INCLUDE_DIR) -O3 $(OMP_FLAGS) -CXX_FLAGS = -std=gnu++0x -I$(TNL_INCLUDE_DIR) -O3 $(OMP_FLAGS) -pg -#CXX_FLAGS = -DHAVE_NOT_CXX11 -I$(TNL_INCLUDE_DIR) -O3 $(OMP_FLAGS) - +#CXX_FLAGS = -std=gnu++0x -I$(TNL_INCLUDE_DIR) -O0 -g -DDEBUG $(OMP_FLAGS) -DTEMPLATE_EXPLICIT_INSTANTIATION +CXX_FLAGS = -std=gnu++0x -I$(TNL_INCLUDE_DIR) -O3 $(OMP_FLAGS) -DTEMPLATE_EXPLICIT_INSTANTIATION +#CXX_FLAGS = -std=gnu++0x -I$(TNL_INCLUDE_DIR) -O3 $(OMP_FLAGS) -pg -DTEMPLATE_EXPLICIT_INSTANTIATION +#CXX_FLAGS = -DHAVE_NOT_CXX11 -I$(TNL_INCLUDE_DIR) -O3 $(OMP_FLAGS) -DTEMPLATE_EXPLICIT_INSTANTIATION + + LD_FLAGS = -L$(TNL_INSTALL_DIR) -ltnl-0.1 -lgomp SOURCES = main.cpp diff --git a/examples/navier-stokes/navierStokesSolver_impl.h b/examples/navier-stokes/navierStokesSolver_impl.h index bcd38a8cec..a2a0380866 100644 --- a/examples/navier-stokes/navierStokesSolver_impl.h +++ b/examples/navier-stokes/navierStokesSolver_impl.h @@ -111,7 +111,7 @@ bool navierStokesSolver< Mesh, EulerScheme > :: init( const tnlParameterContaine cerr << "Error: y-size must be positive integer number! It is " << meshes. y() << " now." << endl; return false; } - this -> mesh. setDimensions( meshes. y(), meshes. x() ); + this -> mesh. setDimensions( meshes. x(), meshes. y() ); RealType hx = this -> mesh. getParametricStep(). x(); RealType hy = this -> mesh. getParametricStep(). y(); mesh. save( tnlString( "mesh.tnl" ) ); @@ -176,7 +176,7 @@ bool navierStokesSolver< Mesh, EulerScheme > :: setInitialCondition( const tnlPa for( IndexType j = 0; j < ySize; j ++ ) for( IndexType i = 0; i < xSize; i ++ ) { - const IndexType c = mesh. getElementIndex( j, i ); + const IndexType c = mesh. getElementIndex( i, j ); const RealType x = i * hx; const RealType y = j * hy; @@ -269,7 +269,7 @@ void navierStokesSolver< Mesh, EulerScheme > :: updatePhysicalQuantities( const for( IndexType j = 0; j < ySize; j ++ ) for( IndexType i = 0; i < xSize; i ++ ) { - IndexType c = mesh. getElementIndex( j, i ); + IndexType c = mesh. getElementIndex( i, j ); u1[ c ] = rho_u1[ c ] / rho[ c ]; u2[ c ] = rho_u2[ c ] / rho[ c ]; p[ c ] = rho[ c ] * this -> R * this -> T; @@ -359,10 +359,10 @@ void navierStokesSolver< Mesh, EulerScheme > :: GetExplicitRHS( const RealType& */ for( IndexType i = 0; i < xSize; i ++ ) { - const IndexType c1 = mesh. getElementIndex( 0, i ); - const IndexType c2 = mesh. getElementIndex( 1, i ); - const IndexType c3 = mesh. getElementIndex( ySize - 1, i ); - const IndexType c4 = mesh. getElementIndex( ySize - 2, i ); + const IndexType c1 = mesh. getElementIndex( i, 0 ); + const IndexType c2 = mesh. getElementIndex( i, 1 ); + const IndexType c3 = mesh. getElementIndex( i, ySize - 1 ); + const IndexType c4 = mesh. getElementIndex( i, ySize - 2 ); RealType x = i * hx / mesh. getProportions(). x(); @@ -388,10 +388,10 @@ void navierStokesSolver< Mesh, EulerScheme > :: GetExplicitRHS( const RealType& } for( IndexType j = 0; j < ySize; j ++ ) { - const IndexType c1 = mesh. getElementIndex( j, 0 ); - const IndexType c2 = mesh. getElementIndex( j, 1 ); - const IndexType c3 = mesh. getElementIndex( j, xSize - 1 ); - const IndexType c4 = mesh. getElementIndex( j, xSize - 2 ); + const IndexType c1 = mesh. getElementIndex( 0, j ); + const IndexType c2 = mesh. getElementIndex( 1, j ); + const IndexType c3 = mesh. getElementIndex( xSize - 1, j ); + const IndexType c4 = mesh. getElementIndex( xSize - 2, j ); RealType y = j * hy / mesh. getProportions(). y(); @@ -422,7 +422,7 @@ void navierStokesSolver< Mesh, EulerScheme > :: GetExplicitRHS( const RealType& for( IndexType j = 0; j < ySize; j ++ ) for( IndexType i = 0; i < xSize; i ++ ) { - IndexType c = mesh. getElementIndex( j, i ); + IndexType c = mesh. getElementIndex( i, j ); if( i == 0 || j == 0 || i == xSize - 1 || j == ySize - 1 ) { diff --git a/src/implementation/mesh/tnlGrid2D_impl.h b/src/implementation/mesh/tnlGrid2D_impl.h index e963d163a0..828920379c 100644 --- a/src/implementation/mesh/tnlGrid2D_impl.h +++ b/src/implementation/mesh/tnlGrid2D_impl.h @@ -58,7 +58,7 @@ template< typename Real, typename Device, typename Index, template< int, typename, typename, typename > class Geometry > -void tnlGrid< 2, Real, Device, Index, Geometry > :: setDimensions( const Index ySize, const Index xSize ) +void tnlGrid< 2, Real, Device, Index, Geometry > :: setDimensions( const Index xSize, const Index ySize ) { tnlAssert( xSize > 0, cerr << "The number of Elements along x-axis must be larger than 0." ); @@ -80,7 +80,7 @@ template< typename Real, template< int, typename, typename, typename > class Geometry > void tnlGrid< 2, Real, Device, Index, Geometry > :: setDimensions( const tnlTuple< 2, Index >& dimensions ) { - this -> setDimensions( dimensions. y(), dimensions. x() ); + this -> setDimensions( dimensions. x(), dimensions. y() ); } template< typename Real, @@ -169,7 +169,7 @@ template< typename Real, typename Device, typename Index, template< int, typename, typename, typename > class Geometry > -Index tnlGrid< 2, Real, Device, Index, Geometry > :: getElementIndex( const Index j, const Index i ) const +Index tnlGrid< 2, Real, Device, Index, Geometry > :: getElementIndex( const Index i, const Index j ) const { tnlAssert( i < dimensions. x(), cerr << "Index i ( " << i @@ -183,6 +183,21 @@ Index tnlGrid< 2, Real, Device, Index, Geometry > :: getElementIndex( const Inde return j * this -> dimensions. x() + i; } +template< typename Real, + typename Device, + typename Index, + template< int, typename, typename, typename > class Geometry > +void tnlGrid< 2, Real, Device, Index, Geometry > :: getElementCoordinates( const Index element, + tnlTuple< 2, Index >& coordinates ) const +{ + tnlAssert( i >= 0 && i < dofs, + cerr << " i = " << i << " dofs = " << dofs + << " in tnlGrid " << this -> getName(); ); + + coordinates. x() = element % this -> dimensions. x(); + coordinates. y() = element / this -> dimensions. x(); +} + template< typename Real, typename Device, typename Index, @@ -201,9 +216,9 @@ Index tnlGrid< 2, Real, Device, Index, Geometry > :: getElementNeighbour( const template< typename Real, - typename Device, - typename Index, - template< int, typename, typename, typename > class Geometry > + typename Device, + typename Index, + template< int, typename, typename, typename > class Geometry > Index tnlGrid< 2, Real, Device, Index, Geometry > :: getDofs() const { return this -> dofs; @@ -213,16 +228,25 @@ template< typename Real, typename Device, typename Index, template< int, typename, typename, typename > class Geometry > -Real tnlGrid< 2, Real, Device, Index, Geometry > :: getElementMeasure( const Index j, - const Index i ) const +void tnlGrid< 2, Real, Device, Index, Geometry > :: getElementCenter( const tnlTuple< 2, Index >& coordinates, + tnlTuple< 2, Real >& center ) const { - return geometry. getElementMeasure( j, i ); + geometry. getElementCenter( origin, coordinates, center ); } template< typename Real, typename Device, typename Index, template< int, typename, typename, typename > class Geometry > +Real tnlGrid< 2, Real, Device, Index, Geometry > :: getElementMeasure( const tnlTuple< 2, Index >& coordinates ) const +{ + return geometry. getElementMeasure( coordinates ); +} + +/*template< typename Real, + typename Device, + typename Index, + template< int, typename, typename, typename > class Geometry > template< Index dy, Index dx > Real tnlGrid< 2, Real, Device, Index, Geometry > :: getElementsDistance( const Index j, const Index i ) const @@ -250,7 +274,7 @@ tnlTuple< 2, Real > tnlGrid< 2, Real, Device, Index, Geometry > :: getEdgeNormal const Index i ) const { return geometry. getEdgeNormal< dy, dx >( j, i ); -} +}*/ template< typename Real, typename Device, @@ -349,7 +373,7 @@ bool tnlGrid< 2, Real, Device, Index, Geometry > :: write( const MeshFunction& f { const RealType x = this -> getOrigin(). x() + i * hx; const RealType y = this -> getOrigin(). y() + j * hy; - file << x << " " << " " << y << " " << function[ this -> getElementIndex( j, i ) ] << endl; + file << x << " " << " " << y << " " << function[ this -> getElementIndex( i, j ) ] << endl; } file << endl; } diff --git a/src/implementation/mesh/tnlGrid3D_impl.h b/src/implementation/mesh/tnlGrid3D_impl.h index 91022ab6c3..bfdac5f7a6 100644 --- a/src/implementation/mesh/tnlGrid3D_impl.h +++ b/src/implementation/mesh/tnlGrid3D_impl.h @@ -55,7 +55,7 @@ template< typename Real, typename Device, typename Index, template< int, typename, typename, typename > class Geometry > -void tnlGrid< 3, Real, Device, Index, Geometry > :: setDimensions( const Index zSize, const Index ySize, const Index xSize ) +void tnlGrid< 3, Real, Device, Index, Geometry > :: setDimensions( const Index xSize, const Index ySize, const Index zSize ) { tnlAssert( xSize > 1, cerr << "The number of Elements along x-axis must be larger than 1." ); @@ -76,15 +76,7 @@ template< typename Real, template< int, typename, typename, typename > class Geometry > void tnlGrid< 3, Real, Device, Index, Geometry > :: setDimensions( const tnlTuple< 3, Index >& dimensions ) { - tnlAssert( dimensions. x() > 1, - cerr << "The number of Elements along x-axis must be larger than 1." ); - tnlAssert( dimensions. y() > 1, - cerr << "The number of Elements along y-axis must be larger than 1." ); - tnlAssert( dimensions. z() > 1, - cerr << "The number of Elements along z-axis must be larger than 1." ); - - this -> dimensions = dimensions; - dofs = this -> dimensions. x() * this -> dimensions. y() * this -> dimensions. z(); + this -> setDimensions( this -> dimensions. x(), this -> dimensions. y(), this -> dimensions. z() ); } template< typename Real, @@ -177,7 +169,7 @@ template< typename Real, typename Device, typename Index, template< int, typename, typename, typename > class Geometry > -Index tnlGrid< 3, Real, Device, Index, Geometry > :: getElementIndex( const Index k, const Index j, const Index i ) const +Index tnlGrid< 3, Real, Device, Index, Geometry > :: getElementIndex( const Index i, const Index j, const Index k ) const { tnlAssert( i < dimensions. x(), cerr << "Index i ( " << i @@ -297,7 +289,7 @@ bool tnlGrid< 3, Real, Device, Index, Geometry > :: write( const MeshFunction& f { const RealType x = this -> getOrigin(). x() + i * hx; const RealType y = this -> getOrigin(). y() + j * hy; - //file << x << " " << " " << y << " " << function[ this -> getElementIndex( j, i ) ] << endl; + //file << x << " " << " " << y << " " << function[ this -> getElementIndex( i, j ) ] << endl; } file << endl; } diff --git a/src/implementation/mesh/tnlIdenticalGridGeometry_impl.h b/src/implementation/mesh/tnlIdenticalGridGeometry_impl.h index db1396ef27..b052df5956 100644 --- a/src/implementation/mesh/tnlIdenticalGridGeometry_impl.h +++ b/src/implementation/mesh/tnlIdenticalGridGeometry_impl.h @@ -45,11 +45,11 @@ const tnlTuple< 1, Real >& tnlIdenticalGridGeometry< 1, Real, Device, Index > :: template< typename Real, typename Device, typename Index > -void tnlIdenticalGridGeometry< 1, Real, Device, Index > :: getElementCoordinates( const Index i, - const tnlTuple< 1, Real >& origin, - tnlTuple< 1, Real >& coordinates ) const +void tnlIdenticalGridGeometry< 1, Real, Device, Index > :: getElementCenter( const tnlTuple< 1, Real >& origin, + const tnlTuple< 1, Index >& coordinates, + tnlTuple< 1, Real >& center ) const { - coordinates. x() = ( i + 0.5 ) * parametricStep. x(); + center. x() = ( coordinates. x() + 0.5 ) * parametricStep. x(); } template< typename Real, @@ -145,20 +145,18 @@ const tnlTuple< 2, Real >& tnlIdenticalGridGeometry< 2, Real, Device, Index > :: template< typename Real, typename Device, typename Index > -void tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementCoordinates( const Index j, - const Index i, - const tnlTuple< 2, Real >& origin, - tnlTuple< 2, Real >& coordinates ) const +void tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementCenter( const tnlTuple< 2, Real >& origin, + const tnlTuple< 2, Index >& coordinates, + tnlTuple< 2, Real >& center ) const { - coordinates. x() = ( i + 0.5 ) * parametricStep. x(); - coordinates. y() = ( j + 0.5 ) * parametricStep. y(); + center. x() = ( coordinates. x() + 0.5 ) * parametricStep. x(); + center. y() = ( coordinates. y() + 0.5 ) * parametricStep. y(); } template< typename Real, typename Device, typename Index > -Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementMeasure( const Index j, - const Index i ) const +Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementMeasure( const tnlTuple< 2, Index >& coordinates ) const { return elementMeasure; } @@ -166,9 +164,9 @@ Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementMeasure( co template< typename Real, typename Device, typename Index > - template< Index dy, Index dx > -Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementsDistance( const Index j, - const Index i ) const + template< Index dx, Index dy > +Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementsDistance( const Index i, + const Index j ) const { if( dy == 0 && dx == 1 ) return parametricStep. x(); @@ -182,9 +180,9 @@ Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getElementsDistance( template< typename Real, typename Device, typename Index > -template< Index dy, Index dx > -void tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeCoordinates( const Index j, - const Index i, +template< Index dx, Index dy > +void tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeCoordinates( const Index i, + const Index j, const tnlTuple< 2, Real >& origin, tnlTuple< 2, Real >& coordinates ) const { @@ -195,9 +193,9 @@ void tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeCoordinates( c template< typename Real, typename Device, typename Index > -template< Index dy, Index dx > -Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeLength( const Index j, - const Index i ) const +template< Index dx, Index dy > +Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeLength( const Index i, + const Index j ) const { if( dy == 0 && dx == 1 ) return parametricStep. y(); @@ -209,9 +207,9 @@ Real tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeLength( const template< typename Real, typename Device, typename Index > -template< Index dy, Index dx > -void tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeNormal( const Index j, - const Index i, +template< Index dx, Index dy > +void tnlIdenticalGridGeometry< 2, Real, Device, Index > :: getEdgeNormal( const Index i, + const Index j, tnlTuple< 2, Real >& normal ) const { tnlAssert( ( dx == 0 || dx == 1 || dx == -1 || diff --git a/src/implementation/schemes/euler/fvm/tnlLaxFridrichs_impl.h b/src/implementation/schemes/euler/fvm/tnlLaxFridrichs_impl.h index a1cae1b79a..573aa1a098 100644 --- a/src/implementation/schemes/euler/fvm/tnlLaxFridrichs_impl.h +++ b/src/implementation/schemes/euler/fvm/tnlLaxFridrichs_impl.h @@ -107,9 +107,9 @@ template< typename Real, typename Index, typename PressureGradient > void tnlLaxFridrichs< tnlGrid< 2, Real, Device, Index >, PressureGradient > :: getExplicitRhs( const IndexType centralVolume, - RealType& rho_t, - RealType& rho_u1_t, - RealType& rho_u2_t ) const + RealType& rho_t, + RealType& rho_u1_t, + RealType& rho_u2_t ) const { tnlAssert( mesh, cerr << "No mesh has been binded with the Lax-Fridrichs scheme." ); tnlAssert( pressureGradient, cerr << "No pressure gradient was set in the the Lax-Fridrichs scheme." ) @@ -156,6 +156,81 @@ void tnlLaxFridrichs< tnlGrid< 2, Real, Device, Index >, PressureGradient > :: - ( rho_u2[ e ] * u1_e - rho_u2[ w ] * u1_w ) / ( 2.0 * hx ) - ( rho_u2[ n ] * u2_n - rho_u2[ s ] * u2_s ) / ( 2.0 * hy ) - p_y; + + /**** + * Scheme for deformed grids + */ + const RealType u1_c = rho_u1[ c ] / regularize( rho[ c ] ); + const RealType u1_n = rho_u1[ n ] / regularize( rho[ n ] ); + const RealType u1_s = rho_u1[ s ] / regularize( rho[ s ] ); + const RealType u2_c = rho_u2[ c ] / regularize( rho[ c ] ); + const RealType u2_e = rho_u2[ e ] / regularize( rho[ e ] ); + const RealType u2_w = rho_u2[ w ] / regularize( rho[ w ] ); + + + /**** + * Get the central volume and its neighbours (east, north, west, south) coordinates + */ + tnlTuple< 2, IndexType > c_coordinates, e_coordinates, n_coordinates, w_coordinates, s_coordinates; + this -> mesh -> getElementCoordinates( c, c_coordinates ); + e_coordinates = n_coordinates = w_coordinates = s_coordinates = c_coordinates; + e_coordinates. x() ++; + w_coordinates. x() --; + n_coordinates. y() ++; + s_coordinates. y() --; + + /**** + * Get the volumes measure + */ + const RealType mu_D_c = this -> mesh -> getElementMeasure( c_coordinates ); + const RealType mu_D_e = this -> mesh -> getElementMeasure( e_coordinates ); + const RealType mu_D_n = this -> mesh -> getElementMeasure( n_coordinates ); + const RealType mu_D_w = this -> mesh -> getElementMeasure( w_coordinates ); + const RealType mu_D_s = this -> mesh -> getElementMeasure( s_coordinates ); + + /**** + * Get the volumes centers of gravity + */ + tnlTuple< 2, RealType > c_center, e_center, w_center, n_center, s_center; + this -> mesh -> getElementCenter( c_coordinates, c_center ); + this -> mesh -> getElementCenter( c_coordinates, e_center ); + this -> mesh -> getElementCenter( c_coordinates, w_center ); + this -> mesh -> getElementCenter( c_coordinates, n_center ); + this -> mesh -> getElementCenter( c_coordinates, s_center ); + + /**** + * Get delta x and delta y between the volumes + */ + const RealType dx_e = e_center. x() - c_center. x(); + const RealType dx_w = w_center. x() - w_center. x(); + const RealType dx_n = n_center. x() - n_center. x(); + const RealType dx_s = s_center. x() - s_center. x(); + const RealType dy_e = e_center. y() - c_center. y(); + const RealType dy_w = w_center. y() - w_center. y(); + const RealType dy_n = n_center. y() - n_center. y(); + const RealType dy_s = s_center. y() - s_center. y(); + + /**** + * Compute the fluxes + */ + const RealType rho_f_e = 0.5 * ( rho[ c ] * u1_c + rho[ e ] * u1_e ); + const RealType rho_f_w = 0.5 * ( rho[ c ] * u1_c + rho[ w ] * u1_w ); + const RealType rho_f_n = 0.5 * ( rho[ c ] * u1_c + rho[ n ] * u1_n ); + const RealType rho_f_s = 0.5 * ( rho[ c ] * u1_c + rho[ s ] * u1_s ); + const RealType rho_g_e = 0.5 * ( rho[ c ] * u2_c + rho[ e ] * u2_e ); + const RealType rho_g_w = 0.5 * ( rho[ c ] * u2_c + rho[ w ] * u2_w ); + const RealType rho_g_n = 0.5 * ( rho[ c ] * u2_c + rho[ n ] * u2_n ); + const RealType rho_g_s = 0.5 * ( rho[ c ] * u2_c + rho[ s ] * u2_s ); + + rho_t = - 1.0 / mu_D_c * ( rho_f_e * dy_e - rho_g_e * dx_e + + rho_f_n * dy_n - rho_g_n * dx_n + + rho_f_w * dy_w - rho_g_w * dx_w + + rho_f_s * dy_s - rho_g_s * dx_s ) + + 1.0 / ( 8.0 * mu_D_c ) * + ( ( mu_D_c + mu_D_e ) * ( rho[ e ] - rho[ c ] ) + + ( mu_D_c + mu_D_n ) * ( rho[ n ] - rho[ c ] ) + + ( mu_D_c + mu_D_w ) * ( rho[ w ] - rho[ c ] ) + + ( mu_D_c + mu_D_s ) * ( rho[ s ] - rho[ c ] ) ); } template< typename Real, diff --git a/src/mesh/tnlGrid.h b/src/mesh/tnlGrid.h index ed1a2444b9..c0bd40569e 100644 --- a/src/mesh/tnlGrid.h +++ b/src/mesh/tnlGrid.h @@ -121,7 +121,7 @@ class tnlGrid< 2, Real, Device, Index, Geometry > : public tnlObject tnlString getType() const; - void setDimensions( const Index ySize, const Index xSize ); + void setDimensions( const Index xSize, const Index ySize ); void setDimensions( const tnlTuple< 2, Index >& ); @@ -139,8 +139,11 @@ class tnlGrid< 2, Real, Device, Index, Geometry > : public tnlObject const tnlTuple< 2, Real >& getParametricStep() const; - Index getElementIndex( const Index j, - const Index i ) const; + Index getElementIndex( const Index i, + const Index j ) const; + + void getElementCoordinates( const Index i, + tnlTuple< 2, Index >& coordinates ) const; Index getElementNeighbour( const Index Element, const Index dy, @@ -148,10 +151,12 @@ class tnlGrid< 2, Real, Device, Index, Geometry > : public tnlObject Index getDofs() const; - Real getElementMeasure( const Index j, - const Index i ) const; + void getElementCenter( const tnlTuple< 2, Index >& coordinates, + tnlTuple< 2, Real >& center ) const; + + Real getElementMeasure( const tnlTuple< 2, Index >& coordinates ) const; - template< Index dy, Index dx > + /*template< Index dy, Index dx > Real getElementsDistance( const Index j, const Index i ) const; @@ -161,7 +166,7 @@ class tnlGrid< 2, Real, Device, Index, Geometry > : public tnlObject template< int dy, int dx > tnlTuple< 2, Real > getEdgeNormal( const Index j, - const Index i ) const; + const Index i ) const;*/ //! Method for saving the object to a file as a binary data @@ -210,7 +215,7 @@ class tnlGrid< 3, Real, Device, Index, Geometry > : public tnlObject tnlString getType() const; - void setDimensions( const Index zSize, const Index ySize, const Index xSize ); + void setDimensions( const Index xSize, const Index ySize, const Index zSize ); void setDimensions( const tnlTuple< 3, Index >& ); @@ -228,7 +233,7 @@ class tnlGrid< 3, Real, Device, Index, Geometry > : public tnlObject tnlTuple< 3, Real > getParametricStep() const; - Index getElementIndex( const Index k, const Index j, const Index i ) const; + Index getElementIndex( const Index i, const Index j, const Index k ) const; Index getDofs() const; diff --git a/src/mesh/tnlIdenticalGridGeometry.h b/src/mesh/tnlIdenticalGridGeometry.h index a3d84f94e0..0fabcbf5ea 100644 --- a/src/mesh/tnlIdenticalGridGeometry.h +++ b/src/mesh/tnlIdenticalGridGeometry.h @@ -43,9 +43,10 @@ class tnlIdenticalGridGeometry< 1, Real, Device, Index > const tnlTuple< 1, Real >& getParametricStep() const; - void getElementCoordinates( const Index i, - const tnlTuple< 1, Real >& origin, - tnlTuple< 1, Real >& coordinates ) const; + void getElementCenter( const tnlTuple< 1, Real >& origin, + const tnlTuple< 1, Index >& coordinates, + tnlTuple< 1, Real >& center ) const; + Real getElementMeasure( const Index i ) const; @@ -94,35 +95,33 @@ class tnlIdenticalGridGeometry< 2, Real, Device, Index > const tnlTuple< 2, Real >& getParametricStep() const; - void getElementCoordinates( const Index j, - const Index i, - const tnlTuple< 2, Real >& origin, - tnlTuple< 2, Real >& coordinates ) const; + void getElementCenter( const tnlTuple< 2, Real >& origin, + const tnlTuple< 2, Index >& coordinates, + tnlTuple< 2, Real >& center ) const; - Real getElementMeasure( const Index j, - const Index i ) const; + Real getElementMeasure( const tnlTuple< 2, Index >& coordinates ) const; - template< Index dy, Index dx > - Real getElementsDistance( const Index j, - const Index i ) const; + template< Index dx, Index dy > + Real getElementsDistance( const Index i, + const Index j ) const; - template< Index dy, Index dx > - void getEdgeCoordinates( const Index j, - const Index i, + template< Index dx, Index dy > + void getEdgeCoordinates( const Index i, + const Index j, const tnlTuple< 2, Real >& origin, tnlTuple< 2, Real >& coordinates ) const; - template< Index dy, Index dx > - Real getEdgeLength( const Index j, - const Index i ) const; + template< Index dx, Index dy > + Real getEdgeLength( const Index i, + const Index j ) const; - template< Index dy, Index dx > - void getEdgeNormal( const Index j, - const Index i, + template< Index dx, Index dy > + void getEdgeNormal( const Index i, + const Index j, tnlTuple< 2, Real >& normal ) const; - void getVertexCoordinates( const Index j, - const Index i, + void getVertexCoordinates( const Index i, + const Index j, const tnlTuple< 2, Real >& origin, tnlTuple< 2, Real >& coordinates ) const; diff --git a/tests/unit-tests/core/tnlArrayTest.cpp b/tests/unit-tests/core/tnlArrayTest.cpp index 7dc6f8e5aa..9782352a30 100644 --- a/tests/unit-tests/core/tnlArrayTest.cpp +++ b/tests/unit-tests/core/tnlArrayTest.cpp @@ -21,7 +21,7 @@ int main( int argc, char* argv[] ) { #ifdef HAVE_CPPUNIT - if( ! tnlUnitTestStarter :: run< tnlArrayTester< char, tnlHost, int > >() || +/* if( ! tnlUnitTestStarter :: run< tnlArrayTester< char, tnlHost, int > >() || ! tnlUnitTestStarter :: run< tnlArrayTester< int, tnlHost, int > >() || ! tnlUnitTestStarter :: run< tnlArrayTester< long int, tnlHost, int > >() || ! tnlUnitTestStarter :: run< tnlArrayTester< float, tnlHost, int > >() || @@ -33,7 +33,7 @@ int main( int argc, char* argv[] ) ! tnlUnitTestStarter :: run< tnlArrayTester< float, tnlHost, long int > >() || ! tnlUnitTestStarter :: run< tnlArrayTester< double, tnlHost, long int > >() || ! tnlUnitTestStarter :: run< tnlArrayTester< long double, tnlHost, long int > >() ) - return EXIT_FAILURE; + return EXIT_FAILURE;*/ return EXIT_SUCCESS; #else return EXIT_FAILURE; -- GitLab