Non classé

OpenFOAM – fixedTabulatedScalarValue

Cette condition limite permet d’imposer au bord du domaine fluide des valeurs scalaires stockées dans une table. Les valeurs scalaires doivent être associées à un point définit par ses coordonnées.

xyzepsilon
-12-4200.00800494
-12-420.2740780.00730994
-12-41.60640.2740780.00792009

Le programme utilise un algorithme des plus proches voisins pour associer à chaque centre de faces du patch sur lequel est imposée la condition limite une valeur d’entrée calculée par la méthode de la pondération inverse à la distance.

Les développements font appel à la librarie nanoflann accessible en cliquant sur le lien suivant. Cette librarie ne requiert aucune compilation préalable, toutefois les liens vers les headers doivent être correctement définis pour pouvoir compiler les sources de fixedTabularScalarValue.

Dans la suite, nous fournissons les sources du .h et du .cpp ainsi que des deux fichiers options et files nécessaires à la compilation. Les chemins sont susceptibles de changer en fonction de l’environnement de travail.

 

  • options
EXE_INC = \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(FOAM_DEV)/libraries/nanoflann/include \
-I$(FOAM_DEV)/libraries/utils/include
EXE_LIBS =<span id="mce_marker" data-mce-type="bookmark" data-mce-fragment="1">​</span>

 

  • files
fvPatchFields = ./
derivedFvPatchFields = $(fvPatchFields)/derived

$(derivedFvPatchFields)/fixedTabulatedScalarValue/fixedTabulatedScalarValueFvPatchScalarField.C

LIB = $(FOAM_USER_LIBBIN)/libArepCustomPatch

 

  • fixedTabulatedScalarValueFvPatchScalarField.H
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
 * Copyright 2018 arep <alexis.sauvageon@arep.fr>
 *
 * Redistribution and use in source and binary forms, with or without
 * modification, are permitted provided that the following conditions are
 * met:
 * 
 * * Redistributions of source code must retain the above copyright
 *   notice, this list of conditions and the following disclaimer.
 * * Redistributions in binary form must reproduce the above
 *   copyright notice, this list of conditions and the following disclaimer
 *   in the documentation and/or other materials provided with the
 *   distribution.
 * * Neither the name of the  nor the names of its
 *   contributors may be used to endorse or promote products derived from
 *   this software without specific prior written permission.
 * 
 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
 * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
 * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
 * A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
 * OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
 * SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
 * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
 * DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
 * THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
 * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
 * OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.

Class
    fixedTabulatedScalarValueFvPatchScalarField

Description
    Compute boundaries conditions from a field stored in table (format .raw)
	[X Y Z Scalar]
    the two first lines of the table are not interpreted.
    The interpollated value is computed thanks to a PID method.

    \heading Patch usage

    \table
        Property    | Description             	 	| Required   	| Default value
        fileDir     | path to the input file  		| yes        	| 
        nbPoint     | number of point for PID 	        | no       	| 4
        wFactor     | p factor for PID      		| no        	| 2.0

    \endtable

    Example of the boundary condition specification:
    \verbatim
    myPatch
    {
        type            fixedTabulatedScalarValue;
        fileDir         "constant/inputTableData/p_boxValue.raw";
        nbPoint              8;
        wFactor            2.0;

    }
    \endverbatim


SourceFiles
    fixedTabulatedScalarValueFvPatchScalarField.C

Author
    Alexis Sauvageon.  All rights reserved

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

#ifndef fixedTabulatedScalarValueFvPatchScalarField_H
#define fixedTabulatedScalarValueFvPatchScalarField_H

#include "fvPatchFields.H"
#include "fixedValueFvPatchFields.H"


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

namespace Foam
{

/*---------------------------------------------------------------------------*\
              Class fixedTabulatedScalarValueFvPatchField Declaration
\*---------------------------------------------------------------------------*/

class fixedTabulatedScalarValueFvPatchScalarField
:
    public fixedValueFvPatchScalarField
{
    // Private data

        //- Peak velocity magnitude
        std::string fileDir_;
	unsigned int nbPoint_;
	scalar wFactor_;

public:

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


    // Constructors

        //- Construct from patch and internal field
        fixedTabulatedScalarValueFvPatchScalarField
        (
            const fvPatch&,
            const DimensionedField<scalar, volMesh>&
        );

        //- Construct from patch, internal field and dictionary
        fixedTabulatedScalarValueFvPatchScalarField
        (
            const fvPatch&,
            const DimensionedField<scalar, volMesh>&,
            const dictionary&
        );

        //- Construct by mapping given fixedTabulatedScalarValueFvPatchScalarField
        //  onto a new patch
        fixedTabulatedScalarValueFvPatchScalarField
        (
            const fixedTabulatedScalarValueFvPatchScalarField&,
            const fvPatch&,
            const DimensionedField<scalar, volMesh>&,
            const fvPatchFieldMapper&
        );

        //- Construct and return a clone
        virtual tmp<fvPatchScalarField> clone() const
        {
            return tmp<fvPatchScalarField>
            (
                new fixedTabulatedScalarValueFvPatchScalarField(*this)
            );
        }

        //- Construct as copy setting internal field reference
        fixedTabulatedScalarValueFvPatchScalarField
        (
            const fixedTabulatedScalarValueFvPatchScalarField&,
            const DimensionedField<scalar, volMesh>&
        );

        //- Construct and return a clone setting internal field reference
        virtual tmp<fvPatchScalarField> clone
        (
            const DimensionedField<scalar, volMesh>& iF
        ) const
        {
            return tmp<fvPatchScalarField>
            (
                new fixedTabulatedScalarValueFvPatchScalarField(*this, iF)
            );
        }


    // Member functions

        //- Return fileDir_
        std::string& fileDir()
        {
            return fileDir_;
        }
        //- Return nbPoint_
        unsigned int& nbPoint()
        {
            return nbPoint_;
        }
        //- Return wFactor_
        scalar& wFactor()
        {
            return wFactor_;
        }


        //- Update coefficients
        virtual void updateCoeffs();

        //- Write
        virtual void write(Ostream&) const;
};


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

} // End namespace Foam

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

#endif

// ************************************************************************* //
<span id="mce_marker" data-mce-type="bookmark" data-mce-fragment="1">​</span>

 

  • fixedTabulatedScalarValueFvPatchScalarField.C
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | foam-extend: Open Source CFD
   \\    /   O peration     | Version:     4.0
    \\  /    A nd           | Web:         http://www.foam-extend.org
     \\/     M anipulation  | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
 * Copyright 2018 arep <alexis.sauvageon@arep.fr>
 *
 * Redistribution and use in source and binary forms, with or without
 * modification, are permitted provided that the following conditions are
 * met:
 * 
 * * Redistributions of source code must retain the above copyright
 *   notice, this list of conditions and the following disclaimer.
 * * Redistributions in binary form must reproduce the above
 *   copyright notice, this list of conditions and the following disclaimer
 *   in the documentation and/or other materials provided with the
 *   distribution.
 * * Neither the name of the  nor the names of its
 *   contributors may be used to endorse or promote products derived from
 *   this software without specific prior written permission.
 * 
 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
 * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
 * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
 * A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
 * OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
 * SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
 * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
 * DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
 * THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
 * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
 * OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.

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

#include "fixedTabulatedScalarValueFvPatchScalarField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "surfaceFields.H"
#include <fstream>
#include <sstream>
#include <iterator>

#include "nanoflann.hpp"
#include "kdtree.h"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //


namespace Foam
{

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

fixedTabulatedScalarValueFvPatchScalarField::fixedTabulatedScalarValueFvPatchScalarField
(
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF
)
:
    fixedValueFvPatchScalarField(p, iF),
    fileDir_("."),
    nbPoint_(4),
    wFactor_(2)
{}


fixedTabulatedScalarValueFvPatchScalarField::fixedTabulatedScalarValueFvPatchScalarField
(
    const fixedTabulatedScalarValueFvPatchScalarField& ptf,
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF,
    const fvPatchFieldMapper& mapper
)
:
    fixedValueFvPatchScalarField(ptf, p, iF, mapper),
    fileDir_(ptf.fileDir_),
    nbPoint_(ptf.nbPoint_),
    wFactor_(ptf.wFactor_)
{}


fixedTabulatedScalarValueFvPatchScalarField::fixedTabulatedScalarValueFvPatchScalarField
(
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF,
    const dictionary& dict
)
:
    fixedValueFvPatchScalarField(p, iF),
    fileDir_(string(dict.lookup("fileDir"))),
    nbPoint_(dict.lookupOrDefault<scalar>("nbPoint",4)),
    wFactor_(dict.lookupOrDefault<scalar>("wFactor",2))
{

// check the existance of the file

    std::ifstream inFile(fileDir_,std::ios::in);  
 
    if(inFile) 
    {       
        inFile.close(); 
    }
    else 
    {
        FatalErrorIn("fixedTabulatedScalarValueFvPatchScalarField(dict)")
            << "path to the raw file is not correct"
            << abort(FatalError);
    }

    if( (nbPoint_ < 1) || (wFactor_ <= 0) ) 
    {
        FatalErrorIn("fixedTabulatedScalarValueFvPatchScalarField(dict)")
            << "the number of interpolation point is less than 1 or the weighting factor is not positive."
            << abort(FatalError);
    }

    evaluate();
}


fixedTabulatedScalarValueFvPatchScalarField::fixedTabulatedScalarValueFvPatchScalarField
(
    const fixedTabulatedScalarValueFvPatchScalarField& fcvpvf,
    const DimensionedField<scalar, volMesh>& iF
)
:
    fixedValueFvPatchScalarField(fcvpvf, iF),
    fileDir_(fcvpvf.fileDir_),
    nbPoint_(fcvpvf.nbPoint_),
    wFactor_(fcvpvf.wFactor_)
{}


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

void fixedTabulatedScalarValueFvPatchScalarField::updateCoeffs()
{
    if (updated())
    {
        return;
    }
// >> open the input data file
	std::ifstream inFile(fileDir_,std::ios::in);
        // count the lines in the file
	inFile.unsetf(std::ios_base::skipws);
    	unsigned nbLines = std::count(
        	std::istream_iterator<char>(inFile),
        	std::istream_iterator<char>(), 
        	'\n');
	nbLines=nbLines-2;
	//...skip two first lines (commentary)
	inFile.close();
	inFile.open(fileDir_);
	std::string line;  
        std::getline(inFile, line);  
        std::getline(inFile, line);  

	//Initialize a plain continous array for the data
        scalar vectorData[nbLines][4];
	unsigned int index_(0);

	PointCloud<scalar> cloud;
	//...for all line until the end, create a scalarField with the computed point
        while(std::getline(inFile, line))  // tant que l'on peut mettre la ligne dans "contenu"
        {
		std::istringstream buffer(line);
		std::string word_("");

	        buffer >> word_;
		scalar x_(atof(word_.c_str()));
	        buffer >> word_;
		scalar y_(atof(word_.c_str()));
	        buffer >> word_;
		scalar z_(atof(word_.c_str()));
	        buffer >> word_;
		scalar val_(atof(word_.c_str()));

		updatePointCloud(cloud, x_, y_, z_); 

		vectorData[index_][0] = x_;
                vectorData[index_][1] = y_;
                vectorData[index_][2] = z_;
                vectorData[index_][3] = val_;
		index_++;
	
        }

	
	// construct a kd-tree index:
	typedef nanoflann::KDTreeSingleIndexAdaptor<
		nanoflann::L2_Simple_Adaptor<scalar, PointCloud<scalar> > ,
		PointCloud<scalar>,
		3
		> my_kd_tree_t;

	my_kd_tree_t   index(3, cloud, nanoflann::KDTreeSingleIndexAdaptorParams(10) );
	index.buildIndex();

#if 0
	// Test resize of dataset and rebuild of index:
	cloud.pts.resize(cloud.pts.size()*0.5);
	index.buildIndex();
#endif


// >> get the face of the current patch & initialize a scalarField of the right size
	const vectorField& faceCenter = patch().Cf();
	scalarField values(faceCenter.size(), scalar(0));

// >> compute the scalar value for all face centers
        forAll(faceCenter, f)
        {
// >>>>>>> get face center point
                const point faceCenterPoint = faceCenter[f];
// >>>>>>> find closest neighbourg
		// define the requested point
		const scalar query_pt[3] = { faceCenterPoint.x(), faceCenterPoint.y(), faceCenterPoint.z()};
		// find the 8 closest points
		size_t num_results = nbPoint_;
		std::vector<size_t>   ret_index(num_results);
		std::vector<scalar> out_dist_sqr(num_results);

		num_results = index.knnSearch(&query_pt[0], num_results, &ret_index[0], &out_dist_sqr[0]);
		
		// In case of less points in the tree than requested:
		ret_index.resize(num_results);
		out_dist_sqr.resize(num_results);

		// weight factors for PID
		std::vector<scalar> out_dist_weight(num_results);
		for (size_t i = 0; i < num_results; i++)
		{
			out_dist_weight[i]=1/std::pow(out_dist_sqr[i],wFactor_);
		}	
		// Distance inverse weighting methode (PID)
		scalar WxScalar(0);
		scalar Wtot(0);
		for (size_t i = 0; i < num_results; i++)
		{
			WxScalar+=out_dist_weight[i]*vectorData[ret_index[i]][3];
			Wtot+=out_dist_weight[i];
		}
		values[f]=WxScalar/Wtot;

		// print the index of the closest point and their distance
		/*cout << faceCenterPoint.x() << "    " << faceCenterPoint.y() << "    " << faceCenterPoint.z()  << std::endl;
		cout << values[f] << std::endl;
		cout << std::endl;

		for (size_t i = 0; i < num_results; i++)
		{
		cout << vectorData[ret_index[i]][0]  << "    " << vectorData[ret_index[i]][1] << "    " << vectorData[ret_index[i]][2] << std::endl;
		cout << "dist : " << out_dist_sqr[i] << " value : " << vectorData[ret_index[i]][4]  << std::endl;
		cout << "weight : " << out_dist_weight[i]   << std::endl;
		}*/

        }



    // Calculate local 1-D coordinate for the parabolic profile
   // scalarField coord = 0; //2*((c - ctr) & y_)/((bb.max() - bb.min()) & y_);

   scalarField::operator=(values);
}


// Write
void fixedTabulatedScalarValueFvPatchScalarField::write(Ostream& os) const
{
    fvPatchScalarField::write(os);
    os.writeKeyword("fileDir")
        << fileDir_ << token::END_STATEMENT << nl;
    writeEntry("value", os);
}


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

makePatchTypeField(fvPatchScalarField, fixedTabulatedScalarValueFvPatchScalarField);

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

} // End namespace Foam

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