Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / solidMechanics / elasticPlasticNonLinTLSolidFoam / createFields.H
blob1f615df600049aa76ca287e71c7b92316b5c835e
1     Info<< "Reading field DU\n" << endl;
2     volVectorField DU
3     (
4         IOobject
5         (
6             "DU",
7             runTime.timeName(),
8             mesh,
9             IOobject::MUST_READ,
10             IOobject::AUTO_WRITE
11         ),
12         mesh
13     );
15     volTensorField gradDU = fvc::grad(DU);
17     volVectorField U
18     (
19         IOobject
20         (
21             "U",
22             runTime.timeName(),
23             mesh,
24             IOobject::READ_IF_PRESENT,
25             IOobject::AUTO_WRITE
26         ),
27         mesh,
28         dimensionedVector("zero", dimLength, vector::zero)
29     );
31     volTensorField gradU
32     (
33         IOobject
34         (
35          "grad(U)",
36          runTime.timeName(),
37          mesh,
38          IOobject::READ_IF_PRESENT,
39          IOobject::AUTO_WRITE
40          ),
41         mesh,
42         dimensionedTensor("zero", dimless, tensor::zero)
43      );
45     //- Increment of Green finite strain tensor
46     volSymmTensorField DEpsilon
47     (
48         IOobject
49         (
50             "DEpsilon",
51             runTime.timeName(),
52             mesh,
53             IOobject::READ_IF_PRESENT,
54             IOobject::AUTO_WRITE
55         ),
56         mesh,
57         dimensionedSymmTensor("zero", dimless, symmTensor::zero)
58     );
60     volSymmTensorField epsilon
61     (
62         IOobject
63         (
64             "epsilon",
65             runTime.timeName(),
66             mesh,
67             IOobject::READ_IF_PRESENT,
68             IOobject::AUTO_WRITE
69         ),
70         mesh,
71         dimensionedSymmTensor("zero", dimless, symmTensor::zero)
72     );
74     //- plastic strain
75     volSymmTensorField epsilonP
76     (
77         IOobject
78         (
79             "epsilonP",
80             runTime.timeName(),
81             mesh,
82             IOobject::READ_IF_PRESENT,
83             IOobject::AUTO_WRITE
84         ),
85         mesh,
86         dimensionedSymmTensor("zero", dimless, symmTensor::zero)
87     );
90     //- Increment of 2nd Piola-Kirchhoff stress tensor
91     volSymmTensorField DSigma
92     (
93         IOobject
94         (
95             "DSigma",
96             runTime.timeName(),
97             mesh,
98             IOobject::READ_IF_PRESENT,
99             IOobject::AUTO_WRITE
100         ),
101         mesh,
102         dimensionedSymmTensor("zero", dimForce/dimArea, symmTensor::zero)
103     );
105     //- 2nd Piola-Kirchhoff stress tensor
106     volSymmTensorField sigma
107     (
108         IOobject
109         (
110             "sigma",
111             runTime.timeName(),
112             mesh,
113             IOobject::READ_IF_PRESENT,
114             IOobject::AUTO_WRITE
115         ),
116         mesh,
117         dimensionedSymmTensor("zero", dimForce/dimArea, symmTensor::zero)
118     );
120       volVectorField divDSigmaExp
121       (
122         IOobject
123         (
124          "divDSigmaExp",
125          runTime.timeName(),
126          mesh,
127          IOobject::NO_READ,
128          IOobject::NO_WRITE
129          ),
130         mesh,
131         dimensionedVector("zero", dimensionSet(1,-2,-2,0,0,0,0), vector::zero)
132        );
134       volVectorField divDSigmaNonLinExp
135       (
136         IOobject
137         (
138          "divDSigmaNonLinExp",
139          runTime.timeName(),
140          mesh,
141          IOobject::NO_READ,
142          IOobject::NO_WRITE
143          ),
144         mesh,
145         dimensionedVector("zero", dimensionSet(1,-2,-2,0,0,0,0), vector::zero)
146        );
148     constitutiveModel rheology(sigma, DU);
149     
150     volScalarField rho = rheology.rho();
152     volScalarField mu = rheology.mu();
153     volScalarField lambda = rheology.lambda();
154     surfaceScalarField muf = fvc::interpolate(mu, "mu");
155     surfaceScalarField lambdaf = fvc::interpolate(lambda, "lambda");
157     surfaceVectorField n = mesh.Sf()/mesh.magSf();
159    // plastic strain increment
160    const volSymmTensorField& DEpsilonP = rheology.DEpsilonP();
162    // for aitken relaxation
163    volVectorField aitkenDelta
164     (
165         IOobject
166         (
167             "aitkenDelta",
168             runTime.timeName(),
169             mesh,
170             IOobject::NO_READ,
171             IOobject::NO_WRITE
172         ),
173         mesh,
174         dimensionedVector("zero", dimLength, vector::zero)
175     );
176 // aitken relaxation factor
177 scalar aitkenInitialRes = 1.0;
178 scalar aitkenTheta = 0.1;
180 //    volVectorField resid
181 //    (
182 //         IOobject
183 //         (
184 //          "resid",
185 //          runTime.timeName(),
186 //          mesh,
187 //          IOobject::NO_READ,
188 //          IOobject::AUTO_WRITE
189 //          ),
190 //         mesh,
191 //         dimensionedVector("zero", dimless, vector::zero)
192 //     );