velocityComponentLaplacianFvMotionSolver.C 4.17 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2012 OpenFOAM Foundation
6
7
8
9
10
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

11
12
13
14
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
15
16
17
18
19
20
21

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    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
27
28
29

\*---------------------------------------------------------------------------*/

#include "velocityComponentLaplacianFvMotionSolver.H"
#include "motionDiffusivity.H"
#include "fvmLaplacian.H"
#include "addToRunTimeSelectionTable.H"
30
#include "volPointInterpolation.H"
31
32
33
34
35
36
37
38
39

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

namespace Foam
{
    defineTypeNameAndDebug(velocityComponentLaplacianFvMotionSolver, 0);

    addToRunTimeSelectionTable
    (
40
        motionSolver,
41
42
43
44
45
46
47
48
49
50
51
52
        velocityComponentLaplacianFvMotionSolver,
        dictionary
    );
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::velocityComponentLaplacianFvMotionSolver::
velocityComponentLaplacianFvMotionSolver
(
    const polyMesh& mesh,
53
    const IOdictionary& dict
54
55
)
:
56
57
    componentVelocityMotionSolver(mesh, dict, typeName),
    fvMotionSolverCore(mesh),
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
    cellMotionU_
    (
        IOobject
        (
            "cellMotionU" + cmptName_,
            mesh.time().timeName(),
            mesh,
            IOobject::READ_IF_PRESENT,
            IOobject::AUTO_WRITE
        ),
        fvMesh_,
        dimensionedScalar
        (
            "cellMotionU",
            pointMotionU_.dimensions(),
            0
        ),
        cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
    ),
    diffusivityPtr_
    (
79
        motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
80
    )
81
{}
82
83
84
85
86
87
88
89
90
91
92
93
94
95


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

Foam::velocityComponentLaplacianFvMotionSolver::
~velocityComponentLaplacianFvMotionSolver()
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

Foam::tmp<Foam::pointField>
Foam::velocityComponentLaplacianFvMotionSolver::curPoints() const
{
96
97
98
99
100
    volPointInterpolation::New(fvMesh_).interpolate
    (
        cellMotionU_,
        pointMotionU_
    );
101
102
103
104
105
106
107

    tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));

    tcurPoints().replace
    (
        cmpt_,
        tcurPoints().component(cmpt_)
graham's avatar
graham committed
108
      + fvMesh_.time().deltaTValue()*pointMotionU_.internalField()
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
    );

    twoDCorrectPoints(tcurPoints());

    return tcurPoints;
}


void Foam::velocityComponentLaplacianFvMotionSolver::solve()
{
    // The points have moved so before interpolation update
    // the fvMotionSolver accordingly
    movePoints(fvMesh_.points());

    diffusivityPtr_->correct();
    pointMotionU_.boundaryField().updateCoeffs();

    Foam::solve
    (
        fvm::laplacian
        (
            diffusivityPtr_->operator()(),
            cellMotionU_,
            "laplacian(diffusivity,cellMotionU)"
        )
    );
}


void Foam::velocityComponentLaplacianFvMotionSolver::updateMesh
(
    const mapPolyMesh& mpm
)
{
143
    componentVelocityMotionSolver::updateMesh(mpm);
144
145
146
147

    // Update diffusivity. Note two stage to make sure old one is de-registered
    // before creating/registering new one.
    diffusivityPtr_.reset(NULL);
148
149
150
151
152
    diffusivityPtr_ = motionDiffusivity::New
    (
        fvMesh_,
        coeffDict().lookup("diffusivity")
    );
153
154
155
156
}


// ************************************************************************* //