Skip to content

Commit

Permalink
Allowed variance to skip reading data files and return 0 (#533)
Browse files Browse the repository at this point in the history
* Used 0 if variance objFunc can't find data files.

* Fixed the res print for pimple.
  • Loading branch information
friedenhe authored Dec 3, 2023
1 parent eb23bb0 commit 6db1720
Show file tree
Hide file tree
Showing 2 changed files with 89 additions and 59 deletions.
5 changes: 4 additions & 1 deletion dafoam/pyDAFoam.py
Original file line number Diff line number Diff line change
Expand Up @@ -2299,7 +2299,7 @@ def solveAdjointUnsteady(self):

# NOTE: this step is critical because we need to compute the residual for
# self.solverAD once to get the proper oldTime level for unsteady adjoint
self.solverAD.calcPrimalResidualStatistics("print".encode())
self.solverAD.calcPrimalResidualStatistics("calc".encode())

# calc the total number of time instances
# we assume the adjoint is for deltaT to endTime
Expand All @@ -2319,6 +2319,9 @@ def solveAdjointUnsteady(self):
# now we can read the variables
self.readStateVars(endTime, deltaT)

# now we can print the residual for the endTime state
self.solverAD.calcPrimalResidualStatistics("print".encode())

# init dRdWTMF
self.solverAD.initializedRdWTMatrixFree(self.xvVec, self.wVec)

Expand Down
143 changes: 85 additions & 58 deletions src/adjoint/DAObjFunc/DAObjFuncVariance.C
Original file line number Diff line number Diff line change
Expand Up @@ -64,14 +64,17 @@ DAObjFuncVariance::DAObjFuncVariance(

if (varType_ == "scalar")
{
scalar greatVal = varUpperBound_ * 10.0;
volScalarField varData(
IOobject(
varName_ + "Data",
Foam::name(endTime),
mesh_,
IOobject::MUST_READ,
IOobject::READ_IF_PRESENT,
IOobject::NO_WRITE),
mesh_);
mesh_,
dimensionedScalar("dummy", dimensionSet(0, 0, 0, 0, 0, 0, 0), greatVal),
"zeroGradient");

forAll(varData, cellI)
{
Expand All @@ -83,14 +86,17 @@ DAObjFuncVariance::DAObjFuncVariance(
}
else if (varType_ == "vector")
{
scalar greatVal = varUpperBound_ * 10.0;
volVectorField varData(
IOobject(
varName_ + "Data",
Foam::name(endTime),
mesh_,
IOobject::MUST_READ,
IOobject::READ_IF_PRESENT,
IOobject::NO_WRITE),
mesh_);
mesh_,
dimensionedVector("dummy", dimensionSet(0, 0, 0, 0, 0, 0, 0), {greatVal, greatVal, greatVal}),
"zeroGradient");

forAll(varData, cellI)
{
Expand All @@ -114,70 +120,88 @@ DAObjFuncVariance::DAObjFuncVariance(
nRefPointsGlobal_ = nRefPoints_;
reduce(nRefPointsGlobal_, sumOp<label>());

refCellIndex_.setSize(nRefPoints_, -1);
refCellComp_.setSize(nRefPoints_, -1);
refValue_.setSize(nTimeSteps);

// second loop, set refValue
if (varType_ == "scalar")
if (nRefPointsGlobal_ == 0)
{
// if we cant find the varData file or there is no valid values in the varData file
// we just create a zero-size refCellIndex_, and the calcObjFunc will return 0
// without calculating the variance. Plus, we will print out an warning.
refCellIndex_.setSize(0);

Info << endl;
Info << "**************************************************************************** " << endl;
Info << "* WARNING! Can't find data files or can't find valid * " << endl;
Info << "* values in the data file for the variance objFunc! * " << endl;
Info << "**************************************************************************** " << endl;
Info << endl;
}
else
{
for (label n = 0; n < nTimeSteps; n++)

refCellIndex_.setSize(nRefPoints_, -1);
refCellComp_.setSize(nRefPoints_, -1);
refValue_.setSize(nTimeSteps);

// second loop, set refValue
if (varType_ == "scalar")
{
scalar t = (n + 1) * deltaT;
word timeName = Foam::name(t);

refValue_[n].setSize(nRefPoints_);

volScalarField varData(
IOobject(
varName_ + "Data",
timeName,
mesh_,
IOobject::MUST_READ,
IOobject::NO_WRITE),
mesh_);

label pointI = 0;
forAll(varData, cellI)
for (label n = 0; n < nTimeSteps; n++)
{
if (fabs(varData[cellI]) < varUpperBound_)
scalar t = (n + 1) * deltaT;
word timeName = Foam::name(t);

refValue_[n].setSize(nRefPoints_);

volScalarField varData(
IOobject(
varName_ + "Data",
timeName,
mesh_,
IOobject::MUST_READ,
IOobject::NO_WRITE),
mesh_);

label pointI = 0;
forAll(varData, cellI)
{
refValue_[n][pointI] = varData[cellI];
refCellIndex_[pointI] = cellI;
pointI++;
if (fabs(varData[cellI]) < varUpperBound_)
{
refValue_[n][pointI] = varData[cellI];
refCellIndex_[pointI] = cellI;
pointI++;
}
}
}
}
}
else if (varType_ == "vector")
{
for (label n = 0; n < nTimeSteps; n++)
else if (varType_ == "vector")
{
scalar t = (n + 1) * deltaT;
word timeName = Foam::name(t);

refValue_[n].setSize(nRefPoints_);

volVectorField varData(
IOobject(
varName_ + "Data",
timeName,
mesh_,
IOobject::MUST_READ,
IOobject::NO_WRITE),
mesh_);

label pointI = 0;
forAll(varData, cellI)
for (label n = 0; n < nTimeSteps; n++)
{
for (label comp = 0; comp < 3; comp++)
scalar t = (n + 1) * deltaT;
word timeName = Foam::name(t);

refValue_[n].setSize(nRefPoints_);

volVectorField varData(
IOobject(
varName_ + "Data",
timeName,
mesh_,
IOobject::MUST_READ,
IOobject::NO_WRITE),
mesh_);

label pointI = 0;
forAll(varData, cellI)
{
if (fabs(varData[cellI][comp]) < varUpperBound_)
for (label comp = 0; comp < 3; comp++)
{
refValue_[n][pointI] = varData[cellI][comp];
refCellIndex_[pointI] = cellI;
refCellComp_[pointI] = comp;
pointI++;
if (fabs(varData[cellI][comp]) < varUpperBound_)
{
refValue_[n][pointI] = varData[cellI][comp];
refCellIndex_[pointI] = cellI;
refCellComp_[pointI] = comp;
pointI++;
}
}
}
}
Expand Down Expand Up @@ -252,7 +276,10 @@ void DAObjFuncVariance::calcObjFunc(
// need to reduce the sum of force across all processors
reduce(objFuncValue, sumOp<scalar>());

objFuncValue /= nRefPointsGlobal_;
if (nRefPointsGlobal_ != 0)
{
objFuncValue /= nRefPointsGlobal_;
}

return;
}
Expand Down

0 comments on commit 6db1720

Please sign in to comment.