ATCUaGradU.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 "ATCUaGradU.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
42(
46);
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52(
53 const fvMesh& mesh,
54 const incompressibleVars& primalVars,
55 const incompressibleAdjointVars& adjointVars,
56 const dictionary& dict
57)
58:
59 ATCModel(mesh, primalVars, adjointVars, dict)
60{}
61
62
63// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64
66{
67 const volVectorField& U = primalVars_.U();
71
72 // Build Ua to go into the ATC term, based on whether to smooth field or not
73 autoPtr<volVectorField> UaForATC(nullptr);
75 {
76 UaForATC.reset(new volVectorField(fvc::reconstruct(phia)));
77 }
78 else
79 {
80 UaForATC.reset(new volVectorField(Ua));
81 }
82
83 // Main ATC term
84 ATC_ = -fvc::grad(UaForATC(), "gradUaATC") & U;
85
86 if (extraConvection_ > 0)
87 {
88 // Implicit part added to increase diagonal dominance
89 // Note: Maybe this needs to be multiplied with the ATClimiter ??
90 UaEqn += extraConvection_*fvm::div(-phi, Ua);
91
92 // correct rhs due to implicitly augmenting the adjoint convection
93 ATC_ += extraConvection_*(fvc::grad(UaForATC(), "gradUaATC")().T() & U);
94 }
95
96 // Zero ATC on cells next to given patch types
97 smoothATC();
98
99 // Actual ATC term
100 UaEqn += ATC_.internalField();
101}
102
103
105{
106 tmp<volTensorField> tvolSDTerm
107 (
109 (
111 (
112 "ATCFISensitivityTerm" + type(),
113 mesh_.time().timeName(),
114 mesh_,
117 ),
118 mesh_,
120 )
121 );
122 volTensorField& volSDTerm = tvolSDTerm.ref();
123
124 const volVectorField& U = primalVars_.U();
125 const volVectorField& Ua = adjointVars_.Ua();
126
127 //const volScalarField& mask = getLimiter();
128
129 volSDTerm -=
130 Ua.component(0) * fvc::grad(U.component(0) * U)
131 + Ua.component(1) * fvc::grad(U.component(1) * U)
132 + Ua.component(2) * fvc::grad(U.component(2) * U);
133
134 return tvolSDTerm;
135}
136
137
138// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140} // End namespace Foam
141
142// ************************************************************************* //
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
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 Conservative form of the Momentum equations.
Definition: ATCUaGradU.H:57
virtual tmp< volTensorField > getFISensitivityTerm() const
Get the FI sensitivity derivatives term coming from the ATC.
Definition: ATCUaGradU.C:104
virtual void addATC(fvVectorMatrix &UaEqn)
Add ATC.
Definition: ATCUaGradU.C:65
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
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 Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const volVectorField & Ua() const
Return const reference to velocity.
const volVectorField & UaInst() const
Return const reference to velocity.
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
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
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict