2013-10-11 13:31:14 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
|
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
|
|
\\ / O peration |
|
|
|
|
\\ / A nd | Copyright (C) 2004-2007 Hrvoje Jasak
|
|
|
|
\\/ 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., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
|
|
|
|
|
|
|
|
Application
|
|
|
|
elasticIncrAcpSolidFoam
|
|
|
|
|
|
|
|
Description
|
|
|
|
Incremental form of elasticAcpSolidFoam
|
|
|
|
arbitrary crack propagation solver
|
|
|
|
|
|
|
|
Author
|
|
|
|
Zeljko Tukovic, FSB Zagreb
|
|
|
|
Declan Carolan UCD
|
|
|
|
Philip Cardiff UCD
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "fvCFD.H"
|
|
|
|
#include "constitutiveModel.H"
|
2013-10-14 08:26:40 +00:00
|
|
|
//#include "componentReferenceList.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
#include "crackerFvMesh.H"
|
|
|
|
#include "processorPolyPatch.H"
|
|
|
|
#include "SortableList.H"
|
|
|
|
#include "solidInterface.H"
|
|
|
|
#include "solidCohesiveFvPatchVectorField.H"
|
|
|
|
#include "solidCohesiveFixedModeMixFvPatchVectorField.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
# include "setRootCase.H"
|
|
|
|
# include "createTime.H"
|
|
|
|
# include "createCrackerMesh.H"
|
|
|
|
# include "createFields.H"
|
|
|
|
# include "createCrack.H"
|
|
|
|
//# include "createReference.H"
|
|
|
|
# include "createHistory.H"
|
|
|
|
# include "readDivDSigmaExpMethod.H"
|
|
|
|
# include "createSolidInterfaceIncrNoModify.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
lduMatrix::debug = 0;
|
|
|
|
|
|
|
|
scalar maxEffTractionFraction = 0;
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
while (runTime.run())
|
|
|
|
{
|
2013-10-14 08:26:40 +00:00
|
|
|
# include "readSolidMechanicsControls.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "setDeltaT.H"
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
runTime++;
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
Info<< "\nTime: " << runTime.timeName() << " s\n" << endl;
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
volScalarField rho = rheology.rho();
|
|
|
|
volScalarField mu = rheology.mu();
|
|
|
|
volScalarField lambda = rheology.lambda();
|
|
|
|
surfaceScalarField muf = fvc::interpolate(mu);
|
|
|
|
surfaceScalarField lambdaf = fvc::interpolate(lambda);
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
solidInterfacePtr->modifyProperties(muf, lambdaf);
|
|
|
|
//# include "waveCourantNo.H"
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
int iCorr = 0;
|
|
|
|
lduMatrix::solverPerformance solverPerf;
|
|
|
|
scalar initialResidual = 0;
|
|
|
|
scalar relativeResidual = 1;
|
|
|
|
//scalar forceResidual = 1;
|
|
|
|
label nFacesToBreak = 0;
|
|
|
|
label nCoupledFacesToBreak = 0;
|
|
|
|
bool topoChange = false;
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
// DU from the previous timestep is usually a good guess
|
|
|
|
// for the next timestep, but it can cause faces to prematurely
|
|
|
|
// crack.
|
|
|
|
// so I will reduce DU here to stop this happening
|
2013-11-07 20:22:30 +00:00
|
|
|
if (!predictor)
|
|
|
|
{
|
|
|
|
DU *= 0.0;
|
|
|
|
}
|
2013-10-11 13:31:14 +00:00
|
|
|
|
|
|
|
do
|
2013-11-07 20:22:30 +00:00
|
|
|
{
|
|
|
|
surfaceVectorField n = mesh.Sf()/mesh.magSf();
|
|
|
|
do
|
|
|
|
{
|
|
|
|
DU.storePrevIter();
|
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "calculateDivDSigmaExp.H"
|
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
fvVectorMatrix DUEqn
|
|
|
|
(
|
|
|
|
rho*fvm::d2dt2(DU)
|
|
|
|
==
|
|
|
|
fvm::laplacian(2*muf + lambdaf, DU, "laplacian(DDU,DU)")
|
|
|
|
+ divDSigmaExp
|
|
|
|
);
|
|
|
|
|
|
|
|
//# include "setReference.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
if(solidInterfacePtr)
|
|
|
|
{
|
|
|
|
solidInterfacePtr->correct(DUEqn);
|
|
|
|
}
|
2013-10-11 13:31:14 +00:00
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
//DUEqn.relax();
|
2013-10-11 13:31:14 +00:00
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
solverPerf = DUEqn.solve();
|
|
|
|
|
|
|
|
if (aitkenRelax)
|
|
|
|
{
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "aitkenRelaxation.H"
|
2013-11-07 20:22:30 +00:00
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
DU.relax();
|
|
|
|
}
|
|
|
|
|
|
|
|
if (iCorr == 0)
|
|
|
|
{
|
|
|
|
initialResidual = solverPerf.initialResidual();
|
|
|
|
aitkenInitialRes = gMax(mag(DU.internalField()));
|
|
|
|
}
|
|
|
|
|
|
|
|
//gradDU = solidInterfacePtr->grad(DU);
|
|
|
|
// use leastSquaresSolidInterface grad scheme
|
|
|
|
gradDU = fvc::grad(DU);
|
2013-10-11 13:31:14 +00:00
|
|
|
|
|
|
|
# include "calculateRelativeResidual.H"
|
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
if (iCorr % infoFrequency == 0)
|
|
|
|
{
|
|
|
|
Info << "\tTime " << runTime.value()
|
|
|
|
<< ", Corr " << iCorr
|
|
|
|
<< ", Solving for " << DU.name()
|
|
|
|
<< " using " << solverPerf.solverName()
|
|
|
|
<< ", res = " << solverPerf.initialResidual()
|
|
|
|
<< ", rel res = " << relativeResidual;
|
|
|
|
if (aitkenRelax)
|
|
|
|
{
|
|
|
|
Info << ", aitken = " << aitkenTheta;
|
|
|
|
}
|
|
|
|
Info << ", inner iters " << solverPerf.nIterations() << endl;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
while
|
|
|
|
(
|
|
|
|
//iCorr++ == 0
|
|
|
|
iCorr++ < 10
|
|
|
|
||
|
|
|
|
(
|
|
|
|
//solverPerf.initialResidual() > convergenceTolerance
|
|
|
|
relativeResidual > convergenceTolerance
|
|
|
|
&&
|
|
|
|
iCorr < nCorr
|
|
|
|
)
|
|
|
|
);
|
|
|
|
|
|
|
|
Info<< "Solving for " << DU.name() << " using "
|
|
|
|
<< solverPerf.solverName()
|
|
|
|
<< ", Initial residual = " << initialResidual
|
|
|
|
<< ", Final residual = " << solverPerf.initialResidual()
|
|
|
|
<< ", No outer iterations " << iCorr
|
|
|
|
<< ", Relative residual " << relativeResidual << endl;
|
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "calculateTraction.H"
|
|
|
|
# include "updateCrack.H"
|
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
Info<< "Max effective traction fraction: "
|
|
|
|
<< maxEffTractionFraction << endl;
|
2013-10-11 13:31:14 +00:00
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
// reset counter if faces want to crack
|
|
|
|
if ((nFacesToBreak > 0) || (nCoupledFacesToBreak > 0)) iCorr = 0;
|
|
|
|
}
|
2013-10-11 13:31:14 +00:00
|
|
|
while( (nFacesToBreak > 0) || (nCoupledFacesToBreak > 0));
|
2013-11-07 20:22:30 +00:00
|
|
|
|
|
|
|
if (cohesivePatchDUPtr)
|
|
|
|
{
|
|
|
|
if (returnReduce(cohesivePatchDUPtr->size(), sumOp<label>()))
|
|
|
|
{
|
|
|
|
cohesivePatchDUPtr->cracking();
|
|
|
|
}
|
|
|
|
}
|
2013-10-11 13:31:14 +00:00
|
|
|
else
|
2013-11-07 20:22:30 +00:00
|
|
|
{
|
|
|
|
if
|
|
|
|
(
|
|
|
|
returnReduce
|
|
|
|
(
|
|
|
|
cohesivePatchDUFixedModePtr->size(),
|
|
|
|
sumOp<label>())
|
|
|
|
)
|
|
|
|
{
|
|
|
|
Pout << "Number of faces in crack: "
|
|
|
|
<< cohesivePatchDUFixedModePtr->size() << endl;
|
|
|
|
cohesivePatchDUFixedModePtr->relativeSeparationDistance();
|
|
|
|
}
|
|
|
|
}
|
2013-10-11 13:31:14 +00:00
|
|
|
|
|
|
|
# include "calculateDEpsilonDSigma.H"
|
|
|
|
|
|
|
|
// update total quantities
|
|
|
|
U += DU;
|
|
|
|
epsilon += DEpsilon;
|
|
|
|
sigma += DSigma;
|
|
|
|
|
|
|
|
# include "writeFields.H"
|
|
|
|
# include "writeHistory.H"
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
2013-11-07 20:22:30 +00:00
|
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s\n\n"
|
|
|
|
<< endl;
|
2013-10-11 13:31:14 +00:00
|
|
|
}
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
Info<< "End\n" << endl;
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
return(0);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|