mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: uniformInterpolatedDisplacement: point boundary condition by interpolating existing motion
This commit is contained in:
@ -38,6 +38,7 @@ $(derivedPoint)/surfaceDisplacement/surfaceDisplacementPointPatchVectorField.C
|
|||||||
$(derivedPoint)/waveDisplacement/waveDisplacementPointPatchVectorField.C
|
$(derivedPoint)/waveDisplacement/waveDisplacementPointPatchVectorField.C
|
||||||
|
|
||||||
$(derivedPoint)/timeVaryingMappedFixedValue/timeVaryingMappedFixedValuePointPatchFields.C
|
$(derivedPoint)/timeVaryingMappedFixedValue/timeVaryingMappedFixedValuePointPatchFields.C
|
||||||
|
$(derivedPoint)/uniformInterpolatedDisplacement/uniformInterpolatedDisplacementPointPatchVectorField.C
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|||||||
@ -0,0 +1,291 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | Copyright (C) 2012 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 "uniformInterpolatedDisplacementPointPatchVectorField.H"
|
||||||
|
#include "pointPatchFields.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
#include "Time.H"
|
||||||
|
#include "polyMesh.H"
|
||||||
|
#include "interpolationWeights.H"
|
||||||
|
#include "uniformInterpolate.H"
|
||||||
|
#include "ReadFields.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField::
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const pointPatch& p,
|
||||||
|
const DimensionedField<vector, pointMesh>& iF
|
||||||
|
)
|
||||||
|
:
|
||||||
|
fixedValuePointPatchField<vector>(p, iF)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField::
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const pointPatch& p,
|
||||||
|
const DimensionedField<vector, pointMesh>& iF,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
:
|
||||||
|
fixedValuePointPatchField<vector>(p, iF, dict),
|
||||||
|
fieldName_(dict.lookup("fieldName")),
|
||||||
|
interpolationScheme_(dict.lookup("interpolationScheme"))
|
||||||
|
{
|
||||||
|
const pointMesh& pMesh = this->dimensionedInternalField().mesh();
|
||||||
|
|
||||||
|
// Read time values
|
||||||
|
instantList allTimes = Time::findTimes(pMesh().time().path());
|
||||||
|
|
||||||
|
// Only keep those that contain the field
|
||||||
|
DynamicList<word> names(allTimes.size());
|
||||||
|
DynamicList<scalar> values(allTimes.size());
|
||||||
|
|
||||||
|
forAll(allTimes, i)
|
||||||
|
{
|
||||||
|
IOobject io
|
||||||
|
(
|
||||||
|
fieldName_,
|
||||||
|
allTimes[i].name(),
|
||||||
|
pMesh(),
|
||||||
|
IOobject::MUST_READ,
|
||||||
|
IOobject::NO_WRITE,
|
||||||
|
false
|
||||||
|
);
|
||||||
|
if (io.headerOk())
|
||||||
|
{
|
||||||
|
names.append(allTimes[i].name());
|
||||||
|
values.append(allTimes[i].value());
|
||||||
|
}
|
||||||
|
}
|
||||||
|
timeNames_.transfer(names);
|
||||||
|
timeVals_.transfer(values);
|
||||||
|
|
||||||
|
Info<< type() << " : found " << fieldName_ << " for times "
|
||||||
|
<< timeNames_ << endl;
|
||||||
|
|
||||||
|
if (timeNames_.size() < 1)
|
||||||
|
{
|
||||||
|
FatalErrorIn
|
||||||
|
(
|
||||||
|
"uniformInterpolatedDisplacementPointPatchVectorField::\n"
|
||||||
|
"uniformInterpolatedDisplacementPointPatchVectorField\n"
|
||||||
|
"(\n"
|
||||||
|
" const pointPatch&,\n"
|
||||||
|
" const DimensionedField<vector, pointMesh>&,\n"
|
||||||
|
" const dictionary&\n"
|
||||||
|
")\n"
|
||||||
|
) << "Did not find any times with " << fieldName_
|
||||||
|
<< exit(FatalError);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
if (!dict.found("value"))
|
||||||
|
{
|
||||||
|
updateCoeffs();
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField::
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const uniformInterpolatedDisplacementPointPatchVectorField& ptf,
|
||||||
|
const pointPatch& p,
|
||||||
|
const DimensionedField<vector, pointMesh>& iF,
|
||||||
|
const pointPatchFieldMapper& mapper
|
||||||
|
)
|
||||||
|
:
|
||||||
|
fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
|
||||||
|
fieldName_(ptf.fieldName_),
|
||||||
|
interpolationScheme_(ptf.interpolationScheme_),
|
||||||
|
timeNames_(ptf.timeNames_),
|
||||||
|
timeVals_(ptf.timeVals_),
|
||||||
|
interpolatorPtr_(ptf.interpolatorPtr_)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField::
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const uniformInterpolatedDisplacementPointPatchVectorField& ptf,
|
||||||
|
const DimensionedField<vector, pointMesh>& iF
|
||||||
|
)
|
||||||
|
:
|
||||||
|
fixedValuePointPatchField<vector>(ptf, iF),
|
||||||
|
fieldName_(ptf.fieldName_),
|
||||||
|
interpolationScheme_(ptf.interpolationScheme_),
|
||||||
|
timeNames_(ptf.timeNames_),
|
||||||
|
timeVals_(ptf.timeVals_),
|
||||||
|
interpolatorPtr_(ptf.interpolatorPtr_)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
void uniformInterpolatedDisplacementPointPatchVectorField::updateCoeffs()
|
||||||
|
{
|
||||||
|
if (this->updated())
|
||||||
|
{
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (!interpolatorPtr_.valid())
|
||||||
|
{
|
||||||
|
interpolatorPtr_ = interpolationWeights::New
|
||||||
|
(
|
||||||
|
interpolationScheme_,
|
||||||
|
timeVals_
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
const pointMesh& pMesh = this->dimensionedInternalField().mesh();
|
||||||
|
const Time& t = pMesh().time();
|
||||||
|
|
||||||
|
// Update indices of times and weights
|
||||||
|
bool timesChanged = interpolatorPtr_->valueWeights
|
||||||
|
(
|
||||||
|
t.timeOutputValue(),
|
||||||
|
currentIndices_,
|
||||||
|
currentWeights_
|
||||||
|
);
|
||||||
|
|
||||||
|
const wordList currentTimeNames
|
||||||
|
(
|
||||||
|
UIndirectList<word>(timeNames_, currentIndices_)
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
// Load if necessary fields for this interpolation
|
||||||
|
if (timesChanged)
|
||||||
|
{
|
||||||
|
objectRegistry& fieldsCache = const_cast<objectRegistry&>
|
||||||
|
(
|
||||||
|
pMesh.thisDb().subRegistry("fieldsCache", true)
|
||||||
|
);
|
||||||
|
// Save old times so we now which ones have been loaded and need
|
||||||
|
// 'correctBoundaryConditions'. Bit messy.
|
||||||
|
HashSet<word> oldTimes(fieldsCache.toc());
|
||||||
|
|
||||||
|
ReadFields<pointVectorField>
|
||||||
|
(
|
||||||
|
fieldName_,
|
||||||
|
pMesh,
|
||||||
|
currentTimeNames
|
||||||
|
);
|
||||||
|
|
||||||
|
forAllConstIter(objectRegistry, fieldsCache, fieldsCacheIter)
|
||||||
|
{
|
||||||
|
if (!oldTimes.found(fieldsCacheIter.key()))
|
||||||
|
{
|
||||||
|
// Newly loaded fields. Make sure the internal
|
||||||
|
// values are consistent with the boundary conditions.
|
||||||
|
// This is quite often not the case since these
|
||||||
|
// fields typically are constructed 'by hand'
|
||||||
|
|
||||||
|
const objectRegistry& timeCache = dynamic_cast
|
||||||
|
<
|
||||||
|
const objectRegistry&
|
||||||
|
>(*fieldsCacheIter());
|
||||||
|
|
||||||
|
|
||||||
|
pointVectorField& d = const_cast<pointVectorField&>
|
||||||
|
(
|
||||||
|
timeCache.lookupObject<pointVectorField>
|
||||||
|
(
|
||||||
|
fieldName_
|
||||||
|
)
|
||||||
|
);
|
||||||
|
d.correctBoundaryConditions();
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// Interpolate the whole field
|
||||||
|
pointVectorField result
|
||||||
|
(
|
||||||
|
uniformInterpolate<pointVectorField>
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
word("uniformInterpolate(")
|
||||||
|
+ this->dimensionedInternalField().name()
|
||||||
|
+ ')',
|
||||||
|
pMesh.time().timeName(),
|
||||||
|
pMesh.thisDb(),
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::AUTO_WRITE
|
||||||
|
),
|
||||||
|
fieldName_,
|
||||||
|
currentTimeNames,
|
||||||
|
currentWeights_
|
||||||
|
)
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
// Extract back from the internal field
|
||||||
|
this->operator==
|
||||||
|
(
|
||||||
|
this->patchInternalField(result.dimensionedInternalField())
|
||||||
|
);
|
||||||
|
|
||||||
|
fixedValuePointPatchField<vector>::updateCoeffs();
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void uniformInterpolatedDisplacementPointPatchVectorField::write(Ostream& os)
|
||||||
|
const
|
||||||
|
{
|
||||||
|
pointPatchField<vector>::write(os);
|
||||||
|
os.writeKeyword("fieldName")
|
||||||
|
<< fieldName_ << token::END_STATEMENT << nl;
|
||||||
|
os.writeKeyword("interpolationScheme")
|
||||||
|
<< interpolationScheme_ << token::END_STATEMENT << nl;
|
||||||
|
writeEntry("value", os);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
makePointPatchTypeField
|
||||||
|
(
|
||||||
|
pointPatchVectorField,
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
);
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,183 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | Copyright (C) 2012 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::uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
|
||||||
|
Description
|
||||||
|
Interpolates pre-specified motion.
|
||||||
|
|
||||||
|
Motion specified as pointVectorFields. E.g.
|
||||||
|
|
||||||
|
walls
|
||||||
|
{
|
||||||
|
type uniformInterpolatedDisplacement;
|
||||||
|
value uniform (0 0 0);
|
||||||
|
fieldName wantedDisplacement;
|
||||||
|
interpolationScheme linear;
|
||||||
|
}
|
||||||
|
|
||||||
|
This will scan the case for 'wantedDisplacement' pointVectorFields
|
||||||
|
and interpolate those in time (using 'linear' interpolation) to
|
||||||
|
obtain the current displacement.
|
||||||
|
The advantage of specifying displacement in this way is that it
|
||||||
|
automatically works through decomposePar.
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef uniformInterpolatedDisplacementPointPatchVectorField_H
|
||||||
|
#define uniformInterpolatedDisplacementPointPatchVectorField_H
|
||||||
|
|
||||||
|
#include "fixedValuePointPatchField.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
class interpolationWeights;
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class uniformInterpolatedDisplacementPointPatchVectorField Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
:
|
||||||
|
public fixedValuePointPatchField<vector>
|
||||||
|
{
|
||||||
|
// Private data
|
||||||
|
|
||||||
|
//- Name of displacement field
|
||||||
|
const word fieldName_;
|
||||||
|
|
||||||
|
const word interpolationScheme_;
|
||||||
|
|
||||||
|
//- Times with pre-specified displacement
|
||||||
|
wordList timeNames_;
|
||||||
|
|
||||||
|
//- Times with pre-specified displacement
|
||||||
|
scalarField timeVals_;
|
||||||
|
|
||||||
|
//- User-specified interpolator
|
||||||
|
autoPtr<interpolationWeights> interpolatorPtr_;
|
||||||
|
|
||||||
|
|
||||||
|
//- Cached interpolation times
|
||||||
|
labelList currentIndices_;
|
||||||
|
|
||||||
|
//- Cached interpolation weights
|
||||||
|
scalarField currentWeights_;
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("uniformInterpolatedDisplacement");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from patch and internal field
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const pointPatch&,
|
||||||
|
const DimensionedField<vector, pointMesh>&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct from patch, internal field and dictionary
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const pointPatch&,
|
||||||
|
const DimensionedField<vector, pointMesh>&,
|
||||||
|
const dictionary&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct by mapping given patchField<vector> onto a new patch
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const uniformInterpolatedDisplacementPointPatchVectorField&,
|
||||||
|
const pointPatch&,
|
||||||
|
const DimensionedField<vector, pointMesh>&,
|
||||||
|
const pointPatchFieldMapper&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct and return a clone
|
||||||
|
virtual autoPtr<pointPatchField<vector> > clone() const
|
||||||
|
{
|
||||||
|
return autoPtr<pointPatchField<vector> >
|
||||||
|
(
|
||||||
|
new uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
*this
|
||||||
|
)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
//- Construct as copy setting internal field reference
|
||||||
|
uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
const uniformInterpolatedDisplacementPointPatchVectorField&,
|
||||||
|
const DimensionedField<vector, pointMesh>&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct and return a clone setting internal field reference
|
||||||
|
virtual autoPtr<pointPatchField<vector> > clone
|
||||||
|
(
|
||||||
|
const DimensionedField<vector, pointMesh>& iF
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
return autoPtr<pointPatchField<vector> >
|
||||||
|
(
|
||||||
|
new uniformInterpolatedDisplacementPointPatchVectorField
|
||||||
|
(
|
||||||
|
*this,
|
||||||
|
iF
|
||||||
|
)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// Member functions
|
||||||
|
|
||||||
|
// Evaluation functions
|
||||||
|
|
||||||
|
//- Update the coefficients associated with the patch field
|
||||||
|
virtual void updateCoeffs();
|
||||||
|
|
||||||
|
|
||||||
|
//- Write
|
||||||
|
virtual void write(Ostream&) const;
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
Reference in New Issue
Block a user