Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Reduce the IO load for pimple #526

Merged
merged 5 commits into from
Nov 22, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
4 changes: 4 additions & 0 deletions dafoam/pyDAFoam.py
Original file line number Diff line number Diff line change
Expand Up @@ -467,6 +467,7 @@ def __init__(self):
"objFuncEndTime": -1.0,
"PCMatPrecomputeInterval": 100,
"PCMatUpdateInterval": 1,
"reduceIO": True,
}

## At which iteration should we start the averaging of objective functions.
Expand Down Expand Up @@ -3126,18 +3127,21 @@ def runColoring(self):
from .pyColoringIncompressible import pyColoringIncompressible

solverArg = "ColoringIncompressible -python " + self.parallelFlag

solver = pyColoringIncompressible(solverArg.encode(), self.options)
elif solverName in self.solverRegistry["Compressible"]:

from .pyColoringCompressible import pyColoringCompressible

solverArg = "ColoringCompressible -python " + self.parallelFlag

solver = pyColoringCompressible(solverArg.encode(), self.options)
elif solverName in self.solverRegistry["Solid"]:

from .pyColoringSolid import pyColoringSolid

solverArg = "ColoringSolid -python " + self.parallelFlag

solver = pyColoringSolid(solverArg.encode(), self.options)
else:
raise Error("pyDAFoam: %s not registered! Check _solverRegistry(self)." % solverName)
Expand Down
18 changes: 14 additions & 4 deletions src/adjoint/DASolver/DAPimpleFoam/DAPimpleFoam.C
Original file line number Diff line number Diff line change
Expand Up @@ -151,6 +151,14 @@ label DAPimpleFoam::solvePrimal(
// reset the unsteady obj func to zeros
this->initUnsteadyObjFuncs();

// we need to reduce the number of files written to the disk to minimize the file IO load
label reduceIO = daOptionPtr_->getAllOptions().subDict("unsteadyAdjoint").getLabel("reduceIO");
if (reduceIO)
{
// set all states and vars to NO_WRITE
this->disableStateAutoWrite();
}

// main loop
while (runTime.run())
{
Expand Down Expand Up @@ -208,11 +216,13 @@ label DAPimpleFoam::solvePrimal(
<< nl << endl;
}

runTime.write();

if (mode_ == "hybrid")
if (reduceIO)
{
this->writeAdjStates();
}
else
{
this->saveTimeInstanceFieldHybrid(timeInstanceI);
runTime.write();
}
}

Expand Down
4 changes: 2 additions & 2 deletions src/adjoint/DASolver/DAPimpleFoam/DAPimpleFoam.H
Original file line number Diff line number Diff line change
Expand Up @@ -95,8 +95,8 @@ protected:
/// pressure reference value
scalar pRefValue_ = 0.0;

/// whether the hybrid adjoint or time accurate adjoint is active
word mode_ = "None";
/// the primal IO interval
label IOInterval_ = 1;

public:
TypeName("DAPimpleFoam");
Expand Down
144 changes: 144 additions & 0 deletions src/adjoint/DASolver/DASolver.C
Original file line number Diff line number Diff line change
Expand Up @@ -8668,6 +8668,150 @@ void DASolver::calcPCMatWithFvMatrix(Mat PCMat)
#endif
}

void DASolver::disableStateAutoWrite()
{
/*
Description:
set state variables to NO_WRITE to reduce the file IO load
*/

// volVector states
forAll(stateInfo_["volVectorStates"], idxI)
{
const word stateName = stateInfo_["volVectorStates"][idxI];
volVectorField& state =
const_cast<volVectorField&>(meshPtr_->thisDb().lookupObject<volVectorField>(stateName));

state.writeOpt() = IOobject::NO_WRITE;
}

// volScalar states
forAll(stateInfo_["volScalarStates"], idxI)
{
const word stateName = stateInfo_["volScalarStates"][idxI];
volScalarField& state =
const_cast<volScalarField&>(meshPtr_->thisDb().lookupObject<volScalarField>(stateName));

state.writeOpt() = IOobject::NO_WRITE;
}

// model states
forAll(stateInfo_["modelStates"], idxI)
{
const word stateName = stateInfo_["modelStates"][idxI];
volScalarField& state =
const_cast<volScalarField&>(meshPtr_->thisDb().lookupObject<volScalarField>(stateName));

state.writeOpt() = IOobject::NO_WRITE;
}

// surfaceScalar states
forAll(stateInfo_["surfaceScalarStates"], idxI)
{
const word stateName = stateInfo_["surfaceScalarStates"][idxI];
surfaceScalarField& state =
const_cast<surfaceScalarField&>(meshPtr_->thisDb().lookupObject<surfaceScalarField>(stateName));

state.writeOpt() = IOobject::NO_WRITE;
}

// some potential extra variables
if (meshPtr_->thisDb().foundObject<volScalarField>("nut"))
{
volScalarField& nut =
const_cast<volScalarField&>(meshPtr_->thisDb().lookupObject<volScalarField>("nut"));
nut.writeOpt() = IOobject::NO_WRITE;
}

if (meshPtr_->thisDb().foundObject<volScalarField>("betaFI"))
{
volScalarField& betaFI =
const_cast<volScalarField&>(meshPtr_->thisDb().lookupObject<volScalarField>("betaFI"));
betaFI.writeOpt() = IOobject::NO_WRITE;
}

if (meshPtr_->thisDb().foundObject<volVectorField>("fvSource"))
{
volVectorField& fvSource =
const_cast<volVectorField&>(meshPtr_->thisDb().lookupObject<volVectorField>("fvSource"));
fvSource.writeOpt() = IOobject::NO_WRITE;
}
}

void DASolver::writeAdjStates()
{
/*
Description:
Write only the adjoint states
*/

// volVector states
forAll(stateInfo_["volVectorStates"], idxI)
{
const word stateName = stateInfo_["volVectorStates"][idxI];
volVectorField& state =
const_cast<volVectorField&>(meshPtr_->thisDb().lookupObject<volVectorField>(stateName));

state.write();
}

// volScalar states
forAll(stateInfo_["volScalarStates"], idxI)
{
const word stateName = stateInfo_["volScalarStates"][idxI];
volScalarField& state =
const_cast<volScalarField&>(meshPtr_->thisDb().lookupObject<volScalarField>(stateName));

state.write();
}

// model states
forAll(stateInfo_["modelStates"], idxI)
{
const word stateName = stateInfo_["modelStates"][idxI];
volScalarField& state =
const_cast<volScalarField&>(meshPtr_->thisDb().lookupObject<volScalarField>(stateName));

state.write();
}

// surfaceScalar states
forAll(stateInfo_["surfaceScalarStates"], idxI)
{
const word stateName = stateInfo_["surfaceScalarStates"][idxI];
surfaceScalarField& state =
const_cast<surfaceScalarField&>(meshPtr_->thisDb().lookupObject<surfaceScalarField>(stateName));

state.write();
}

scalar endTime = runTimePtr_->endTime().value();
scalar deltaT = runTimePtr_->deltaT().value();
label nInstances = round(endTime / deltaT);

// write these extra suppressed variables for the last time step
if (runTimePtr_->timeIndex() == nInstances)
{
if (meshPtr_->thisDb().foundObject<volScalarField>("nut"))
{
const volScalarField& nut = meshPtr_->thisDb().lookupObject<volScalarField>("nut");
nut.write();
}

if (meshPtr_->thisDb().foundObject<volScalarField>("betaFI"))
{
const volScalarField& betaFI = meshPtr_->thisDb().lookupObject<volScalarField>("betaFI");
betaFI.write();
}

if (meshPtr_->thisDb().foundObject<volVectorField>("fvSource"))
{
const volVectorField& fvSource = meshPtr_->thisDb().lookupObject<volVectorField>("fvSource");
fvSource.write();
}
}
}

void DASolver::readStateVars(
scalar timeVal,
label oldTimeLevel)
Expand Down
6 changes: 6 additions & 0 deletions src/adjoint/DASolver/DASolver.H
Original file line number Diff line number Diff line change
Expand Up @@ -1134,6 +1134,12 @@ public:
DAUtility::pyCalcBetaJacVecProd = jacVecProd;
}

/// set state variables to NO_WRITE
void disableStateAutoWrite();

/// write state variables that are NO_WRITE to disk
void writeAdjStates();

#ifdef CODI_AD_REVERSE

/// global tape for reverse-mode AD
Expand Down
Loading