solidDisplacementFoam: Added support for fvOptions in both D and T equations
This commit is contained in:
@ -6,4 +6,5 @@ EXE_INC = \
|
||||
|
||||
EXE_LIBS = \
|
||||
-lfiniteVolume \
|
||||
-lfvOptions \
|
||||
-lmeshTools
|
||||
|
||||
@ -1,6 +1,6 @@
|
||||
#include "readMechanicalProperties.H"
|
||||
#include "readThermalProperties.H"
|
||||
|
||||
mesh.Sf();
|
||||
Info<< "Reading field D\n" << endl;
|
||||
volVectorField D
|
||||
(
|
||||
@ -77,3 +77,5 @@ else
|
||||
}
|
||||
|
||||
mesh.setFluxRequired(D.name());
|
||||
|
||||
#include "createFvOptions.H"
|
||||
|
||||
@ -36,6 +36,7 @@ Description
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "fvCFD.H"
|
||||
#include "fvOptions.H"
|
||||
#include "Switch.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
@ -68,10 +69,16 @@ int main(int argc, char *argv[])
|
||||
if (thermalStress)
|
||||
{
|
||||
volScalarField& T = Tptr();
|
||||
solve
|
||||
fvScalarMatrix TEqn
|
||||
(
|
||||
fvm::ddt(T) == fvm::laplacian(DT, T)
|
||||
fvm::ddt(T) == fvm::laplacian(DT, T) + fvOptions(T)
|
||||
);
|
||||
|
||||
fvOptions.constrain(TEqn);
|
||||
|
||||
TEqn.solve();
|
||||
|
||||
fvOptions.correct(T);
|
||||
}
|
||||
|
||||
{
|
||||
@ -81,6 +88,7 @@ int main(int argc, char *argv[])
|
||||
==
|
||||
fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)")
|
||||
+ divSigmaExp
|
||||
+ fvOptions.d2dt2(D)
|
||||
);
|
||||
|
||||
if (thermalStress)
|
||||
@ -89,8 +97,7 @@ int main(int argc, char *argv[])
|
||||
DEqn += fvc::grad(threeKalpha*T);
|
||||
}
|
||||
|
||||
// DEqn.setComponentReference(1, 0, vector::X, 0);
|
||||
// DEqn.setComponentReference(1, 0, vector::Z, 0);
|
||||
fvOptions.constrain(DEqn);
|
||||
|
||||
initialResidual = DEqn.solve().max().initialResidual();
|
||||
|
||||
|
||||
Reference in New Issue
Block a user