Skip to content

Commit

Permalink
Adding in flag for calculating the heat flux per unit area or for tot…
Browse files Browse the repository at this point in the history
…al surface
  • Loading branch information
ChrisPsenica committed Nov 14, 2024
1 parent bd386fd commit 1649d8f
Show file tree
Hide file tree
Showing 2 changed files with 43 additions and 10 deletions.
50 changes: 40 additions & 10 deletions src/adjoint/DAObjFunc/DAObjFuncWallHeatFlux.C
Original file line number Diff line number Diff line change
Expand Up @@ -74,11 +74,23 @@ DAObjFuncWallHeatFlux::DAObjFuncWallHeatFlux(
"calculated")
#endif
{

// Assign type, this is common for all objectives
objFuncDict_.readEntry<word>("type", objFuncType_);

objFuncDict_.readEntry<scalar>("scale", scale_);

// Heat flux can be calculated by either per unit area or over entire surface. Default value is byUnitArea
if (objFuncDict_.found("scheme"))
{
objFuncDict_.readEntry<word>("scheme", calcMode_);
}
else
{
calcMode_ = "byUnitArea";
}


#ifdef CompressibleFlow

// setup the connectivity for heat flux, this is needed in Foam::DAJacCondFdW
Expand Down Expand Up @@ -162,17 +174,20 @@ void DAObjFuncWallHeatFlux::calcObjFunc(
objFuncValue: the sum of objective, reduced across all processors and scaled by "scale"
*/

// always calculate the area of all the heat flux patches
areaSum_ = 0.0;
forAll(objFuncFaceSources, idxI)
// if calcMode is per unit area then calculate the area of all the heat flux patches
if (calcMode_ == "byUnitArea")
{
const label& objFuncFaceI = objFuncFaceSources[idxI];
label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
const label patchI = daIndex_.bFacePatchI[bFaceI];
const label faceI = daIndex_.bFaceFaceI[bFaceI];
areaSum_ += mesh_.magSf().boundaryField()[patchI][faceI];
areaSum_ = 0.0;
forAll(objFuncFaceSources, idxI)
{
const label& objFuncFaceI = objFuncFaceSources[idxI];
label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
const label patchI = daIndex_.bFacePatchI[bFaceI];
const label faceI = daIndex_.bFaceFaceI[bFaceI];
areaSum_ += mesh_.magSf().boundaryField()[patchI][faceI];
}
reduce(areaSum_, sumOp<scalar>());
}
reduce(areaSum_, sumOp<scalar>());

// initialize faceValues to zero
forAll(objFuncFaceValues, idxI)
Expand Down Expand Up @@ -238,7 +253,22 @@ void DAObjFuncWallHeatFlux::calcObjFunc(
const label faceI = daIndex_.bFaceFaceI[bFaceI];

scalar area = mesh_.magSf().boundaryField()[patchI][faceI];
objFuncFaceValues[idxI] = scale_ * wallHeatFluxBf[patchI][faceI] * area / areaSum_;

if (calcMode_ == "byUnitArea")
{
objFuncFaceValues[idxI] = scale_ * wallHeatFluxBf[patchI][faceI] * area / areaSum_;
}
else if (calcMode_ == "total")
{
objFuncFaceValues[idxI] = scale_ * wallHeatFluxBf[patchI][faceI] * area;
}
else
{
FatalErrorIn(" ") << "mode for "
<< objFuncName_ << " " << objFuncPart_ << " not valid!"
<< "Options: byUnitArea (default value), total."
<< abort(FatalError);
}

objFuncValue += objFuncFaceValues[idxI];
}
Expand Down
3 changes: 3 additions & 0 deletions src/adjoint/DAObjFunc/DAObjFuncWallHeatFlux.H
Original file line number Diff line number Diff line change
Expand Up @@ -50,6 +50,9 @@ protected:
/// the area of all heat flux patches
scalar areaSum_ = -9999.0;

/// if calculating flux per unit area or total, which mode to use
word calcMode_;

public:
TypeName("wallHeatFlux");
// Constructors
Expand Down

0 comments on commit 1649d8f

Please sign in to comment.