This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/stressAnalysis/icoFsiFoam/tractionDisplacement/tractionDisplacementFvPatchVectorField.C

204 lines
5.7 KiB
C
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright held by original author
\\/ 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 "tractionDisplacementFvPatchVectorField.H"
#include "addToRunTimeSelectionTable.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
tractionDisplacementFvPatchVectorField::
tractionDisplacementFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF
)
:
fixedGradientFvPatchVectorField(p, iF),
traction_(p.size(), vector::zero),
pressure_(p.size(), 0.0)
{
fvPatchVectorField::operator=(patchInternalField());
gradient() = vector::zero;
}
tractionDisplacementFvPatchVectorField::
tractionDisplacementFvPatchVectorField
(
const tractionDisplacementFvPatchVectorField& tdpvf,
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedGradientFvPatchVectorField(tdpvf, p, iF, mapper),
traction_(tdpvf.traction_, mapper),
pressure_(tdpvf.pressure_, mapper)
{}
tractionDisplacementFvPatchVectorField::
tractionDisplacementFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const dictionary& dict
)
:
fixedGradientFvPatchVectorField(p, iF),
traction_("traction", dict, p.size()),
pressure_("pressure", dict, p.size())
{
fvPatchVectorField::operator=(patchInternalField());
gradient() = vector::zero;
}
tractionDisplacementFvPatchVectorField::
tractionDisplacementFvPatchVectorField
(
const tractionDisplacementFvPatchVectorField& tdpvf
)
:
fixedGradientFvPatchVectorField(tdpvf),
traction_(tdpvf.traction_),
pressure_(tdpvf.pressure_)
{}
tractionDisplacementFvPatchVectorField::
tractionDisplacementFvPatchVectorField
(
const tractionDisplacementFvPatchVectorField& tdpvf,
const DimensionedField<vector, volMesh>& iF
)
:
fixedGradientFvPatchVectorField(tdpvf, iF),
traction_(tdpvf.traction_),
pressure_(tdpvf.pressure_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void tractionDisplacementFvPatchVectorField::autoMap
(
const fvPatchFieldMapper& m
)
{
fixedGradientFvPatchVectorField::autoMap(m);
traction_.autoMap(m);
pressure_.autoMap(m);
}
// Reverse-map the given fvPatchField onto this fvPatchField
void tractionDisplacementFvPatchVectorField::rmap
(
const fvPatchVectorField& ptf,
const labelList& addr
)
{
fixedGradientFvPatchVectorField::rmap(ptf, addr);
const tractionDisplacementFvPatchVectorField& dmptf =
refCast<const tractionDisplacementFvPatchVectorField>(ptf);
traction_.rmap(dmptf.traction_, addr);
pressure_.rmap(dmptf.pressure_, addr);
}
// Update the coefficients associated with the patch field
void tractionDisplacementFvPatchVectorField::updateCoeffs()
{
if (updated())
{
return;
}
const dictionary& mechanicalProperties =
db().lookupObject<IOdictionary>("mechanicalProperties");
dimensionedScalar rho(mechanicalProperties.lookup("rho"));
dimensionedScalar rhoE(mechanicalProperties.lookup("E"));
dimensionedScalar nu(mechanicalProperties.lookup("nu"));
dimensionedScalar E = rhoE/rho;
dimensionedScalar mu = E/(2.0*(1.0 + nu));
dimensionedScalar lambda = nu*E/((1.0 + nu)*(1.0 - 2.0*nu));
dimensionedScalar threeK = E/(1.0 - 2.0*nu);
Switch planeStress(mechanicalProperties.lookup("planeStress"));
if (planeStress)
{
lambda = nu*E/((1.0 + nu)*(1.0 - nu));
threeK = E/(1.0 - nu);
}
vectorField n = patch().nf();
const fvPatchField<tensor>& gradU =
patch().lookupPatchField<volTensorField, tensor>("grad(U)");
gradient() =
(
(traction_ - pressure_*n)/rho.value()
- (n & (mu.value()*gradU.T() - (mu + lambda).value()*gradU))
- n*tr(gradU)*lambda.value()
)/(2.0*mu + lambda).value();
fixedGradientFvPatchVectorField::updateCoeffs();
}
// Write
void tractionDisplacementFvPatchVectorField::write(Ostream& os) const
{
fvPatchVectorField::write(os);
traction_.writeEntry("traction", os);
pressure_.writeEntry("pressure", os);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePatchTypeField(fvPatchVectorField, tractionDisplacementFvPatchVectorField);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// ************************************************************************* //