ATCstandard.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) 2007-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "ATCstandard.H"
32#include "wallFvPatch.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
43(
47);
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
53(
54 const fvMesh& mesh,
55 const incompressibleVars& primalVars,
56 const incompressibleAdjointVars& adjointVars,
57 const dictionary& dict
58)
59:
60 ATCModel(mesh, primalVars, adjointVars, dict),
61 gradU_
62 (
64 (
65 "gradUATC",
66 mesh_.time().timeName(),
67 mesh_,
68 IOobject::NO_READ,
69 IOobject::NO_WRITE
70 ),
71 mesh_,
73 )
74{}
75
76
77// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78
80{
81 addProfiling(ATCstandard, "ATCstandard::addATC");
82 const volVectorField& U = primalVars_.U();
85
86
87 // Main ATC term
88 ATC_ = gradU_ & Ua;
89
90 if (extraConvection_ > 0)
91 {
92 // Implicit part added to increase diagonal dominance
93 // Note: Maybe this needs to be multiplied with the ATClimiter ??
94 UaEqn += extraConvection_*fvm::div(-phi, Ua);
95
96 // correct rhs due to implicitly augmenting the adjoint convection
97 ATC_ += extraConvection_*(fvc::grad(Ua, "gradUaATC")().T() & U);
98 }
99
100 //zero ATC on cells next to given patch types
101 smoothATC();
102
103 // actual ATC term
104 UaEqn += ATC_.internalField();
105}
106
107
109{
110 const volVectorField& U = primalVars_.U();
111 const volVectorField& Ua = adjointVars_.Ua();
112
113 // We only need to modify the boundaryField of gradU locally.
114 // If grad(U) is cached then
115 // a. The .ref() call fails since the tmp is initialised from a
116 // const ref
117 // b. we would be changing grad(U) for all other places in the code
118 // that need it
119 // So, always allocate new memory and avoid registering the new field
120 tmp<volTensorField> tgradU(volTensorField::New("gradULocal", fvc::grad(U)));
121 volTensorField::Boundary& gradUbf = tgradU.ref().boundaryFieldRef();
122
123 // Explicitly correct the boundary gradient to get rid of
124 // the tangential component
125 forAll(mesh_.boundary(), patchI)
126 {
127 const fvPatch& patch = mesh_.boundary()[patchI];
128 if (isA<wallFvPatch>(patch))
129 {
130 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
131 gradUbf[patchI] = tnf*U.boundaryField()[patchI].snGrad();
132 }
133 }
134
135 const volScalarField& mask = getLimiter();
136
137 return
139 (
140 "ATCFISensitivityTerm" + type(),
141 - (tgradU & Ua)*U*mask
142 );
143}
144
145
147{
148 const volVectorField& U = primalVars_.U();
150 // Build U to go into the ATC term, based on whether to smooth field or not
151 autoPtr<volVectorField> UForATC(nullptr);
153 {
154 gradU_ = fvc::grad(fvc::reconstruct(phi), "gradUATC");
155 }
156 else
157 {
158 gradU_ = fvc::grad(U, "gradUATC");
159 }
160}
161
162
163// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164
165} // End namespace Foam
166
167// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition: ATCModel.H:63
const incompressibleAdjointVars & adjointVars_
Definition: ATCModel.H:81
const fvMesh & mesh_
Definition: ATCModel.H:79
const scalar extraConvection_
Definition: ATCModel.H:85
volVectorField ATC_
Definition: ATCModel.H:92
const bool reconstructGradients_
Definition: ATCModel.H:88
const volScalarField & getLimiter() const
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:180
void smoothATC()
Limit ATC field using ATClimiter_.
Definition: ATCModel.C:51
const incompressibleVars & primalVars_
Definition: ATCModel.H:80
The ATC formualtion resulting by differentiating the Non-conservative form of the momentum equations.
Definition: ATCstandard.H:57
virtual void updatePrimalBasedQuantities()
Update quantities related with the primal fields.
Definition: ATCstandard.C:146
virtual tmp< volTensorField > getFISensitivityTerm() const
Get the FI sensitivity derivatives term coming from the ATC.
Definition: ATCstandard.C:108
virtual void addATC(fvVectorMatrix &UaEqn)
Add ATC.
Definition: ATCstandard.C:79
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const volVectorField & Ua() const
Return const reference to velocity.
const volVectorField & UaInst() const
Return const reference to velocity.
Class including all adjoint fields for incompressible flows.
Base class for solution control classes.
const volVectorField & U() const
Return const reference to velocity.
const surfaceScalarField & phi() const
Return const reference to volume flux.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const volScalarField & T
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333