DAInputStateVar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAInputStateVar.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAInputStateVar, 0);
16 addToRunTimeSelectionTable(DAInput, DAInputStateVar, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const word inputName,
21  const word inputType,
22  fvMesh& mesh,
23  const DAOption& daOption,
24  const DAModel& daModel,
25  const DAIndex& daIndex)
26  : DAInput(
27  inputName,
28  inputType,
29  mesh,
30  daOption,
31  daModel,
32  daIndex)
33 {
34 }
35 
36 void DAInputStateVar::run(const scalarList& input)
37 {
38  /*
39  Description:
40  Assign the input array to OF's state variables
41  */
42 
43 #ifndef CODI_ADR
44  Info << "DAInputStateVar. " << endl;
45  Info << "Setting state variables. " << endl;
46 #endif
47 
48  forAll(stateInfo_["volVectorStates"], idxI)
49  {
50  const word stateName = stateInfo_["volVectorStates"][idxI];
51  volVectorField& state = const_cast<volVectorField&>(
52  mesh_.thisDb().lookupObject<volVectorField>(stateName));
53 
54  forAll(mesh_.cells(), cellI)
55  {
56  for (label i = 0; i < 3; i++)
57  {
58  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, cellI, i);
59  state[cellI][i] = input[localIdx];
60  }
61  }
62  state.correctBoundaryConditions();
63  }
64 
65  forAll(stateInfo_["volScalarStates"], idxI)
66  {
67  const word stateName = stateInfo_["volScalarStates"][idxI];
68  volScalarField& state = const_cast<volScalarField&>(
69  mesh_.thisDb().lookupObject<volScalarField>(stateName));
70 
71  forAll(mesh_.cells(), cellI)
72  {
73  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, cellI);
74  state[cellI] = input[localIdx];
75  }
76  state.correctBoundaryConditions();
77  }
78 
79  forAll(stateInfo_["modelStates"], idxI)
80  {
81  const word stateName = stateInfo_["modelStates"][idxI];
82  volScalarField& state = const_cast<volScalarField&>(
83  mesh_.thisDb().lookupObject<volScalarField>(stateName));
84 
85  forAll(mesh_.cells(), cellI)
86  {
87  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, cellI);
88  state[cellI] = input[localIdx];
89  }
90  state.correctBoundaryConditions();
91  }
92 
93  forAll(stateInfo_["surfaceScalarStates"], idxI)
94  {
95  const word stateName = stateInfo_["surfaceScalarStates"][idxI];
96  surfaceScalarField& state = const_cast<surfaceScalarField&>(
97  mesh_.thisDb().lookupObject<surfaceScalarField>(stateName));
98 
99  forAll(mesh_.faces(), faceI)
100  {
101  label localIdx = daIndex_.getLocalAdjointStateIndex(stateName, faceI);
102 
103  if (faceI < daIndex_.nLocalInternalFaces)
104  {
105  state[faceI] = input[localIdx];
106  }
107  else
108  {
109  label relIdx = faceI - daIndex_.nLocalInternalFaces;
110  label patchIdx = daIndex_.bFacePatchI[relIdx];
111  label faceIdx = daIndex_.bFaceFaceI[relIdx];
112  state.boundaryFieldRef()[patchIdx][faceIdx] = input[localIdx];
113  }
114  }
115  }
116 }
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 } // End namespace Foam
121 
122 // ************************************************************************* //
Foam::DAInput::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAInput.H:57
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAOption
Definition: DAOption.H:29
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DAInputStateVar::run
virtual void run(const scalarList &input)
Definition: DAInputStateVar.C:36
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAInput::stateInfo_
HashTable< wordList > stateInfo_
the stateInfo_ list from DAStateInfo object
Definition: DAInput.H:60
Foam::DAModel
Definition: DAModel.H:57
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
DAInputStateVar.H
Foam::DAInputStateVar::DAInputStateVar
DAInputStateVar(const word inputName, const word inputType, fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAInputStateVar.C:19
Foam::DAInput
Definition: DAInput.H:30
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
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::DAInput::mesh_
fvMesh & mesh_
fvMesh
Definition: DAInput.H:48
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116