mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Added InflowBoundaryModel submodel with NoInflow and FreeStream models. Free stream inserts particles at a specified patch. Currently inserting particles half way along the face centre to cell centre line - will randomise. Fixed label/scalar inconsistency between declaration and accessors of nParticles.
This commit is contained in:
@ -11,5 +11,6 @@ clouds/derived/dsmcCloud/dsmcCloud.C
|
||||
parcels/derived/dsmcParcel/defineDsmcParcel.C
|
||||
parcels/derived/dsmcParcel/makeDsmcParcelBinaryCollisionModels.C
|
||||
parcels/derived/dsmcParcel/makeDsmcParcelWallInteractionModels.C
|
||||
parcels/derived/dsmcParcel/makeDsmcParcelInflowBoundaryModels.C
|
||||
|
||||
LIB = $(FOAM_LIBBIN)/libdsmc
|
||||
|
||||
@ -27,6 +27,7 @@ License
|
||||
#include "DsmcCloud.H"
|
||||
#include "BinaryCollisionModel.H"
|
||||
#include "WallInteractionModel.H"
|
||||
#include "InflowBoundaryModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
@ -234,8 +235,8 @@ void Foam::DsmcCloud<ParcelType>::collisions()
|
||||
scalar sigmaTcRMax = sigmaTcRMax_[celli];
|
||||
|
||||
scalar selectedPairs = collisionSelectionRemainder_[celli]
|
||||
+ 0.5*nC*(nC-1)*nParticle_*sigmaTcRMax*deltaT
|
||||
/mesh_.cellVolumes()[celli];
|
||||
+ 0.5*nC*(nC-1)*nParticle_*sigmaTcRMax*deltaT
|
||||
/mesh_.cellVolumes()[celli];
|
||||
|
||||
label nCandidates(selectedPairs);
|
||||
|
||||
@ -403,6 +404,14 @@ Foam::DsmcCloud<ParcelType>::DsmcCloud
|
||||
particleProperties_,
|
||||
*this
|
||||
)
|
||||
),
|
||||
inflowBoundaryModel_
|
||||
(
|
||||
InflowBoundaryModel<DsmcCloud<ParcelType> >::New
|
||||
(
|
||||
particleProperties_,
|
||||
*this
|
||||
)
|
||||
)
|
||||
{
|
||||
buildConstProps();
|
||||
@ -490,7 +499,8 @@ Foam::DsmcCloud<ParcelType>::DsmcCloud
|
||||
)
|
||||
),
|
||||
binaryCollisionModel_(),
|
||||
wallInteractionModel_()
|
||||
wallInteractionModel_(),
|
||||
inflowBoundaryModel_()
|
||||
{
|
||||
clear();
|
||||
|
||||
@ -526,13 +536,14 @@ void Foam::DsmcCloud<ParcelType>::evolve()
|
||||
{
|
||||
typename ParcelType::trackData td(*this);
|
||||
|
||||
//this->injection().inject(td);
|
||||
|
||||
if (debug)
|
||||
{
|
||||
this->dumpParticlePositions();
|
||||
}
|
||||
|
||||
// Insert new particles from the inflow boundary
|
||||
this->inflowBoundary().inflow();
|
||||
|
||||
// Move the particles ballistically with their current velocities
|
||||
Cloud<ParcelType>::move(td);
|
||||
|
||||
|
||||
@ -59,6 +59,9 @@ class BinaryCollisionModel;
|
||||
template<class CloudType>
|
||||
class WallInteractionModel;
|
||||
|
||||
template<class CloudType>
|
||||
class InflowBoundaryModel;
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class DsmcCloud Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -127,6 +130,10 @@ class DsmcCloud
|
||||
autoPtr<WallInteractionModel<DsmcCloud<ParcelType> > >
|
||||
wallInteractionModel_;
|
||||
|
||||
//- Inflow boundary model
|
||||
autoPtr<InflowBoundaryModel<DsmcCloud<ParcelType> > >
|
||||
inflowBoundaryModel_;
|
||||
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
@ -205,7 +212,7 @@ public:
|
||||
|
||||
//- Return the number of real particles represented by one
|
||||
// parcel
|
||||
inline label nParticle() const;
|
||||
inline scalar nParticle() const;
|
||||
|
||||
//- Return the cell occupancy addressing
|
||||
inline const List<DynamicList<ParcelType*> >&
|
||||
@ -292,6 +299,14 @@ public:
|
||||
inline WallInteractionModel<DsmcCloud<ParcelType> >&
|
||||
wallInteraction();
|
||||
|
||||
//- Return reference to wall interaction model
|
||||
inline const InflowBoundaryModel<DsmcCloud<ParcelType> >&
|
||||
inflowBoundary() const;
|
||||
|
||||
//- Return non-const reference to wall interaction model
|
||||
inline InflowBoundaryModel<DsmcCloud<ParcelType> >&
|
||||
inflowBoundary();
|
||||
|
||||
|
||||
// Check
|
||||
|
||||
|
||||
@ -57,7 +57,7 @@ Foam::DsmcCloud<ParcelType>::typeIdList() const
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline Foam::label Foam::DsmcCloud<ParcelType>::nParticle() const
|
||||
inline Foam::scalar Foam::DsmcCloud<ParcelType>::nParticle() const
|
||||
{
|
||||
return nParticle_;
|
||||
}
|
||||
@ -166,6 +166,22 @@ Foam::DsmcCloud<ParcelType>::wallInteraction()
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline const Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
|
||||
Foam::DsmcCloud<ParcelType>::inflowBoundary() const
|
||||
{
|
||||
return inflowBoundaryModel_;
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline Foam::InflowBoundaryModel<Foam::DsmcCloud<ParcelType> >&
|
||||
Foam::DsmcCloud<ParcelType>::inflowBoundary()
|
||||
{
|
||||
return inflowBoundaryModel_();
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline Foam::scalar Foam::DsmcCloud<ParcelType>::massInSystem() const
|
||||
{
|
||||
|
||||
@ -0,0 +1,52 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "dsmcParcel.H"
|
||||
#include "DsmcCloud.H"
|
||||
#include "FreeStream.H"
|
||||
#include "NoInflow.H"
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
makeInflowBoundaryModel(DsmcCloud<dsmcParcel>);
|
||||
|
||||
// Add instances of inflow boundary model to the table
|
||||
makeInflowBoundaryModelType
|
||||
(
|
||||
FreeStream,
|
||||
DsmcCloud,
|
||||
dsmcParcel
|
||||
);
|
||||
makeInflowBoundaryModelType
|
||||
(
|
||||
NoInflow,
|
||||
DsmcCloud,
|
||||
dsmcParcel
|
||||
);
|
||||
};
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -109,7 +109,8 @@ Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::LarsenBorgnakkeVariableHardS
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
template <class CloudType>
|
||||
Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::~LarsenBorgnakkeVariableHardSphere()
|
||||
Foam::LarsenBorgnakkeVariableHardSphere<CloudType>::
|
||||
~LarsenBorgnakkeVariableHardSphere()
|
||||
{}
|
||||
|
||||
|
||||
|
||||
@ -50,7 +50,7 @@ class LarsenBorgnakkeVariableHardSphere
|
||||
{
|
||||
// Private data
|
||||
|
||||
//- Temperature
|
||||
//- Relaxation collision number
|
||||
const scalar relaxationCollisionNumber_;
|
||||
|
||||
|
||||
|
||||
@ -0,0 +1,218 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "FreeStream.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
template <class CloudType>
|
||||
Foam::FreeStream<CloudType>::FreeStream
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& cloud
|
||||
)
|
||||
:
|
||||
InflowBoundaryModel<CloudType>(dict, cloud, typeName),
|
||||
patchIndex_(),
|
||||
temperature_(readScalar(this->coeffDict().lookup("temperature"))),
|
||||
velocity_(this->coeffDict().lookup("velocity")),
|
||||
moleculeTypeIds_(),
|
||||
numberDensities_(),
|
||||
particleFluxAccumulators_()
|
||||
{
|
||||
word patchName = this->coeffDict().lookup("patch");
|
||||
|
||||
patchIndex_ = cloud.mesh().boundaryMesh().findPatchID(patchName);
|
||||
|
||||
const polyPatch& patch = cloud.mesh().boundaryMesh()[patchIndex_];
|
||||
|
||||
if (patchIndex_ == -1)
|
||||
{
|
||||
FatalErrorIn("Foam::DsmcCloud<ParcelType>::initialise")
|
||||
<< "patch " << patchName << " not found." << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
const dictionary& numberDensitiesDict
|
||||
(
|
||||
this->coeffDict().subDict("numberDensities")
|
||||
);
|
||||
|
||||
List<word> molecules(numberDensitiesDict.toc());
|
||||
|
||||
numberDensities_.setSize(molecules.size());
|
||||
|
||||
moleculeTypeIds_.setSize(molecules.size());
|
||||
|
||||
forAll(molecules, i)
|
||||
{
|
||||
numberDensities_[i] = readScalar
|
||||
(
|
||||
numberDensitiesDict.lookup(molecules[i])
|
||||
);
|
||||
|
||||
moleculeTypeIds_[i] = findIndex(cloud.typeIdList(), molecules[i]);
|
||||
|
||||
if (moleculeTypeIds_[i] == -1)
|
||||
{
|
||||
FatalErrorIn("Foam::DsmcCloud<ParcelType>::initialise")
|
||||
<< "typeId " << molecules[i] << "not defined in cloud." << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
numberDensities_ /= cloud.nParticle();
|
||||
|
||||
particleFluxAccumulators_.setSize
|
||||
(
|
||||
molecules.size(),
|
||||
Field<scalar>(patch.size(), 0)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
template <class CloudType>
|
||||
Foam::FreeStream<CloudType>::~FreeStream()
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template <class CloudType>
|
||||
void Foam::FreeStream<CloudType>::inflow()
|
||||
{
|
||||
CloudType& cloud(this->owner());
|
||||
|
||||
const polyMesh& mesh(cloud.mesh());
|
||||
|
||||
const scalar deltaT = mesh.time().deltaT().value();
|
||||
|
||||
Random& rndGen(cloud.rndGen());
|
||||
|
||||
const polyPatch& patch = mesh.boundaryMesh()[patchIndex_];
|
||||
|
||||
label particlesInserted = 0;
|
||||
|
||||
// Add mass to the accumulators. negative face area dotted with the
|
||||
// velocity to point flux into the domain.
|
||||
|
||||
forAll(particleFluxAccumulators_, i)
|
||||
{
|
||||
particleFluxAccumulators_[i] +=
|
||||
-patch.faceAreas() & (velocity_*numberDensities_[i]*deltaT);
|
||||
|
||||
forAll(particleFluxAccumulators_[i], f)
|
||||
{
|
||||
scalar& faceAccumulator = particleFluxAccumulators_[i][f];
|
||||
|
||||
// Number of particles to insert
|
||||
label nI = max(label(faceAccumulator), 0);
|
||||
|
||||
faceAccumulator -= nI;
|
||||
|
||||
label typeId = moleculeTypeIds_[i];
|
||||
|
||||
scalar mass = cloud.constProps(typeId).mass();
|
||||
|
||||
labelList faceVertices = patch[f];
|
||||
|
||||
label globalFaceIndex = f + patch.start();
|
||||
|
||||
label cell = mesh.faceOwner()[globalFaceIndex];
|
||||
|
||||
const vector& fC = patch.faceCentres()[f];
|
||||
|
||||
for (label n = 0; n < nI; n++)
|
||||
{
|
||||
// Temporarily insert particles half way between the face and
|
||||
// cell centres
|
||||
vector p = 0.5*(fC + mesh.cellCentres()[cell]);
|
||||
|
||||
// Normal unit vector *negative* so normal is pointing into the
|
||||
// domain
|
||||
vector nw = patch.faceAreas()[f];
|
||||
nw /= -mag(nw);
|
||||
|
||||
// Wall tangential unit vector. Use the direction between the
|
||||
// face centre and the first vertex in the list
|
||||
vector tw1 = fC - (mesh.points()[faceVertices[0]]);
|
||||
tw1 /= mag(tw1);
|
||||
|
||||
// Other tangential unit vector. Rescaling in case face is not
|
||||
// flat and nw and tw1 aren't perfectly orthogonal
|
||||
vector tw2 = nw ^ tw1;
|
||||
tw2 /= mag(tw2);
|
||||
|
||||
scalar C = sqrt(CloudType::kb*temperature_/mass);
|
||||
|
||||
vector U =
|
||||
C
|
||||
*(
|
||||
rndGen.GaussNormal()*tw1
|
||||
+ rndGen.GaussNormal()*tw2
|
||||
- sqrt(-2.0*log(max(1 - rndGen.scalar01(),VSMALL)))*nw
|
||||
);
|
||||
|
||||
U += velocity_;
|
||||
|
||||
scalar Ei =
|
||||
0.5*cloud.constProps(typeId).internalDegreesOfFreedom()
|
||||
*CloudType::kb*temperature_;
|
||||
|
||||
cloud.addNewParcel
|
||||
(
|
||||
p,
|
||||
U,
|
||||
Ei,
|
||||
cell,
|
||||
typeId
|
||||
);
|
||||
|
||||
particlesInserted++;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
reduce(particlesInserted, sumOp<label>());
|
||||
|
||||
Info<< " Particles inserted = "
|
||||
<< particlesInserted << endl;
|
||||
|
||||
// Info<< "insert particles now! " << nl
|
||||
// << temperature_ << nl
|
||||
// << velocity_ << nl
|
||||
// << moleculeTypeIds_ << nl
|
||||
// << numberDensities_ << nl
|
||||
// << particleFluxAccumulators_ << nl
|
||||
// << patch
|
||||
// << endl;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,117 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
Class
|
||||
Foam::FreeStream
|
||||
|
||||
Description
|
||||
Inserting new particles across the faces of a specified patch for a free
|
||||
stream. Uniform values of temperature, velocity and number densities
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef FreeStream_H
|
||||
#define FreeStream_H
|
||||
|
||||
#include "InflowBoundaryModel.H"
|
||||
#include "polyMesh.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class FreeStream Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
template<class CloudType>
|
||||
class FreeStream
|
||||
:
|
||||
public InflowBoundaryModel<CloudType>
|
||||
{
|
||||
|
||||
// Private data
|
||||
|
||||
//- Index of patch to introduce particles across
|
||||
label patchIndex_;
|
||||
|
||||
//- Temperature of the free stream
|
||||
scalar temperature_;
|
||||
|
||||
//- Velocity of the free stream
|
||||
vector velocity_;
|
||||
|
||||
//- The molecule types to be introduced
|
||||
List<label> moleculeTypeIds_;
|
||||
|
||||
//- The number density of the species in the inflow
|
||||
Field<scalar> numberDensities_;
|
||||
|
||||
//- A List of Fields, one Field for every species to be introduced, each
|
||||
// field entry corresponding to a face on the patch to be injected
|
||||
// across.
|
||||
List<Field<scalar> > particleFluxAccumulators_;
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("FreeStream");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from dictionary
|
||||
FreeStream
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& cloud
|
||||
);
|
||||
|
||||
|
||||
// Destructor
|
||||
virtual ~FreeStream();
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Introduce particles
|
||||
virtual void inflow();
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#ifdef NoRepository
|
||||
# include "FreeStream.C"
|
||||
#endif
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,86 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "InflowBoundaryModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
template<class CloudType>
|
||||
Foam::InflowBoundaryModel<CloudType>::InflowBoundaryModel
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& owner,
|
||||
const word& type
|
||||
)
|
||||
: dict_(dict),
|
||||
owner_(owner),
|
||||
coeffDict_(dict.subDict(type + "Coeffs"))
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
template<class CloudType>
|
||||
Foam::InflowBoundaryModel<CloudType>::~InflowBoundaryModel()
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class CloudType>
|
||||
const CloudType& Foam::InflowBoundaryModel<CloudType>::owner() const
|
||||
{
|
||||
return owner_;
|
||||
}
|
||||
|
||||
|
||||
template<class CloudType>
|
||||
CloudType& Foam::InflowBoundaryModel<CloudType>::owner()
|
||||
{
|
||||
return owner_;
|
||||
}
|
||||
|
||||
|
||||
template<class CloudType>
|
||||
const Foam::dictionary& Foam::InflowBoundaryModel<CloudType>::dict() const
|
||||
{
|
||||
return dict_;
|
||||
}
|
||||
|
||||
|
||||
template<class CloudType>
|
||||
const Foam::dictionary& Foam::InflowBoundaryModel<CloudType>::coeffDict() const
|
||||
{
|
||||
return coeffDict_;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#include "NewInflowBoundaryModel.C"
|
||||
|
||||
// ************************************************************************* //
|
||||
|
||||
@ -0,0 +1,167 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
Class
|
||||
Foam::InflowBoundaryModel
|
||||
|
||||
|
||||
Description
|
||||
Templated inflow boundary model class
|
||||
|
||||
SourceFiles
|
||||
InflowBoundaryModel.C
|
||||
NewInflowBoundaryModel.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef InflowBoundaryModel_H
|
||||
#define InflowBoundaryModel_H
|
||||
|
||||
#include "IOdictionary.H"
|
||||
#include "autoPtr.H"
|
||||
#include "runTimeSelectionTables.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class InflowBoundaryModel Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
template<class CloudType>
|
||||
class InflowBoundaryModel
|
||||
{
|
||||
// Private data
|
||||
|
||||
//- The cloud dictionary
|
||||
const dictionary& dict_;
|
||||
|
||||
// Reference to the owner cloud class
|
||||
CloudType& owner_;
|
||||
|
||||
//- The coefficients dictionary
|
||||
const dictionary coeffDict_;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("InflowBoundaryModel");
|
||||
|
||||
//- Declare runtime constructor selection table
|
||||
declareRunTimeSelectionTable
|
||||
(
|
||||
autoPtr,
|
||||
InflowBoundaryModel,
|
||||
dictionary,
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& owner
|
||||
),
|
||||
(dict, owner)
|
||||
);
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from dictionary
|
||||
InflowBoundaryModel
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& owner,
|
||||
const word& type
|
||||
);
|
||||
|
||||
|
||||
//- Destructor
|
||||
virtual ~InflowBoundaryModel();
|
||||
|
||||
|
||||
//- Selector
|
||||
static autoPtr<InflowBoundaryModel<CloudType> > New
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& owner
|
||||
);
|
||||
|
||||
|
||||
// Access
|
||||
|
||||
//- Return const access the owner cloud object
|
||||
inline const CloudType& owner() const;
|
||||
|
||||
//- Return non-const access the owner cloud object for manipulation
|
||||
inline CloudType& owner();
|
||||
|
||||
//- Return the owner cloud dictionary
|
||||
inline const dictionary& dict() const;
|
||||
|
||||
//- Return the coefficients dictionary
|
||||
inline const dictionary& coeffDict() const;
|
||||
|
||||
//- Introduce particles
|
||||
virtual void inflow() = 0;
|
||||
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#define makeInflowBoundaryModel(CloudType) \
|
||||
\
|
||||
defineNamedTemplateTypeNameAndDebug(InflowBoundaryModel<CloudType>, 0); \
|
||||
\
|
||||
defineTemplateRunTimeSelectionTable \
|
||||
( \
|
||||
InflowBoundaryModel<CloudType>, \
|
||||
dictionary \
|
||||
);
|
||||
|
||||
|
||||
#define makeInflowBoundaryModelType(SS, CloudType, ParcelType) \
|
||||
\
|
||||
defineNamedTemplateTypeNameAndDebug(SS<CloudType<ParcelType> >, 0); \
|
||||
\
|
||||
InflowBoundaryModel<CloudType<ParcelType> >:: \
|
||||
adddictionaryConstructorToTable<SS<CloudType<ParcelType> > > \
|
||||
add##SS##CloudType##ParcelType##ConstructorToTable_;
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#ifdef NoRepository
|
||||
# include "InflowBoundaryModel.C"
|
||||
#endif
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,66 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "InflowBoundaryModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
template<class CloudType>
|
||||
Foam::autoPtr<Foam::InflowBoundaryModel<CloudType> >
|
||||
Foam::InflowBoundaryModel<CloudType>::New
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& owner
|
||||
)
|
||||
{
|
||||
word InflowBoundaryModelType
|
||||
(
|
||||
dict.lookup("InflowBoundaryModel")
|
||||
);
|
||||
|
||||
Info<< "Selecting InflowBoundaryModel " << InflowBoundaryModelType << endl;
|
||||
|
||||
typename dictionaryConstructorTable::iterator cstrIter =
|
||||
dictionaryConstructorTablePtr_->find(InflowBoundaryModelType);
|
||||
|
||||
if (cstrIter == dictionaryConstructorTablePtr_->end())
|
||||
{
|
||||
FatalErrorIn
|
||||
(
|
||||
"InflowBoundaryModel<CloudType>::New"
|
||||
"(const dictionary&, CloudType&)"
|
||||
) << "Unknown InflowBoundaryModelType type "
|
||||
<< InflowBoundaryModelType
|
||||
<< ", constructor not in hash table" << nl << nl
|
||||
<< " Valid InflowBoundaryModel types are :" << nl
|
||||
<< dictionaryConstructorTablePtr_->toc() << exit(FatalError);
|
||||
}
|
||||
|
||||
return autoPtr<InflowBoundaryModel<CloudType> >(cstrIter()(dict, owner));
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,56 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "NoInflow.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
template <class CloudType>
|
||||
Foam::NoInflow<CloudType>::NoInflow
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& cloud
|
||||
)
|
||||
:
|
||||
InflowBoundaryModel<CloudType>(dict, cloud, typeName)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
template <class CloudType>
|
||||
Foam::NoInflow<CloudType>::~NoInflow()
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template <class CloudType>
|
||||
void Foam::NoInflow<CloudType>::inflow()
|
||||
{}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,92 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
Class
|
||||
Foam::NoInflow
|
||||
|
||||
Description
|
||||
Not inserting any particles
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef NoInflow_H
|
||||
#define NoInflow_H
|
||||
|
||||
#include "InflowBoundaryModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class NoInflow Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
template<class CloudType>
|
||||
class NoInflow
|
||||
:
|
||||
public InflowBoundaryModel<CloudType>
|
||||
{
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("NoInflow");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from dictionary
|
||||
NoInflow
|
||||
(
|
||||
const dictionary& dict,
|
||||
CloudType& cloud
|
||||
);
|
||||
|
||||
|
||||
// Destructor
|
||||
virtual ~NoInflow();
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Introduce particles (none in this case)
|
||||
virtual void inflow();
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#ifdef NoRepository
|
||||
# include "NoInflow.C"
|
||||
#endif
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -48,13 +48,6 @@ Foam::MaxwellianThermal<CloudType>::~MaxwellianThermal()
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class CloudType>
|
||||
bool Foam::MaxwellianThermal<CloudType>::active() const
|
||||
{
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
template <class CloudType>
|
||||
void Foam::MaxwellianThermal<CloudType>::correct
|
||||
(
|
||||
|
||||
@ -72,9 +72,6 @@ public:
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Flag to indicate whether model activates heat transfer model
|
||||
bool active() const;
|
||||
|
||||
//- Apply wall correction
|
||||
virtual void correct
|
||||
(
|
||||
|
||||
@ -50,13 +50,6 @@ Foam::SpecularReflection<CloudType>::~SpecularReflection()
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class CloudType>
|
||||
bool Foam::SpecularReflection<CloudType>::active() const
|
||||
{
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
template <class CloudType>
|
||||
void Foam::SpecularReflection<CloudType>::correct
|
||||
(
|
||||
|
||||
@ -70,9 +70,6 @@ public:
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Flag to indicate whether model activates heat transfer model
|
||||
bool active() const;
|
||||
|
||||
//- Apply wall correction
|
||||
virtual void correct
|
||||
(
|
||||
|
||||
@ -124,9 +124,6 @@ public:
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Flag to indicate whether model activates heat transfer model
|
||||
virtual bool active() const = 0;
|
||||
|
||||
//- Apply wall correction
|
||||
virtual void correct
|
||||
(
|
||||
|
||||
Reference in New Issue
Block a user