Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
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.
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
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "copiedFixedValueFvPatchScalarField.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::copiedFixedValueFvPatchScalarField::copiedFixedValueFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchScalarField(p, iF),
sourceFieldName_("default")
{}
Foam::copiedFixedValueFvPatchScalarField::copiedFixedValueFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchScalarField(p, iF, dict),
Andrew Heather
committed
sourceFieldName_(dict.lookup("sourceField"))
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
{}
Foam::copiedFixedValueFvPatchScalarField::copiedFixedValueFvPatchScalarField
(
const copiedFixedValueFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchScalarField(ptf, p, iF, mapper),
sourceFieldName_(ptf.sourceFieldName_)
{}
Foam::copiedFixedValueFvPatchScalarField::copiedFixedValueFvPatchScalarField
(
const copiedFixedValueFvPatchScalarField& awfpsf
)
:
fixedValueFvPatchScalarField(awfpsf),
sourceFieldName_(awfpsf.sourceFieldName_)
{}
Foam::copiedFixedValueFvPatchScalarField::copiedFixedValueFvPatchScalarField
(
const copiedFixedValueFvPatchScalarField& awfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchScalarField(awfpsf, iF),
sourceFieldName_(awfpsf.sourceFieldName_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::copiedFixedValueFvPatchScalarField::updateCoeffs()
{
if (this->updated())
{
return;
}
operator==
(
patch().lookupPatchField<volScalarField, scalar>(sourceFieldName_)
);
fixedValueFvPatchScalarField::updateCoeffs();
}
void Foam::copiedFixedValueFvPatchScalarField::write(Ostream& os) const
{
fvPatchField<scalar>::write(os);
os.writeEntry("sourceField", sourceFieldName_);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
copiedFixedValueFvPatchScalarField
);
}
// ************************************************************************* //