fvMeshDistributors::loadBalancer: Prototype general CPU load balancer

used in conjunction with the new loadBalancing option in constant/chemistryProperties:

    loadBalancing   on;

which enables per-cell CPU time caching used by the loadBalancer to redistribute
the mesh.  Currently this option is only provided for chemistry integration but
the implementation is general and in future options will be provided to balance
other local cell loads, in particular Lagrangian particles.

The loadBalancer in enabled by specifying a distributor entry in
constant/dynamicMeshDict, e.g.

distributor
{
    type            loadBalancer;

    libs            ("libfvMeshDistributors.so");

    multiConstraint true;

    // How often to redistribute
    redistributionInterval  10;

    // Maximum fractional cell distribution imbalance
    // before rebalancing
    maxImbalance    0.1;
}

with which the mesh is checked for more than 10% load-imbalance every 10
time-steps and redistributed using a multi-constraint method, i.e. separate CPU
load weights are provided for each of the loads, currently that is the chemistry
integration load and the CPU time taken for the rest of the simulation,
transport equations solution etc.

The fvMeshDistributors::loadBalancer uses the distributor specified in
system/decomposeParDict to redistribute the mesh based on the cell CPU loads,
e.g. to use the Zoltan RCB method specify:

distributor     zoltan;
libs            ("libzoltanDecomp.so");

zoltanCoeffs
{
    lb_method   rcb;
}

Unfortunately only a few available redistribution methods support
multi-constraints: Zoltan::RCB, MeTiS, parMeTiS and xtraPuLP, of these only
Zoltan::RCB is currently available in OpenFOAM.  Load-balancing is possible
without using a multi-constraint method (i.e. using any of the other
decomposition methods provided with OpenFOAM and Zoltan) by summing the various
CPU loads which is selected by setting:

    multiConstraint false;

but the load-balancing is likely to be a lot less effective with this option.

Due to the licencing issues with parMeTiS interfacing to xtraPuLP might be the
best option for further work on load-balancing in OpenFOAM, or MeTiS could be
used in parallel by first agglomerating the distribution graph on the master
processor and redistributing the result; this pseudo-parallel option is already
provided for the Scotch method.
This commit is contained in:
Henry Weller
2022-01-17 11:31:12 +00:00
parent 3c353761ed
commit 472ce5ace6
29 changed files with 533 additions and 63 deletions

View File

@ -1,3 +1,4 @@
distributor/fvMeshDistributorsDistributor.C
loadBalancer/fvMeshDistributorsLoadBalancer.C
LIB = $(FOAM_LIBBIN)/libfvMeshDistributors

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -60,16 +60,13 @@ void Foam::fvMeshDistributors::distributor::readDict()
}
void Foam::fvMeshDistributors::distributor::distribute()
void Foam::fvMeshDistributors::distributor::distribute
(
const labelList& distribution
)
{
fvMesh& mesh = this->mesh();
// Create new decomposition distribution
labelList distribution
(
distributor_->decompose(mesh, mesh.cellCentres())
);
// Mesh distribution engine
fvMeshDistribute distributor(mesh);
@ -114,22 +111,26 @@ Foam::fvMeshDistributors::distributor::~distributor()
bool Foam::fvMeshDistributors::distributor::update()
{
const fvMesh& mesh = this->mesh();
bool redistributed = false;
if
(
Pstream::nProcs() > 1
&& mesh().time().timeIndex() > 1
&& timeIndex_ != mesh().time().timeIndex()
&& mesh().time().timeIndex() % redistributionInterval_ == 0
&& mesh.time().timeIndex() > 1
&& timeIndex_ != mesh.time().timeIndex()
&& mesh.time().timeIndex() % redistributionInterval_ == 0
)
{
timeIndex_ = mesh().time().timeIndex();
timeIndex_ = mesh.time().timeIndex();
const scalar idealNCells =
mesh().globalData().nTotalCells()/Pstream::nProcs();
mesh.globalData().nTotalCells()/Pstream::nProcs();
const scalar imbalance = returnReduce
(
mag(1 - mesh().nCells()/idealNCells),
mag(1 - mesh.nCells()/idealNCells),
maxOp<scalar>()
);
@ -137,19 +138,19 @@ bool Foam::fvMeshDistributors::distributor::update()
{
Info<< "Redistributing mesh with imbalance " << imbalance << endl;
distribute();
// Create new decomposition distribution
const labelList distribution
(
distributor_->decompose(mesh, mesh.cellCentres())
);
return true;
}
else
{
return false;
distribute(distribution);
redistributed = true;
}
}
else
{
return false;
}
return redistributed;
}

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -73,7 +73,9 @@ class distributor
:
public fvMeshDistributor
{
// Private Member Data
protected:
// Protected Member Data
//- Cache the decomposer/distributor
autoPtr<decompositionMethod> distributor_;
@ -89,13 +91,13 @@ class distributor
label timeIndex_;
// Private Member Functions
// Protected Member Functions
//- Read the projection parameters from dictionary
void readDict();
//- Distribute the mesh and mesh data
void distribute();
void distribute(const labelList& distribution);
public:

View File

@ -0,0 +1,181 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "fvMeshDistributorsLoadBalancer.H"
#include "decompositionMethod.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace fvMeshDistributors
{
defineTypeNameAndDebug(loadBalancer, 0);
addToRunTimeSelectionTable
(
fvMeshDistributor,
loadBalancer,
fvMesh
);
}
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
void Foam::fvMeshDistributors::loadBalancer::readDict()
{
distributor::readDict();
const dictionary& distributorDict(dict());
multiConstraint_ =
distributorDict.lookupOrDefault<Switch>("multiConstraint", true);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::fvMeshDistributors::loadBalancer::loadBalancer(fvMesh& mesh)
:
distributor(mesh)
{
readDict();
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::fvMeshDistributors::loadBalancer::~loadBalancer()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::fvMeshDistributors::loadBalancer::update()
{
const fvMesh& mesh = this->mesh();
bool redistributed = false;
if
(
Pstream::nProcs() > 1
&& mesh.time().timeIndex() > 1
&& timeIndex_ != mesh.time().timeIndex()
)
{
timeIndex_ = mesh.time().timeIndex();
const scalar timeStepCpuTime = cpuTime_.cpuTimeIncrement();
// Chemistry CPU load per cell
volScalarField::Internal& chemistryCpuTimeReg =
mesh.lookupObjectRef<volScalarField::Internal>
(
"chemistryCpuTime"
);
const scalarField& chemistryCpuTime = chemistryCpuTimeReg.field();
if (mesh.time().timeIndex() % redistributionInterval_ == 0)
{
timeIndex_ = mesh.time().timeIndex();
const scalar sumCpuLoad(sum(chemistryCpuTime));
const scalar cellCFDCpuTime = returnReduce
(
(timeStepCpuTime - sumCpuLoad)/mesh.nCells(),
minOp<scalar>()
);
// Total CPU time for this processor
const scalar processorCpuTime =
mesh.nCells()*cellCFDCpuTime + sumCpuLoad;
// Average processor CPU time
const scalar averageProcessorCpuTime =
returnReduce(processorCpuTime, sumOp<scalar>())
/Pstream::nProcs();
Pout<< "imbalance "
<< " " << sumCpuLoad
<< " " << mesh.nCells()*cellCFDCpuTime
<< " " << processorCpuTime
<< " " << averageProcessorCpuTime << endl;
const scalar imbalance = returnReduce
(
mag(1 - processorCpuTime/averageProcessorCpuTime),
maxOp<scalar>()
);
scalarField weights;
if (multiConstraint_)
{
const int nWeights = 2;
weights.setSize(nWeights*mesh.nCells());
forAll(chemistryCpuTime, i)
{
weights[nWeights*i] = cellCFDCpuTime;
weights[nWeights*i + 1] = chemistryCpuTime[i];
}
}
else
{
weights = chemistryCpuTime + cellCFDCpuTime;
}
if (imbalance > maxImbalance_)
{
Info<< "Redistributing mesh with imbalance "
<< imbalance << endl;
// Create new decomposition distribution
const labelList distribution
(
distributor_->decompose(mesh, mesh.cellCentres(), weights)
);
distribute(distribution);
redistributed = true;
}
}
chemistryCpuTimeReg.checkOut();
}
return redistributed;
}
// ************************************************************************* //

View File

@ -0,0 +1,126 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::fvMeshDistributors::loadBalancer
Description
Dynamic mesh redistribution using the distributor specified in
decomposeParDict
Usage
Example of single field based refinement in all cells:
\verbatim
distributor
{
type loadBalancer;
libs ("libfvMeshDistributors.so");
// How often to redistribute
redistributionInterval 10;
// Maximum fractional cell distribution imbalance
// before rebalancing
maxImbalance 0.1;
}
\endverbatim
SourceFiles
fvMeshDistributorsloadBalancer.C
\*---------------------------------------------------------------------------*/
#ifndef fvMeshDistributorsLoadBalancer_H
#define fvMeshDistributorsLoadBalancer_H
#include "fvMeshDistributorsDistributor.H"
#include "cpuTime.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace fvMeshDistributors
{
/*---------------------------------------------------------------------------*\
Class loadBalancer Declaration
\*---------------------------------------------------------------------------*/
class loadBalancer
:
public distributor
{
// Private Member Data
//- CPU time consumed during the time-step
cpuTime cpuTime_;
//- Enable multi-constraint load-balancing in which separate weights
// are provided to the distrubutor for each of the CPU loads.
// When disabled the CPU loads are summed and a single weight per cell
// is provided to the distrubutor.
// Defaults to true.
Switch multiConstraint_;
// Private Member Functions
//- Read the projection parameters from dictionary
void readDict();
public:
//- Runtime type information
TypeName("loadBalancer");
// Constructors
//- Construct from fvMesh
explicit loadBalancer(fvMesh& mesh);
//- Destructor
virtual ~loadBalancer();
// Member Functions
//- Distribute the
virtual bool update();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace fvMeshDistributors
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -92,13 +92,23 @@ static void get_vertex_list
const Foam::pointField& points = vertexData.first();
const Foam::scalarField& weights = vertexData.second();
if (wgt_dim != weights.size()/points.size())
{
*ierr = ZOLTAN_FATAL;
return;
}
Foam::globalIndex globalMap(points.size());
for (Foam::label i=0; i<points.size(); i++)
{
localIDs[i] = i;
globalIDs[i] = globalMap.toGlobal(i);
obj_wgts[i] = weights[i];
for(int j=0; j<wgt_dim; j++)
{
obj_wgts[wgt_dim*i + j] = weights[wgt_dim*i + j];
}
}
*ierr = ZOLTAN_OK;
@ -214,7 +224,6 @@ static void get_edge_list
(
(nGID != 1)
|| (nLID != 1)
|| (wgt_dim != 0)
)
{
*ierr = ZOLTAN_FATAL;
@ -260,9 +269,11 @@ Foam::label Foam::zoltanDecomp::decompose
struct Zoltan_Struct *zz = Zoltan_Create(PstreamGlobals::MPI_COMM_FOAM);
const int nWeights = pWeights.size()/points.size();
// Set internal parameters
Zoltan_Set_Param(zz, "return_lists", "export");
Zoltan_Set_Param(zz, "obj_weight_dim", "1");
Zoltan_Set_Param(zz, "obj_weight_dim", name(nWeights).c_str());
Zoltan_Set_Param(zz, "edge_weight_dim", "0");
// General default paramaters

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2016-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2016-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -26,7 +26,7 @@ License
#include "chemistryModel.H"
#include "UniformField.H"
#include "localEulerDdtScheme.H"
#include "clockTime.H"
#include "cpuTime.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
@ -39,6 +39,7 @@ Foam::chemistryModel<ThermoType>::chemistryModel
basicChemistryModel(thermo),
ODESystem(),
log_(this->lookupOrDefault("log", false)),
loadBalancing_(this->lookupOrDefault("loadBalancing", false)),
jacobianType_
(
this->found("jacobian")
@ -709,14 +710,43 @@ Foam::scalar Foam::chemistryModel<ThermoType>::solve
const DeltaTType& deltaT
)
{
if (loadBalancing_)
{
if
(
!this->mesh().objectRegistry::template
foundObject<volScalarField::Internal>("chemistryCpuTime")
)
{
regIOobject::store
(
volScalarField::Internal::New
(
"chemistryCpuTime",
this->mesh(),
dimensionedScalar(dimTime, 0)
).ptr()
);
}
}
volScalarField::Internal& chemistryCpuTime =
loadBalancing_
? this->mesh().objectRegistry::template
lookupObjectRef<volScalarField::Internal>
(
"chemistryCpuTime"
)
: const_cast<volScalarField::Internal&>(volScalarField::Internal::null());
tabulation_.reset();
const basicSpecieMixture& composition = this->thermo().composition();
// CPU time analysis
const clockTime clockTime_ = clockTime();
clockTime_.timeIncrement();
scalar solveChemistryCpuTime_ = 0;
cpuTime cpuTime_;
cpuTime solveCpuTime_;
scalar totalSolveCpuTime_ = 0;
basicChemistryModel::correct();
@ -810,8 +840,8 @@ Foam::scalar Foam::chemistryModel<ThermoType>::solve
if (log_)
{
// Reset the time
clockTime_.timeIncrement();
// Reset the solve time
solveCpuTime_.cpuTimeIncrement();
}
// Calculate the chemical source terms
@ -845,7 +875,7 @@ Foam::scalar Foam::chemistryModel<ThermoType>::solve
if (log_)
{
solveChemistryCpuTime_ += clockTime_.timeIncrement();
totalSolveCpuTime_ += solveCpuTime_.cpuTimeIncrement();
}
// If tabulation is used, we add the information computed here to
@ -880,13 +910,18 @@ Foam::scalar Foam::chemistryModel<ThermoType>::solve
{
RR_[i][celli] = (Y_[i]*rho - Y0[i]*rho0)/deltaT[celli];
}
if (loadBalancing_)
{
chemistryCpuTime[celli] += cpuTime_.cpuTimeIncrement();
}
}
if (log_)
{
cpuSolveFile_()
<< this->time().userTimeValue()
<< " " << solveChemistryCpuTime_ << endl;
<< " " << totalSolveCpuTime_ << endl;
}
mechRed_.update();

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2016-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2016-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -130,6 +130,9 @@ class chemistryModel
//- Switch to select performance logging
Switch log_;
//- Switch to enable loadBalancing performance logging
Switch loadBalancing_;
//- Type of the Jacobian to be calculated
const jacobianType jacobianType_;

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -42,7 +42,6 @@ Foam::chemistryReductionMethod<ThermoType>::chemistryReductionMethod
activeSpecies_(chemistry.nSpecie(), true),
log_(false),
tolerance_(NaN),
clockTime_(clockTime()),
sumnActiveSpecies_(0),
sumn_(0),
reduceMechCpuTime_(0)
@ -64,7 +63,6 @@ Foam::chemistryReductionMethod<ThermoType>::chemistryReductionMethod
activeSpecies_(chemistry.nSpecie(), false),
log_(coeffsDict_.lookupOrDefault<Switch>("log", false)),
tolerance_(coeffsDict_.lookupOrDefault<scalar>("tolerance", 1e-4)),
clockTime_(clockTime()),
sumnActiveSpecies_(0),
sumn_(0),
reduceMechCpuTime_(0)
@ -91,7 +89,7 @@ void Foam::chemistryReductionMethod<ThermoType>::initReduceMechanism()
{
if (log_)
{
clockTime_.timeIncrement();
cpuTime_.cpuTimeIncrement();
}
}
@ -162,7 +160,7 @@ void Foam::chemistryReductionMethod<ThermoType>::endReduceMechanism
{
sumnActiveSpecies_ += nActiveSpecies_;
sumn_++;
reduceMechCpuTime_ += clockTime_.timeIncrement();
reduceMechCpuTime_ += cpuTime_.cpuTimeIncrement();
}
}

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -41,7 +41,7 @@ SourceFiles
#include "IOdictionary.H"
#include "DynamicField.H"
#include "Switch.H"
#include "clockTime.H"
#include "cpuTime.H"
#include "OFstream.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -101,8 +101,8 @@ private:
//- Tolerance for the mechanism reduction algorithm
scalar tolerance_;
//- Clock time for logging
const clockTime clockTime_;
//- CPU time for logging
cpuTime cpuTime_;
//- ...
int64_t sumnActiveSpecies_;

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2016-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2016-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -83,7 +83,6 @@ Foam::chemistryTabulationMethods::ISAT<ThermoType>::ISAT
addNewLeafCpuTime_(0),
growCpuTime_(0),
searchISATCpuTime_(0),
clockTime_(clockTime()),
tabulationResults_
(
IOobject
@ -411,7 +410,7 @@ bool Foam::chemistryTabulationMethods::ISAT<ThermoType>::retrieve
{
if (log_)
{
clockTime_.timeIncrement();
cpuTime_.cpuTimeIncrement();
}
bool retrieved(false);
@ -480,7 +479,7 @@ bool Foam::chemistryTabulationMethods::ISAT<ThermoType>::retrieve
if (log_)
{
searchISATCpuTime_ += clockTime_.timeIncrement();
searchISATCpuTime_ += cpuTime_.cpuTimeIncrement();
}
return retrieved;
@ -499,7 +498,7 @@ Foam::label Foam::chemistryTabulationMethods::ISAT<ThermoType>::add
{
if (log_)
{
clockTime_.timeIncrement();
cpuTime_.cpuTimeIncrement();
}
label growthOrAddFlag = 1;
@ -518,7 +517,7 @@ Foam::label Foam::chemistryTabulationMethods::ISAT<ThermoType>::add
if (log_)
{
growCpuTime_ += clockTime_.timeIncrement();
growCpuTime_ += cpuTime_.cpuTimeIncrement();
}
// the structure of the tree is not modified, return false
@ -607,7 +606,7 @@ Foam::label Foam::chemistryTabulationMethods::ISAT<ThermoType>::add
if (log_)
{
addNewLeafCpuTime_ += clockTime_.timeIncrement();
addNewLeafCpuTime_ += cpuTime_.cpuTimeIncrement();
}
return growthOrAddFlag;

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2016-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2016-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -119,7 +119,7 @@ class ISAT
scalar growCpuTime_;
scalar searchISATCpuTime_;
const clockTime clockTime_;
cpuTime cpuTime_;
autoPtr<OFstream> nRetrievedFile_;
autoPtr<OFstream> nGrowthFile_;

View File

@ -20,6 +20,8 @@ internalField uniform 0.0;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 0;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 0;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 1;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 0;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 2000;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform (0 0 0);
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 0.0;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 0;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type fixedValue;

View File

@ -20,6 +20,8 @@ internalField uniform 1e5;
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
fuel
{
type zeroGradient;

View File

@ -0,0 +1,11 @@
#!/bin/sh
cd ${0%/*} || exit 1 # Run from this directory
# Source tutorial run functions
. $WM_PROJECT_DIR/bin/tools/RunFunctions
runApplication blockMesh
runApplication $(getApplication)
#------------------------------------------------------------------------------

View File

@ -0,0 +1,12 @@
#!/bin/sh
cd ${0%/*} || exit 1 # Run from this directory
# Source tutorial run functions
. $WM_PROJECT_DIR/bin/tools/RunFunctions
runApplication blockMesh
runApplication decomposePar
runParallel $(getApplication)
#------------------------------------------------------------------------------

View File

@ -21,6 +21,8 @@ chemistryType
chemistry on;
loadBalancing on;
initialChemicalTimeStep 1e-07;
odeCoeffs

View File

@ -0,0 +1,27 @@
/*--------------------------------*- C++ -*----------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Version: dev
\\/ M anipulation |
\*---------------------------------------------------------------------------*/
FoamFile
{
format ascii;
class dictionary;
object dynamicMeshDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
distributor
{
type loadBalancer;
libs ("libfvMeshDistributors.so");
multiConstraint true;
redistributionInterval 10;
}
// ************************************************************************* //

View File

@ -0,0 +1,27 @@
/*--------------------------------*- C++ -*----------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Version: dev
\\/ M anipulation |
\*---------------------------------------------------------------------------*/
FoamFile
{
format ascii;
class dictionary;
object dynamicMeshDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
distributor
{
type loadBalancer;
libs ("libfvMeshDistributors.so");
multiConstraint true;
redistributionInterval 10;
}
// ************************************************************************* //

View File

@ -68,6 +68,12 @@ boundary
(0 3 2 1)
);
}
internal
{
type internal;
faces ();
}
);

View File

@ -14,14 +14,22 @@ FoamFile
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
numberOfSubdomains 4;
numberOfSubdomains 12;
method hierarchical;
decomposer hierarchical;
distributor zoltan;
libs ("libzoltanDecomp.so");
hierarchicalCoeffs
{
n (2 2 1);
n (6 2 1);
order xyz;
}
zoltanCoeffs
{
lb_method rcb;
}
// ************************************************************************* //