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-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 defineTypeNameAndDebug(ATCUaGradU, 0);
42 (
43  ATCModel,
44  ATCUaGradU,
46 );
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 ATCUaGradU::ATCUaGradU
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();
68  const volVectorField& Ua = adjointVars_.UaInst();
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
90 
91  // Correct rhs due to implicitly augmenting the adjoint convection
92  ATC_ += extraConvection_*(fvc::grad(UaForATC(), "gradUaATC")().T() & U);
93  }
94 
95  // Zero ATC on cells next to given patch types
97 
98  // Actual ATC term
99  UaEqn += ATC_.internalField();
100 }
101 
102 
104 {
105  tmp<volTensorField> tvolSDTerm
106  (
107  new volTensorField
108  (
109  IOobject
110  (
111  "ATCFISensitivityTerm" + type(),
112  mesh_.time().timeName(),
113  mesh_,
116  ),
117  mesh_,
119  )
120  );
121  volTensorField& volSDTerm = tvolSDTerm.ref();
122 
123  const volVectorField& U = primalVars_.U();
124  const volVectorField& Ua = adjointVars_.Ua();
125 
126  //const volScalarField& mask = getLimiter();
127 
128  volSDTerm -=
129  Ua.component(0) * fvc::grad(U.component(0) * U)
130  + Ua.component(1) * fvc::grad(U.component(1) * U)
131  + Ua.component(2) * fvc::grad(U.component(2) * U);
132 
133  return tvolSDTerm;
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138 
139 } // End namespace Foam
140 
141 // ************************************************************************* //
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
volScalarField ATClimiter_
Definition: ATCModel.H:91
void smoothATC()
Limit ATC field using ATClimiter_.
Definition: ATCModel.C:44
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
virtual tmp< volTensorField > getFISensitivityTerm() const
Get the FI sensitivity derivatives term coming from the ATC.
Definition: ATCUaGradU.C:96
const incompressibleVars & primalVars_
Definition: ATCModel.H:80
Macros for easy insertion into run-time selection tables.
Class including all adjoint fields for incompressible flows.
const incompressibleAdjointVars & adjointVars_
Definition: ATCModel.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const surfaceScalarField & phi() const
Return const reference to volume flux.
Base class for solution control classes.
dynamicFvMesh & mesh
const fvMesh & mesh_
Definition: ATCModel.H:79
const volVectorField & Ua() const
Return const reference to velocity.
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition: ATCModel.H:56
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
volVectorField ATC_
Definition: ATCModel.H:92
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
const scalar extraConvection_
Definition: ATCModel.H:85
const volVectorField & U() const
Return const reference to velocity.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const bool reconstructGradients_
Definition: ATCModel.H:88
const volVectorField & UaInst() const
Return const reference to velocity.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual void addATC(fvVectorMatrix &UaEqn)
Add ATC.
Definition: ATCUaGradU.C:58
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127