DAResidual.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  Description:
7  Compute the residual for a given solver
8 
9 \*---------------------------------------------------------------------------*/
10 
11 #ifndef DAResidual_H
12 #define DAResidual_H
13 
14 #include "runTimeSelectionTables.H"
15 #include "fvOptions.H"
16 #include "surfaceFields.H"
17 #include "DAOption.H"
18 #include "DAModel.H"
19 #include "DAMacroFunctions.H"
20 #include "DAUtility.H"
21 #include "DAIndex.H"
22 #include "DAField.H"
23 #include "DAFvSource.H"
24 #include "IOMRFZoneListDF.H"
25 #include "constrainHbyA.H"
26 
27 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
28 
29 namespace Foam
30 {
31 
32 /*---------------------------------------------------------------------------*\
33  Class DAResidual Declaration
34 \*---------------------------------------------------------------------------*/
35 
37  : public regIOobject
38 {
39 
40 private:
42  DAResidual(const DAResidual&);
43 
45  void operator=(const DAResidual&);
46 
47 protected:
49  const fvMesh& mesh_;
50 
53 
55  const DAModel& daModel_;
56 
58  const DAIndex& daIndex_;
59 
62 
63 public:
65  TypeName("DAResidual");
66 
67  // Declare run-time constructor selection table
69  autoPtr,
70  DAResidual,
71  dictionary,
72  (const word modelType,
73  const fvMesh& mesh,
74  const DAOption& daOption,
75  const DAModel& daModel,
76  const DAIndex& daIndex),
77  (modelType, mesh, daOption, daModel, daIndex));
78 
79  // Constructors
80 
81  //- Construct from components
82  DAResidual(
83  const word modelType,
84  const fvMesh& mesh,
85  const DAOption& daOption,
86  const DAModel& daModel,
87  const DAIndex& daIndex);
88 
89  // Selectors
90 
91  //- Return a reference to the selected model
92  static autoPtr<DAResidual> New(
93  const word modelType,
94  const fvMesh& mesh,
95  const DAOption& daOption,
96  const DAModel& daModel,
97  const DAIndex& daIndex);
98 
99  //- Destructor
100  virtual ~DAResidual()
101  {
102  }
103 
104  // Members
105 
107  virtual void clear() = 0;
108 
110  virtual void calcResiduals(const dictionary& options) = 0;
111 
113  virtual void updateIntermediateVariables() = 0;
114 
116  virtual void correctBoundaryConditions() = 0;
117 
119  void masterFunction(
120  const dictionary& options,
121  const Vec xvVec,
122  const Vec wVec,
123  Vec resVec);
124 
126  virtual void calcPCMatWithFvMatrix(Mat PCMat);
127 
129  bool writeData(Ostream& os) const
130  {
131  return true;
132  }
133 };
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 } // End namespace Foam
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 #endif
142 
143 // ************************************************************************* //
Foam::DAResidual::correctBoundaryConditions
virtual void correctBoundaryConditions()=0
update the boundary condition for all the states in the selected solver
DAOption.H
Foam::DAResidual::writeData
bool writeData(Ostream &os) const
virtual function for regIOobject
Definition: DAResidual.H:129
Foam::DAResidual::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAResidual.H:49
DAIndex.H
IOMRFZoneListDF.H
Foam::DAResidual::calcPCMatWithFvMatrix
virtual void calcPCMatWithFvMatrix(Mat PCMat)
calculating the adjoint preconditioner matrix using fvMatrix
Definition: DAResidual.C:169
Foam::DAOption
Definition: DAOption.H:29
DAUtility.H
Foam::DAResidual::updateIntermediateVariables
virtual void updateIntermediateVariables()=0
update any intermdiate variables that are dependent on state variables and are used in calcResiduals
DAFvSource.H
DAMacroFunctions.H
Foam::DAResidual::daIndex_
const DAIndex & daIndex_
DAIndex.
Definition: DAResidual.H:58
Foam::DAResidual::daField_
DAField daField_
DAField object.
Definition: DAResidual.H:61
Foam::DAResidual::clear
virtual void clear()=0
clear the members
Foam::DAResidual::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAResidual.H:52
DAModel.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAResidual::New
static autoPtr< DAResidual > New(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidual.C:47
Foam::DAField
Definition: DAField.H:36
Foam::DAModel
Definition: DAModel.H:57
Foam::DAResidual::~DAResidual
virtual ~DAResidual()
Definition: DAResidual.H:100
Foam
Definition: checkGeometry.C:32
Foam::DAResidual::TypeName
TypeName("DAResidual")
Runtime type information.
Foam::DAResidual
Definition: DAResidual.H:36
Foam::DAResidual::daModel_
const DAModel & daModel_
DAModel object.
Definition: DAResidual.H:55
Foam::DAResidual::masterFunction
void masterFunction(const dictionary &options, const Vec xvVec, const Vec wVec, Vec resVec)
the master function that compute the residual vector given the state and point vectors
Definition: DAResidual.C:90
Foam::DAResidual::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, DAResidual, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex),(modelType, mesh, daOption, daModel, daIndex))
DAField.H
Foam::DAResidual::calcResiduals
virtual void calcResiduals(const dictionary &options)=0
compute residuals