DAOutputResidual.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAOutputResidual.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAOutputResidual, 0);
16 addToRunTimeSelectionTable(DAOutput, DAOutputResidual, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const word outputName,
21  const word outputType,
22  fvMesh& mesh,
23  const DAOption& daOption,
24  DAModel& daModel,
25  const DAIndex& daIndex,
26  DAResidual& daResidual,
27  UPtrList<DAFunction>& daFunctionList)
28  : DAOutput(
29  outputName,
30  outputType,
31  mesh,
32  daOption,
33  daModel,
34  daIndex,
35  daResidual,
36  daFunctionList)
37 {
38 }
39 
40 void DAOutputResidual::run(scalarList& output)
41 {
42  /*
43  Description:
44  Compute OF's residual variables and then assign them to the output array
45  */
46 
47  label isPC = 0;
48 
49  // calculate the residual
50  dictionary resOptions;
51  resOptions.set("isPC", isPC);
52  daResidual_.calcResiduals(resOptions);
53  daModel_.calcResiduals(resOptions);
54 
55  // assign the calculated residuals to output
56  forAll(stateInfo_["volVectorStates"], idxI)
57  {
58  const word stateName = stateInfo_["volVectorStates"][idxI];
59  const word stateResName = stateName + "Res";
60  volVectorField& stateRes = const_cast<volVectorField&>(
61  mesh_.thisDb().lookupObject<volVectorField>(stateResName));
62 
63  forAll(mesh_.cells(), cellI)
64  {
65  for (label i = 0; i < 3; i++)
66  {
67  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, cellI, i);
68  output[localIdx] = stateRes[cellI][i];
69  }
70  }
71  }
72 
73  forAll(stateInfo_["volScalarStates"], idxI)
74  {
75  const word stateName = stateInfo_["volScalarStates"][idxI];
76  const word stateResName = stateName + "Res";
77  volScalarField& stateRes = const_cast<volScalarField&>(
78  mesh_.thisDb().lookupObject<volScalarField>(stateResName));
79 
80  forAll(mesh_.cells(), cellI)
81  {
82  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, cellI);
83  output[localIdx] = stateRes[cellI];
84  }
85  }
86 
87  forAll(stateInfo_["modelStates"], idxI)
88  {
89  const word stateName = stateInfo_["modelStates"][idxI];
90  const word stateResName = stateName + "Res";
91  volScalarField& stateRes = const_cast<volScalarField&>(
92  mesh_.thisDb().lookupObject<volScalarField>(stateResName));
93 
94  forAll(mesh_.cells(), cellI)
95  {
96  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, cellI);
97  output[localIdx] = stateRes[cellI];
98  }
99  }
100 
101  forAll(stateInfo_["surfaceScalarStates"], idxI)
102  {
103  const word stateName = stateInfo_["surfaceScalarStates"][idxI];
104  const word stateResName = stateName + "Res";
105  surfaceScalarField& stateRes = const_cast<surfaceScalarField&>(
106  mesh_.thisDb().lookupObject<surfaceScalarField>(stateResName));
107 
108  forAll(mesh_.faces(), faceI)
109  {
110  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, faceI);
111 
112  if (faceI < daIndex_.nLocalInternalFaces)
113  {
114  output[localIdx] = stateRes[faceI];
115  }
116  else
117  {
118  label relIdx = faceI - daIndex_.nLocalInternalFaces;
119  label patchIdx = daIndex_.bFacePatchI[relIdx];
120  label faceIdx = daIndex_.bFaceFaceI[relIdx];
121  output[localIdx] = stateRes.boundaryFieldRef()[patchIdx][faceIdx];
122  }
123  }
124  }
125 }
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 } // End namespace Foam
130 
131 // ************************************************************************* //
Foam::DAModel::calcResiduals
void calcResiduals(const dictionary &options)
calculate the residuals for model state variables
Definition: DAModel.C:190
Foam::DAOutput::daResidual_
DAResidual & daResidual_
Definition: DAOutput.H:61
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAOption
Definition: DAOption.H:29
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DAOutput::stateInfo_
HashTable< wordList > stateInfo_
the stateInfo_ list from DAStateInfo object
Definition: DAOutput.H:66
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAOutput::daModel_
DAModel & daModel_
DAIndex object.
Definition: DAOutput.H:56
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAOutput
Definition: DAOutput.H:32
Foam::DAModel
Definition: DAModel.H:57
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAResidual
Definition: DAResidual.H:36
DAOutputResidual.H
Foam::DAIndex::getLocalAdjointStateIndex
label getLocalAdjointStateIndex(const word stateName, const label idxI, const label comp=-1) const
get local adjoint index for a given state name, cell/face indxI and its component (optional,...
Definition: DAIndex.C:516
Foam::DAOutputResidual::DAOutputResidual
DAOutputResidual(const word outputName, const word outputType, fvMesh &mesh, const DAOption &daOption, DAModel &daModel, const DAIndex &daIndex, DAResidual &daResidual, UPtrList< DAFunction > &daFunctionList)
Definition: DAOutputResidual.C:19
Foam::DAOutput::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAOutput.H:59
Foam::DAOutputResidual::run
virtual void run(scalarList &output)
Definition: DAOutputResidual.C:40
Foam::DAResidual::calcResiduals
virtual void calcResiduals(const dictionary &options)=0
compute residuals
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116
Foam::DAOutput::mesh_
fvMesh & mesh_
fvMesh
Definition: DAOutput.H:50