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/solidMechanics/elasticPlasticNonLinULSolidFoam/createFields.H

181 lines
4 KiB
C
Raw Normal View History

Info<< "Reading incremental displacement field DU\n" << endl;
2012-09-11 15:42:55 +00:00
volVectorField DU
(
IOobject
(
"DU",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
volTensorField gradDU = fvc::grad(DU);
Info<< "Reading accumulated displacement field U\n" << endl;
2012-09-11 15:42:55 +00:00
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedVector("zero", dimLength, vector::zero)
);
Info << "Reading accumulated strain field epsilon\n" << endl;
volSymmTensorField epsilon
2012-09-11 15:42:55 +00:00
(
IOobject
(
"epsilon",
2012-09-11 15:42:55 +00:00
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedSymmTensor("zero", dimless, symmTensor::zero)
2012-09-11 15:42:55 +00:00
);
// total plastic strain
volSymmTensorField epsilonP
2012-09-11 15:42:55 +00:00
(
IOobject
(
"epsilonP",
2012-09-11 15:42:55 +00:00
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedSymmTensor("zero", dimless, symmTensor::zero)
);
2012-09-11 15:42:55 +00:00
volSymmTensorField DEpsilon
2012-09-11 15:42:55 +00:00
(
IOobject
(
"DEpsilon",
2012-09-11 15:42:55 +00:00
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedSymmTensor("zero", dimless, symmTensor::zero)
);
Info << "Reading accumulated stress field sigma\n" << endl;
volSymmTensorField sigma
2012-09-11 15:42:55 +00:00
(
IOobject
(
"sigma",
2012-09-11 15:42:55 +00:00
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedSymmTensor("zero", dimForce/dimArea, symmTensor::zero)
2012-09-11 15:42:55 +00:00
);
Info << "Reading incremental stress field DSigma\n" << endl;
volSymmTensorField DSigma
2012-09-11 15:42:55 +00:00
(
IOobject
(
"DSigma",
2012-09-11 15:42:55 +00:00
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedSymmTensor("zero", dimForce/dimArea, symmTensor::zero)
2012-09-11 15:42:55 +00:00
);
//- explicit terms in the momentum equation
volVectorField divDSigmaExp
(
IOobject
(
"divDSigmaExp",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionedVector("zero", dimensionSet(1, -2, -2, 0, 0, 0, 0), vector::zero)
);
volVectorField divDSigmaNonLinExp
(
2012-09-11 15:42:55 +00:00
IOobject
(
"divDSigmaNonLinExp",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
2012-09-11 15:42:55 +00:00
mesh,
dimensionedVector("zero", dimensionSet(1,-2,-2,0,0,0,0), vector::zero)
);
2012-09-11 15:42:55 +00:00
constitutiveModel rheology(sigma, DU);
2012-09-11 15:42:55 +00:00
volScalarField rho = rheology.rho();
volScalarField mu = rheology.mu();
volScalarField lambda = rheology.lambda();
surfaceScalarField muf = fvc::interpolate(rheology.mu());
surfaceScalarField lambdaf = fvc::interpolate(rheology.lambda());
2012-09-11 15:42:55 +00:00
surfaceVectorField n = mesh.Sf()/mesh.magSf();
// plastic strain increment
const volSymmTensorField& DEpsilonP = rheology.DEpsilonP();
// for aitken relaxation
volVectorField aitkenDelta
(
2012-09-11 15:42:55 +00:00
IOobject
(
"aitkenDelta",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionedVector("zero", dimLength, vector::zero)
);
// aitken relaxation factor
scalar aitkenInitialRes = 1.0;
scalar aitkenTheta = 0.1;
// volVectorField resid
// (
// IOobject
// (
// "resid",
// runTime.timeName(),
// mesh,
// IOobject::NO_READ,
// IOobject::AUTO_WRITE
// ),
// mesh,
// dimensionedVector("zero", dimless, vector::zero)
// );