Skip to content
Snippets Groups Projects
Commit 36fc8595 authored by Henry Weller's avatar Henry Weller
Browse files

reactionThermo/functionObjects/moleFractions/moleFractions: New functionObject...

reactionThermo/functionObjects/moleFractions/moleFractions: New functionObject to calculate and write mole-fraction fields

    This function object calculates mole-fraction fields from the mass-fraction
    fields of the psi/rhoReactionThermo and caches them for output and further
    post-processing.

    The names of the mole-fraction fields are obtained from the corresponding
    mass-fraction fields prepended by "X_"

    Example of function object specification:

        moleFractions
        {
            type psiReactionThermoMoleFractions;
        }

     or

        moleFractions
        {
            type rhoReactionThermoMoleFractions;
        }

    depending on the thermodynamics package used in the solver.
parent a6ae4f19
Branches
Tags
No related merge requests found
......@@ -19,4 +19,6 @@ derivedFvPatchFields/fixedUnburntEnthalpy/fixedUnburntEnthalpyFvPatchScalarField
derivedFvPatchFields/gradientUnburntEnthalpy/gradientUnburntEnthalpyFvPatchScalarField.C
derivedFvPatchFields/mixedUnburntEnthalpy/mixedUnburntEnthalpyFvPatchScalarField.C
functionObjects/moleFractions/moleFractionsFunctionObjects.C
LIB = $(FOAM_LIBBIN)/libreactionThermophysicalModels
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 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 "moleFractions.H"
#include "basicThermo.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class ThermoType>
void Foam::moleFractions<ThermoType>::calculateMoleFractions()
{
const ThermoType& thermo =
mesh_.lookupObject<ThermoType>(basicThermo::dictName);
const PtrList<volScalarField>& Y = thermo.composition().Y();
const volScalarField W(thermo.composition().W());
forAll(Y, i)
{
X_[i] = W*Y[i]/thermo.composition().W(i);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class ThermoType>
Foam::moleFractions<ThermoType>::moleFractions
(
const word& name,
const objectRegistry& obr,
const dictionary& dict,
const bool loadFromFiles
)
:
functionObjectFile(obr, name),
name_(name),
mesh_(refCast<const fvMesh>(obr))
{
if (mesh_.foundObject<ThermoType>(basicThermo::dictName))
{
const ThermoType& thermo =
mesh_.lookupObject<ThermoType>(basicThermo::dictName);
const PtrList<volScalarField>& Y = thermo.composition().Y();
X_.setSize(Y.size());
forAll(Y, i)
{
X_.set
(
i,
new volScalarField
(
IOobject
(
"X_" + Y[i].name(),
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh_,
dimensionedScalar("X", dimless, 0)
)
);
}
calculateMoleFractions();
}
else
{
FatalErrorInFunction
<< "Cannot find thermodynamics model of type "
<< ThermoType::typeName
<< exit(FatalError);
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template<class ThermoType>
Foam::moleFractions<ThermoType>::~moleFractions()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class ThermoType>
void Foam::moleFractions<ThermoType>::read
(
const dictionary& dict
)
{}
template<class ThermoType>
void Foam::moleFractions<ThermoType>::execute()
{
calculateMoleFractions();
}
template<class ThermoType>
void Foam::moleFractions<ThermoType>::end()
{}
template<class ThermoType>
void Foam::moleFractions<ThermoType>::timeSet()
{}
template<class ThermoType>
void Foam::moleFractions<ThermoType>::write()
{}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 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::moleFractions
Description
This function object calculates mole-fraction fields from the mass-fraction
fields of the psi/rhoReactionThermo and caches them for output and further
post-processing.
The names of the mole-fraction fields are obtained from the corresponding
mass-fraction fields prepended by "X_"
Example of function object specification:
\verbatim
moleFractions
{
type psiReactionThermoMoleFractions;
}
\endverbatim
or
\verbatim
moleFractions
{
type rhoReactionThermoMoleFractions;
}
\endverbatim
depending on the thermodynamics package used in the solver.
SeeAlso
Foam::functionObject
Foam::OutputFilterFunctionObject
SourceFiles
moleFractions.C
\*---------------------------------------------------------------------------*/
#ifndef moleFractions_H
#define moleFractions_H
#include "functionObjectFile.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class moleFractions Declaration
\*---------------------------------------------------------------------------*/
template<class ThermoType>
class moleFractions
:
public functionObjectFile
{
// Private data
//- Name of moleFractions functionObject
word name_;
//- Reference to the mesh
const fvMesh& mesh_;
//- Species mole fractions
PtrList<volScalarField> X_;
// Private Member Functions
//- Calculate the mole fraction fields
virtual void calculateMoleFractions();
//- Disallow default bitwise copy construct
moleFractions(const moleFractions&);
//- Disallow default bitwise assignment
void operator=(const moleFractions&);
public:
//- Runtime type information
TypeName("moleFractions");
// Constructors
//- Construct for given objectRegistry and dictionary.
// Allow the possibility to load fields from files
moleFractions
(
const word& name,
const objectRegistry&,
const dictionary&,
const bool loadFromFiles = false
);
//- Destructor
virtual ~moleFractions();
// Member Functions
//- Return name of the moleFractions functionObject
virtual const word& name() const
{
return name_;
}
//- Read the moleFractions data
virtual void read(const dictionary&);
//- Execute, currently does nothing
virtual void execute();
//- Execute at the final time-loop, currently does nothing
virtual void end();
//- Called when time was set at the end of the Time::operator++
virtual void timeSet();
//- Calculate the moleFractions and write
virtual void write();
//- Update for changes of mesh
virtual void updateMesh(const mapPolyMesh&)
{}
//- Update for changes of mesh
virtual void movePoints(const polyMesh&)
{}
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
# include "moleFractions.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 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 "moleFractionsFunctionObjects.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTemplateTypeNameAndDebugWithName
(
moleFractions<psiReactionThermo>,
"psiReactionThermoMoleFractions",
0
);
defineTemplateTypeNameAndDebugWithName
(
psiReactionThermoMoleFractionsFunctionObject,
"psiReactionThermoMoleFractions",
0
);
addToRunTimeSelectionTable
(
functionObject,
psiReactionThermoMoleFractionsFunctionObject,
dictionary
);
defineTemplateTypeNameAndDebugWithName
(
moleFractions<rhoReactionThermo>,
"rhoReactionThermoMoleFractions",
0
);
defineTemplateTypeNameAndDebugWithName
(
rhoReactionThermoMoleFractionsFunctionObject,
"rhoReactionThermoMoleFractions",
0
);
addToRunTimeSelectionTable
(
functionObject,
rhoReactionThermoMoleFractionsFunctionObject,
dictionary
);
}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 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/>.
Typedef
Foam::moleFractionsFunctionObject
Description
FunctionObject wrapper around moleFractions to allow
it to be created via the functions entry within controlDict.
SourceFiles
moleFractionsFunctionObjects.C
\*---------------------------------------------------------------------------*/
#ifndef moleFractionsFunctionObjects_H
#define moleFractionsFunctionObjects_H
#include "OutputFilterFunctionObject.H"
#include "moleFractions.H"
#include "psiReactionThermo.H"
#include "rhoReactionThermo.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
typedef OutputFilterFunctionObject<moleFractions<psiReactionThermo>>
psiReactionThermoMoleFractionsFunctionObject;
typedef OutputFilterFunctionObject<moleFractions<rhoReactionThermo>>
rhoReactionThermoMoleFractionsFunctionObject;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment