Fixed URL for libccmio-2.6.1 (bug report #5 by Thomas Oliveira)
[foam-extend-3.2.git] / applications / solvers / solidMechanics / elasticPlasticNonLinTLSolidFoam / createFields.H
blob4f7bcc9426260c67a92c3b6c4f6de82b7a231d83
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     Info<< "Creating field U\n" << endl;
18     volVectorField U
19     (
20         IOobject
21         (
22             "U",
23             runTime.timeName(),
24             mesh,
25             IOobject::READ_IF_PRESENT,
26             IOobject::AUTO_WRITE
27         ),
28         mesh,
29         dimensionedVector("zero", dimLength, vector::zero)
30     );
32     volTensorField gradU
33     (
34         IOobject
35         (
36             "grad(U)",
37             runTime.timeName(),
38             mesh,
39             IOobject::READ_IF_PRESENT,
40             IOobject::AUTO_WRITE
41         ),
42         mesh,
43         dimensionedTensor("zero", dimless, tensor::zero)
44     );
46     //- Increment of Green finite strain tensor
47     volSymmTensorField DEpsilon
48     (
49         IOobject
50         (
51             "DEpsilon",
52             runTime.timeName(),
53             mesh,
54             IOobject::READ_IF_PRESENT,
55             IOobject::AUTO_WRITE
56         ),
57         mesh,
58         dimensionedSymmTensor("zero", dimless, symmTensor::zero)
59     );
61     volSymmTensorField epsilon
62     (
63         IOobject
64         (
65             "epsilon",
66             runTime.timeName(),
67             mesh,
68             IOobject::READ_IF_PRESENT,
69             IOobject::AUTO_WRITE
70         ),
71         mesh,
72         dimensionedSymmTensor("zero", dimless, symmTensor::zero)
73     );
75     //- plastic strain
76     volSymmTensorField epsilonP
77     (
78         IOobject
79         (
80             "epsilonP",
81             runTime.timeName(),
82             mesh,
83             IOobject::READ_IF_PRESENT,
84             IOobject::AUTO_WRITE
85         ),
86         mesh,
87         dimensionedSymmTensor("zero", dimless, symmTensor::zero)
88     );
91     //- Increment of 2nd Piola-Kirchhoff stress tensor
92     volSymmTensorField DSigma
93     (
94         IOobject
95         (
96             "DSigma",
97             runTime.timeName(),
98             mesh,
99             IOobject::READ_IF_PRESENT,
100             IOobject::AUTO_WRITE
101         ),
102         mesh,
103         dimensionedSymmTensor("zero", dimForce/dimArea, symmTensor::zero)
104     );
106     //- 2nd Piola-Kirchhoff stress tensor
107     volSymmTensorField sigma
108     (
109         IOobject
110         (
111             "sigma",
112             runTime.timeName(),
113             mesh,
114             IOobject::READ_IF_PRESENT,
115             IOobject::AUTO_WRITE
116         ),
117         mesh,
118         dimensionedSymmTensor("zero", dimForce/dimArea, symmTensor::zero)
119     );
121     volVectorField divDSigmaExp
122     (
123         IOobject
124         (
125             "divDSigmaExp",
126             runTime.timeName(),
127             mesh,
128             IOobject::NO_READ,
129             IOobject::NO_WRITE
130         ),
131         mesh,
132         dimensionedVector("zero", dimensionSet(1, -2, -2, 0, 0, 0, 0), vector::zero)
133     );
135     volVectorField divDSigmaNonLinExp
136     (
137         IOobject
138         (
139             "divDSigmaNonLinExp",
140             runTime.timeName(),
141             mesh,
142             IOobject::NO_READ,
143             IOobject::NO_WRITE
144         ),
145         mesh,
146         dimensionedVector("zero", dimensionSet(1,-2,-2,0,0,0,0), vector::zero)
147     );
149     constitutiveModel rheology(sigma, DU);
151     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     );
177     // aitken relaxation factor
178     scalar aitkenInitialRes = 1.0;
179     scalar aitkenTheta = 0.1;
181 //    volVectorField resid
182 //    (
183 //        IOobject
184 //        (
185 //            "resid",
186 //            runTime.timeName(),
187 //            mesh,
188 //            IOobject::NO_READ,
189 //            IOobject::AUTO_WRITE
190 //        ),
191 //        mesh,
192 //        dimensionedVector("zero", dimless, vector::zero)
193 //    );