motionDirectionalDiffusivity.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2012 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "surfaceInterpolate.H"
32 #include "velocityMotionSolver.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(motionDirectionalDiffusivity, 0);
39 
41  (
42  motionDiffusivity,
43  motionDirectionalDiffusivity,
44  Istream
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 Foam::motionDirectionalDiffusivity::motionDirectionalDiffusivity
52 (
53  const fvMesh& mesh,
54  Istream& mdData
55 )
56 :
58  diffusivityVector_(mdData)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 {
72  static bool first = true;
73 
74  if (!first)
75  {
76  const volVectorField& cellMotionU =
77  mesh().lookupObject<volVectorField>("cellMotionU");
78 
80  (
81  IOobject
82  (
83  "D",
84  mesh().time().timeName(),
85  mesh()
86  ),
87  diffusivityVector_.y()*vector::one
88  + (diffusivityVector_.x() - diffusivityVector_.y())*cellMotionU
89  /(mag(cellMotionU) + dimensionedScalar("small", dimVelocity, SMALL)),
90 
92  );
93  D.correctBoundaryConditions();
94 
95  const surfaceVectorField n(mesh().Sf()/mesh().magSf());
96  faceDiffusivity_ == (n & cmptMultiply(fvc::interpolate(D), n));
97  }
98  else
99  {
100  first = false;
101 
102  const velocityMotionSolver& mSolver =
103  mesh().lookupObject<velocityMotionSolver>("dynamicMeshDict");
104 
105  const_cast<velocityMotionSolver&>(mSolver).solve();
106  correct();
107  }
108 }
109 
110 
111 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
Uniform uniform finite volume mesh motion diffusivity.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
virtual void correct()
Correct the motion diffusivity.
Macros for easy insertion into run-time selection tables.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
defineTypeNameAndDebug(combustionModel, 0)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label n
const dimensionedScalar & D
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
const dimensionSet dimVelocity