calculateAutoCorrelationFunctions.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 if (mesh.time().timeIndex() % vacf.sampleSteps() == 0)
27 {
28  Field<vector> uVals(molecules.size());
29 
30  label uV = 0;
31 
32  forAllConstIter(IDLList<molecule>, molecules, mol)
33  {
34  uVals[uV++] = mol().U();
35  }
36 
37  vacf.calculateCorrelationFunction(uVals);
38 }
39 
40 if (mesh.time().timeIndex() % pacf.sampleSteps() == 0)
41 {
42  vector p = Zero;
43 
44  forAllConstIter(IDLList<molecule>, molecules, mol)
45  {
46  p.x() +=
47  mol().mass() * mol().U().y() * mol().U().z()
48  + 0.5*mol().rf().yz();
49 
50  p.y() +=
51  mol().mass() * mol().U().z() * mol().U().x()
52  + 0.5*mol().rf().zx();
53 
54  p.z() +=
55  mol().mass() * mol().U().x() * mol().U().y()
56  + 0.5*mol().rf().xy();
57  }
58 
59  pacf.calculateCorrelationFunction(p);
60 }
61 
62 if (mesh.time().timeIndex() % hfacf.sampleSteps() == 0)
63 {
64  vector s = Zero;
65 
66  forAllConstIter(IDLList<molecule>, molecules, mol)
67  {
68  s +=
69  (
70  0.5*mol().mass()*magSqr(mol().U())
71  + mol().potentialEnergy()
72  )*mol().U()
73  + 0.5*(mol().rf() & mol().U());
74  }
75 
76  hfacf.calculateCorrelationFunction(s);
77 }
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Info<< tab<< "pressure"<< endl;const dictionary &pressureACFDict(autocorrelationFunctionDict.subDict("pressure"));correlationFunction< vector > pacf(mesh, pressureACFDict, 1)
Pressure autocorrelation function.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
Info<< tab<< "heat flux"<< endl;const dictionary &heatFluxACFDict(autocorrelationFunctionDict.subDict("heatFlux"));correlationFunction< vector > hfacf(mesh, heatFluxACFDict, 1)
Heat flux autocorrelation function.
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Info<< nl<< "Creating autocorrelation functions."<< endl;IOdictionary mdTransportPropertiesDict(IOobject("mdTransportPropertiesDict", mesh.time().system(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE, false));const dictionary &autocorrelationFunctionDict(mdTransportPropertiesDict.subDict("autocorrelationFunctions"));Info<< tab<< "velocity"<< endl;const dictionary &velocityACFDict(autocorrelationFunctionDict.subDict("velocity"));correlationFunction< vector > vacf(mesh, velocityACFDict, molecules.size())
dimensioned< scalar > magSqr(const dimensioned< Type > &)
U
Definition: pEqn.H:72
volScalarField & p