View Issue Details

IDProjectCategoryView StatusLast Update
0001319OpenFOAMBugpublic2014-06-06 16:09
Reporterdkxls Assigned Tohenry  
PrioritynormalSeverityfeatureReproducibilityalways
Status resolvedResolutionfixed 
PlatformLinux x86_64OSopenSUSEOS Version12.3
Summary0001319: [blockMesh]: B-splines not selectable for edges
DescriptionIt would be a nice little feature to have B-splines for edges in blockMesh.

The basic code for B-splines is already available in src/mesh/blockMesh/curvedEdges/BSpline.H, the missing piece is a class similar to splineEdge for Catmull-Rom splines.

The implementation is trivial, but to make things easier, I attached two files that implement this feature.
TagsNo tags attached.

Activities

dkxls

2014-06-05 19:14

reporter  

bsplineEdge.H (2,987 bytes)   
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Class
    Foam::bsplineEdge

Description
    A curvedEdge interface for B-splines.

SourceFiles
    bsplineEdge.C

\*---------------------------------------------------------------------------*/

#ifndef bsplineEdge_H
#define bsplineEdge_H

#include "curvedEdge.H"
#include "BSpline.H"

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{

/*---------------------------------------------------------------------------*\
                      Class bsplineEdge Declaration
\*---------------------------------------------------------------------------*/

class bsplineEdge
:
    public curvedEdge,
    public BSpline
{
    // Private Member Functions

        //- Disallow default bitwise copy construct
        bsplineEdge(const bsplineEdge&);

        //- Disallow default bitwise assignment
        void operator=(const bsplineEdge&);


public:

    //- Runtime type information
    TypeName("bspline");


    // Constructors

        //- Construct from components
        bsplineEdge
        (
            const pointField&,
            const label start,
            const label end,
            const pointField& internalPoints
        );

        //- Construct from Istream, setting pointsList
        bsplineEdge(const pointField&, Istream&);


    //- Destructor
    virtual ~bsplineEdge();


    // Member Functions

        //- Return the point position corresponding to the curve parameter
        //  0 <= lambda <= 1
        virtual point position(const scalar) const;

        //- Return the length of the spline curve (not implemented)
        virtual scalar length() const;
};


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

} // End namespace Foam

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#endif

// ************************************************************************* //
bsplineEdge.H (2,987 bytes)   

dkxls

2014-06-05 19:14

reporter  

bsplineEdge.C (2,637 bytes)   
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

#include "bsplineEdge.H"
#include "addToRunTimeSelectionTable.H"


// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

namespace Foam
{
    defineTypeNameAndDebug(bsplineEdge, 0);

    addToRunTimeSelectionTable
    (
        curvedEdge,
        bsplineEdge,
        Istream
    );
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::bsplineEdge::bsplineEdge
(
    const pointField& points,
    const label start,
    const label end,
    const pointField& internalPoints
)
:
    curvedEdge(points, start, end),
    BSpline(appendEndPoints(points, start, end, internalPoints))
{}


Foam::bsplineEdge::bsplineEdge(const pointField& points, Istream& is)
:
    curvedEdge(points, is),
    BSpline(appendEndPoints(points, start_, end_, pointField(is)))
{
    token t(is);
    is.putBack(t);

    // discard unused start/end tangents
    if (t == token::BEGIN_LIST)
    {
        vector tangent0Ignored(is);
        vector tangent1Ignored(is);
    }
}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

Foam::bsplineEdge::~bsplineEdge()
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

Foam::point Foam::bsplineEdge::position(const scalar mu) const
{
    return BSpline::position(mu);
}


Foam::scalar Foam::bsplineEdge::length() const
{
    return BSpline::length();
}


// ************************************************************************* //
bsplineEdge.C (2,637 bytes)   

henry

2014-06-06 16:09

manager   ~0003111

Resolved by commit fc54b9a40f91c9e6ddebff284df53b42f8a239a9

Issue History

Date Modified Username Field Change
2014-06-05 19:14 dkxls New Issue
2014-06-05 19:14 dkxls File Added: bsplineEdge.H
2014-06-05 19:14 dkxls File Added: bsplineEdge.C
2014-06-06 16:09 henry Note Added: 0003111
2014-06-06 16:09 henry Status new => resolved
2014-06-06 16:09 henry Resolution open => fixed
2014-06-06 16:09 henry Assigned To => henry