calcForcePerS.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  Description:
7  Calculating force per surface area. This will be used to plot the
8  spanwise force distribution
9 
10 \*---------------------------------------------------------------------------*/
11 
12 #include "fvCFD.H"
13 #include "argList.H"
14 #include "autoPtr.H"
15 #include "Time.H"
16 #include "timeSelector.H"
17 #include "TimePaths.H"
18 #include "ListOps.H"
19 #include "fvMesh.H"
20 #include "OFstream.H"
21 #include "simpleControl.H"
22 #include "fvOptions.H"
23 #include "fluidThermo.H"
24 
25 using namespace Foam;
26 
27 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
28 
29 int main(int argc, char* argv[])
30 {
31  Info << "Computing forces...." << endl;
32 
33  argList::addOption(
34  "patchNames",
35  "'(wing)'",
36  "List of patch names to compute");
37 
38  argList::addOption(
39  "time",
40  "1000",
41  "Time instance to compute, if not provided runs all times");
42 
43 #include "setRootCase.H"
44 #include "createTime.H"
45 
46  if (!args.optionFound("time"))
47  {
48  Info << "Time not set, running all times." << endl;
49  }
50 
51  // Create the processor databases
52  autoPtr<TimePaths> timePaths;
53  timePaths = autoPtr<TimePaths>::New(args.rootPath(), args.caseName());
54 
55  const instantList timeDirs(timeSelector::select(timePaths->times(), args));
56 
57  List<wordRe> patchNames;
58  if (args.optionFound("patchNames"))
59  {
60  patchNames = wordReList(args.optionLookup("patchNames")());
61  }
62  else
63  {
64  Info << "patchNames not set! Exit." << endl;
65  return 1;
66  }
67 
68  forAll(timeDirs, iTime)
69  {
70  runTime.setTime(timeDirs[iTime].value(), 0);
71 
72 #include "createMesh.H"
73 #include "createFields.H"
74 
75  // Initialize surface field for face-centered forces
76  volVectorField forcePerS(
77  IOobject(
78  "forcePerS",
79  runTime.timeName(),
80  mesh,
81  IOobject::NO_READ,
82  IOobject::NO_WRITE),
83  mesh,
84  dimensionedVector("forcePerS", dimensionSet(1, -1, -2, 0, 0, 0, 0), vector::zero),
85  fixedValueFvPatchScalarField::typeName);
86 
87  // this code is pulled from:
88  // src/functionObjects/forcces/forces.C
89  // modified slightly
90  vector forces(vector::zero);
91 
92  const surfaceVectorField::Boundary& Sfb = mesh.Sf().boundaryField();
93  const surfaceScalarField::Boundary& magSfb = mesh.magSf().boundaryField();
94 
95  volSymmTensorField devRhoReff(
96  IOobject(
97  IOobject::groupName("devRhoReff", U.group()),
98  mesh.time().timeName(),
99  mesh,
100  IOobject::NO_READ,
101  IOobject::NO_WRITE),
102  (-rho * nuEff) * dev(twoSymm(fvc::grad(U))));
103 
104  const volSymmTensorField::Boundary& devRhoReffb = devRhoReff.boundaryField();
105 
106  vector totalForce(vector::zero);
107  forAll(patchNames, cI)
108  {
109  // get the patch id label
110  label patchI = mesh.boundaryMesh().findPatchID(patchNames[cI]);
111  if (patchI < 0)
112  {
113  Info << "ERROR: Patch name " << patchNames[cI] << " not found in constant/polyMesh/boundary! Exit!" << endl;
114  return 1;
115  }
116  // create a shorter handle for the boundary patch
117  const fvPatch& patch = mesh.boundary()[patchI];
118  // normal force
119  vectorField fN(Sfb[patchI] * p.boundaryField()[patchI]);
120  // tangential force
121  vectorField fT(Sfb[patchI] & devRhoReffb[patchI]);
122  // sum them up
123  forAll(patch, faceI)
124  {
125  // compute forces
126  forces.x() = fN[faceI].x() + fT[faceI].x();
127  forces.y() = fN[faceI].y() + fT[faceI].y();
128  forces.z() = fN[faceI].z() + fT[faceI].z();
129 
130  // project force direction
131  forcePerS.boundaryFieldRef()[patchI][faceI] = forces / magSfb[patchI][faceI];
132 
133  totalForce.x() += forces.x();
134  totalForce.y() += forces.y();
135  totalForce.z() += forces.z();
136  }
137  }
138  forcePerS.write();
139 
140  Info << "Total force: " << totalForce << endl;
141 
142  Info << "Computing forces.... Completed!" << endl;
143  }
144  return 0;
145 }
146 
147 // ************************************************************************* //
U
volVectorField & U
Definition: createRefsPimpleDyM.H:7
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
p
volScalarField & p
Definition: createRefsPimpleDyM.H:6
main
int main(int argc, char *argv[])
Definition: calcForcePerS.C:29
Foam
Definition: checkGeometry.C:32
rho
volScalarField & rho
Definition: createRefsRhoPimple.H:8
argc
int argc
Definition: setArgs.H:18
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
args
Foam::argList & args
Definition: setRootCasePython.H:42
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
createFields.H
argv
char ** argv
Definition: setArgs.H:26