oversetFvPatchField.C 6.86 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2016-2018 OpenCFD Ltd.
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
     \\/     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 "volFields.H"
27
#include "cellCellStencil.H"
28
#include "cellCellStencilObject.H"
29
#include "dynamicOversetFvMesh.H"
30
31
32
33
34
35
36
37
38
39

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

template<class Type>
Foam::oversetFvPatchField<Type>::oversetFvPatchField
(
    const fvPatch& p,
    const DimensionedField<Type, volMesh>& iF
)
:
40
41
    zeroGradientFvPatchField<Type>(p, iF),
    oversetPatch_(refCast<const oversetFvPatch>(p))
42
43
44
45
46
47
48
49
50
51
52
53
{}


template<class Type>
Foam::oversetFvPatchField<Type>::oversetFvPatchField
(
    const oversetFvPatchField<Type>& ptf,
    const fvPatch& p,
    const DimensionedField<Type, volMesh>& iF,
    const fvPatchFieldMapper& mapper
)
:
54
55
    zeroGradientFvPatchField<Type>(ptf, p, iF, mapper),
    oversetPatch_(refCast<const oversetFvPatch>(p))
56
{}
57
58
59
60
61
62
63
64
65
66


template<class Type>
Foam::oversetFvPatchField<Type>::oversetFvPatchField
(
    const fvPatch& p,
    const DimensionedField<Type, volMesh>& iF,
    const dictionary& dict
)
:
67
68
    zeroGradientFvPatchField<Type>(p, iF, dict),
    oversetPatch_(refCast<const oversetFvPatch>(p, dict))
69
{}
70
71
72
73
74
75
76
77


template<class Type>
Foam::oversetFvPatchField<Type>::oversetFvPatchField
(
    const oversetFvPatchField<Type>& ptf
)
:
78
79
    zeroGradientFvPatchField<Type>(ptf),
    oversetPatch_(ptf.oversetPatch_)
80
81
82
83
84
85
86
87
88
89
{}


template<class Type>
Foam::oversetFvPatchField<Type>::oversetFvPatchField
(
    const oversetFvPatchField<Type>& ptf,
    const DimensionedField<Type, volMesh>& iF
)
:
90
91
    zeroGradientFvPatchField<Type>(ptf, iF),
    oversetPatch_(ptf.oversetPatch_)
92
93
94
95
96
97
98
99
100
101
102
{}


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

template<class Type>
void Foam::oversetFvPatchField<Type>::initEvaluate
(
    const Pstream::commsTypes commsType
)
{
103
    if (this->oversetPatch_.master())
104
105
106
107
108
109
110
111
    {
        // Trigger interpolation
        const fvMesh& mesh = this->internalField().mesh();
        const dictionary& fvSchemes = mesh.schemesDict();
        const word& fldName = this->internalField().name();

        if (&mesh.lduAddr() != &mesh.fvMesh::lduAddr())
        {
112
113
114
            // Running extended addressing. Flux correction already done
            // in the linear solver (through the initUpdateInterfaceMatrix
            // method below)
115
116
            if (debug)
            {
117
118
                Info<< "Skipping overset interpolation for solved-for field "
                    << fldName << endl;
119
120
            }
        }
121
        else if (!fvSchemes.found("oversetInterpolation"))
122
123
        {
            IOWarningInFunction(fvSchemes)
124
125
                << "Missing required dictionary entry"
                << " 'oversetInterpolation'"
126
127
128
                << ". Skipping overset interpolation for field "
                << fldName << endl;
        }
129
        else if (fvSchemes.found("oversetInterpolationRequired"))
130
        {
131
132
            // Backwards compatibility mode: only interpolate what is
            // explicitly mentioned
133

134
            if (fvSchemes.found("oversetInterpolationSuppressed"))
135
            {
136
137
138
139
140
                FatalIOErrorInFunction(fvSchemes)
                    << "Cannot have both dictionary entry"
                    << " 'oversetInterpolationSuppresed' and "
                    << " 'oversetInterpolationRequired' for field "
                    << fldName << exit(FatalIOError);
141
            }
142
143
144
145
146
            const dictionary& intDict = fvSchemes.subDict
            (
                "oversetInterpolationRequired"
            );
            if (intDict.found(fldName))
147
148
149
            {
                if (debug)
                {
150
                    Info<< "Interpolating field " << fldName << endl;
151
                }
152
153
154
155
156
157

                // Interpolate without bc update (since would trigger infinite
                // recursion back to oversetFvPatchField<Type>::evaluate)
                // The only problem is bcs that use the new cell values
                // (e.g. zeroGradient, processor). These need to appear -after-
                // the 'overset' bc.
158
159
160
161
162
163
164
165
166
167
                mesh.interpolate
                (
                    const_cast<Field<Type>&>
                    (
                        this->primitiveField()
                    )
                );
            }
            else if (debug)
            {
168
                Info<< "Skipping overset interpolation for field "
169
170
171
                    << fldName << endl;
            }
        }
172
173
        else
        {
174
            const dictionary* dictPtr
175
            (
176
                fvSchemes.findDict("oversetInterpolationSuppressed")
177
            );
178
179
180
181
182
183
184
185
186
187
188
189
190
191

            const wordHashSet& suppress =
                Stencil::New(mesh).nonInterpolatedFields();

            bool skipInterpolate = suppress.found(fldName);

            if (dictPtr)
            {
                skipInterpolate =
                    skipInterpolate
                 || dictPtr->found(fldName);
            }

            if (skipInterpolate)
192
193
194
            {
                if (debug)
                {
195
196
197
198
199
200
201
202
203
204
                    Info<< "Skipping suppressed overset interpolation"
                        << " for field " << fldName << endl;
                }
            }
            else
            {
                if (debug)
                {
                    Info<< "Interpolating non-suppressed field " << fldName
                        << endl;
205
206
207
208
209
210
211
212
213
214
215
216
217
218
                }
                mesh.interpolate
                (
                    const_cast<Field<Type>&>
                    (
                        this->primitiveField()
                    )
                );
            }
        }
    }
}


219
220
221
222
223
224
225
226
227
template<class Type>
void Foam::oversetFvPatchField<Type>::write(Ostream& os) const
{
    zeroGradientFvPatchField<Type>::write(os);
    // Make sure to write the value for ease of postprocessing.
    this->writeEntry("value", os);
}


228
// ************************************************************************* //