View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0002383 | OpenFOAM | Patch | public | 2016-12-09 15:20 | 2016-12-09 17:00 |
Reporter | MattijsJ | Assigned To | henry | ||
Priority | high | Severity | major | Reproducibility | always |
Status | resolved | Resolution | fixed | ||
Platform | GNU/Linux | OS | OpenSuSE | OS Version | 13.2 |
Product Version | dev | ||||
Fixed in Version | dev | ||||
Summary | 0002383: GAMG solver does not run with processorAgglomerator | ||||
Description | The processorAgglomerator does not work anymore since using PBiCGStab p { solver GAMG; tolerance 1e-08; relTol 0.1; smoother GaussSeidel; nCellsInCoarsestLevel 20; processorAgglomerator masterCoarsest; cacheAgglomeration true; } | ||||
Additional Information | Attached patch. Was using the worldComm inside the coarsest level. | ||||
Tags | No tags attached. | ||||
|
PBiCGStab.C (7,391 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2016 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 "PBiCGStab.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(PBiCGStab, 0); lduMatrix::solver::addasymMatrixConstructorToTable<PBiCGStab> addPBiCGStabAsymMatrixConstructorToTable_; } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::PBiCGStab::PBiCGStab ( const word& fieldName, const lduMatrix& matrix, const FieldField<Field, scalar>& interfaceBouCoeffs, const FieldField<Field, scalar>& interfaceIntCoeffs, const lduInterfaceFieldPtrsList& interfaces, const dictionary& solverControls ) : lduMatrix::solver ( fieldName, matrix, interfaceBouCoeffs, interfaceIntCoeffs, interfaces, solverControls ) {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // Foam::solverPerformance Foam::PBiCGStab::solve ( scalarField& psi, const scalarField& source, const direction cmpt ) const { // --- Setup class containing solver performance data solverPerformance solverPerf ( lduMatrix::preconditioner::getName(controlDict_) + typeName, fieldName_ ); const label nCells = psi.size(); scalar* __restrict__ psiPtr = psi.begin(); scalarField pA(nCells); scalar* __restrict__ pAPtr = pA.begin(); scalarField yA(nCells); scalar* __restrict__ yAPtr = yA.begin(); // --- Calculate A.psi matrix_.Amul(yA, psi, interfaceBouCoeffs_, interfaces_, cmpt); // --- Calculate initial residual field scalarField rA(source - yA); scalar* __restrict__ rAPtr = rA.begin(); // --- Calculate normalisation factor const scalar normFactor = this->normFactor(psi, source, yA, pA); if (lduMatrix::debug >= 2) { Info<< " Normalisation factor = " << normFactor << endl; } // --- Calculate normalised residual norm solverPerf.initialResidual() = gSumMag(rA, matrix().mesh().comm()) /normFactor; solverPerf.finalResidual() = solverPerf.initialResidual(); // --- Check convergence, solve if not converged if ( minIter_ > 0 || !solverPerf.checkConvergence(tolerance_, relTol_) ) { scalarField AyA(nCells); scalar* __restrict__ AyAPtr = AyA.begin(); scalarField sA(nCells); scalar* __restrict__ sAPtr = sA.begin(); scalarField zA(nCells); scalar* __restrict__ zAPtr = zA.begin(); scalarField tA(nCells); scalar* __restrict__ tAPtr = tA.begin(); // --- Store initial residual const scalarField rA0(rA); // --- Initial values not used scalar rA0rA = 0; scalar alpha = 0; scalar omega = 0; // --- Select and construct the preconditioner autoPtr<lduMatrix::preconditioner> preconPtr = lduMatrix::preconditioner::New ( *this, controlDict_ ); // --- Solver iteration do { // --- Store previous rA0rA const scalar rA0rAold = rA0rA; rA0rA = gSumProd(rA0, rA, matrix().mesh().comm()); // --- Test for singularity if (solverPerf.checkSingularity(mag(rA0rA))) { break; } // --- Update pA if (solverPerf.nIterations() == 0) { for (label cell=0; cell<nCells; cell++) { pAPtr[cell] = rAPtr[cell]; } } else { // --- Test for singularity if (solverPerf.checkSingularity(mag(omega))) { break; } const scalar beta = (rA0rA/rA0rAold)*(alpha/omega); for (label cell=0; cell<nCells; cell++) { pAPtr[cell] = rAPtr[cell] + beta*(pAPtr[cell] - omega*AyAPtr[cell]); } } // --- Precondition pA preconPtr->precondition(yA, pA, cmpt); // --- Calculate AyA matrix_.Amul(AyA, yA, interfaceBouCoeffs_, interfaces_, cmpt); const scalar rA0AyA = gSumProd(rA0, AyA, matrix().mesh().comm()); alpha = rA0rA/rA0AyA; // --- Calculate sA for (label cell=0; cell<nCells; cell++) { sAPtr[cell] = rAPtr[cell] - alpha*AyAPtr[cell]; } // --- Test sA for convergence solverPerf.finalResidual() = gSumMag(sA, matrix().mesh().comm())/normFactor; if (solverPerf.checkConvergence(tolerance_, relTol_)) { for (label cell=0; cell<nCells; cell++) { psiPtr[cell] += alpha*yAPtr[cell]; } solverPerf.nIterations()++; return solverPerf; } // --- Precondition sA preconPtr->precondition(zA, sA, cmpt); // --- Calculate tA matrix_.Amul(tA, zA, interfaceBouCoeffs_, interfaces_, cmpt); const scalar tAtA = gSumSqr(tA, matrix().mesh().comm()); // --- Calculate omega from tA and sA // (cheaper than using zA with preconditioned tA) omega = gSumProd(tA, sA, matrix().mesh().comm())/tAtA; // --- Update solution and residual for (label cell=0; cell<nCells; cell++) { psiPtr[cell] += alpha*yAPtr[cell] + omega*zAPtr[cell]; rAPtr[cell] = sAPtr[cell] - omega*tAPtr[cell]; } solverPerf.finalResidual() = gSumMag(rA, matrix().mesh().comm()) /normFactor; } while ( ( solverPerf.nIterations()++ < maxIter_ && !solverPerf.checkConvergence(tolerance_, relTol_) ) || solverPerf.nIterations() < minIter_ ); } return solverPerf; } // ************************************************************************* // |
|
Resolved by commit 9eb6736ef68ed7c3ffd7e68cace51e2ab4f5b066 |
Date Modified | Username | Field | Change |
---|---|---|---|
2016-12-09 15:20 | MattijsJ | New Issue | |
2016-12-09 15:20 | MattijsJ | File Added: PBiCGStab.C | |
2016-12-09 17:00 | henry | Assigned To | => henry |
2016-12-09 17:00 | henry | Status | new => resolved |
2016-12-09 17:00 | henry | Resolution | open => fixed |
2016-12-09 17:00 | henry | Fixed in Version | => dev |
2016-12-09 17:00 | henry | Note Added: 0007450 |