62 this->
normFactor(psi, tsource(), Apsi, finestCorrection);
130 (scratch1.
size() ? scratch1 : Apsi),
131 (scratch2.
size() ? scratch2 : finestCorrection),
140 finestResidual = tsource();
141 finestResidual -= Apsi;
174void Foam::GAMGSolver::Vcycle
193 const label coarsestLevel = matrixLevels_.size() - 1;
196 agglomeration_.restrictField(coarseSources[0], finestResidual, 0,
true);
198 if (nPreSweeps_ && ((log_ >= 2) || (
debug >= 2)))
200 Pout<<
"Pre-smoothing scaling factors: ";
205 for (label leveli = 0; leveli < coarsestLevel; leveli++)
207 if (coarseSources.
set(leveli + 1))
213 coarseCorrFields[leveli] = 0.0;
215 smoothers[leveli + 1].scalarSmooth
217 coarseCorrFields[leveli],
218 coarseSources[leveli],
222 nPreSweeps_ + preSweepsLevelMultiplier_*leveli,
229 if (scaleCorrection_ && leveli < coarsestLevel - 1)
234 coarseCorrFields[leveli].size()
239 coarseCorrFields[leveli],
244 matrixLevels_[leveli],
245 interfaceLevelsBouCoeffs_[leveli],
246 interfaceLevels_[leveli],
247 coarseSources[leveli],
254 matrixLevels_[leveli].residual
256 coarseSources[leveli],
257 coarseCorrFields[leveli],
260 coarseSources[leveli]
262 interfaceLevelsBouCoeffs_[leveli],
263 interfaceLevels_[leveli],
269 agglomeration_.restrictField
271 coarseSources[leveli + 1],
272 coarseSources[leveli],
279 if (nPreSweeps_ && ((log_ >= 2) || (
debug >= 2)))
286 if (coarseCorrFields.
set(coarsestLevel))
290 coarseCorrFields[coarsestLevel],
291 coarseSources[coarsestLevel]
295 if ((log_ >= 2) || (
debug >= 2))
297 Pout<<
"Post-smoothing scaling factors: ";
308 for (label leveli = coarsestLevel - 1; leveli >= 0; leveli--)
310 if (coarseCorrFields.
set(leveli))
318 coarseCorrFields[leveli].size()
325 preSmoothedCoarseCorrField = coarseCorrFields[leveli];
330 const auto& cf = agglomeration_.prolongField
332 coarseCorrFields[leveli],
335 coarseCorrFields.
set(leveli + 1)
336 ? coarseCorrFields[leveli + 1]
347 coarseCorrFields[leveli].size()
354 if (interpolateCorrection_)
360 coarseCorrFields[leveli],
362 matrixLevels_[leveli],
363 interfaceLevelsBouCoeffs_[leveli],
364 interfaceLevels_[leveli],
365 agglomeration_.restrictAddressing(leveli + 1),
376 && (interpolateCorrection_ || leveli < coarsestLevel - 1)
381 coarseCorrFields[leveli],
383 matrixLevels_[leveli],
384 interfaceLevelsBouCoeffs_[leveli],
385 interfaceLevels_[leveli],
386 coarseSources[leveli],
395 coarseCorrFields[leveli] += preSmoothedCoarseCorrField;
398 smoothers[leveli + 1].scalarSmooth
400 coarseCorrFields[leveli],
401 coarseSources[leveli],
405 nPostSweeps_ + postSweepsLevelMultiplier_*leveli,
413 agglomeration_.prolongField
421 if (interpolateCorrection_)
430 agglomeration_.restrictAddressing(0),
436 if (scaleCorrection_)
453 psi[i] += finestCorrection[i];
466void Foam::GAMGSolver::initVcycle
475 label maxSize = matrix_.diag().size();
477 coarseCorrFields.setSize(matrixLevels_.size());
478 coarseSources.setSize(matrixLevels_.size());
479 smoothers.setSize(matrixLevels_.size() + 1);
496 forAll(matrixLevels_, leveli)
498 if (agglomeration_.nCells(leveli) >= 0)
500 label nCoarseCells = agglomeration_.nCells(leveli);
505 if (matrixLevels_.set(leveli))
507 const lduMatrix& mat = matrixLevels_[leveli];
509 label nCoarseCells = mat.
diag().
size();
511 maxSize =
max(maxSize, nCoarseCells);
521 matrixLevels_[leveli],
522 interfaceLevelsBouCoeffs_[leveli],
523 interfaceLevelsIntCoeffs_[leveli],
524 interfaceLevels_[leveli],
531 if (maxSize > matrix_.diag().size())
534 scratch1.setSize(maxSize);
535 scratch2.setSize(maxSize);
540Foam::dictionary Foam::GAMGSolver::PCGsolverDict
554Foam::dictionary Foam::GAMGSolver::PBiCGStabSolverDict
568void Foam::GAMGSolver::solveCoarsestLevel
574 const label coarsestLevel = matrixLevels_.size() - 1;
576 const label coarseComm = matrixLevels_[coarsestLevel].mesh().comm();
578 if (directSolveCoarsest_)
582 coarsestLUMatrixPtr_->solve
689 coarsestCorrField = 0;
692 coarsestSolverPtr_->scalarSolve
699 if ((log_ >= 2) ||
debug)
701 coarseSolverPerf.print(
Info.masterStream(coarseComm));
A const Field/List wrapper with possible data conversion.
SubField< solveScalar > subField
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve.
Input from string buffer, using a ISstream. Always UNCOMPRESSED.
A non-const Field/List wrapper with possible data conversion.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
static autoPtr< smoother > New(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new smoother.
const lduMatrix & matrix_
label maxIter_
Maximum number of iterations in the solver.
scalar tolerance_
Final convergence tolerance.
lduInterfaceFieldPtrsList interfaces_
const FieldField< Field, scalar > & interfaceBouCoeffs_
label minIter_
Minimum number of iterations in the solver.
const lduMatrix & matrix() const noexcept
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField, const lduMatrix::normTypes normType) const
Return the matrix norm using the specified norm method.
int log_
Verbosity level for solver output statements.
scalar relTol_
Convergence tolerance relative to the initial.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
const scalarField & diag() const
void setResidualField(const scalarField &residual, const word &fieldName, const bool initial) const
Set the residual field using an IOField on the object registry if it exists.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
const volScalarField & psi
Namespace for handling debugging switches.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< solveScalar > solveScalarField
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
#define forAll(list, i)
Loop across all elements in list.