Loading...
Searching...
No Matches
GAMGSolverSolve.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2021,2023 OpenCFD Ltd.
10 Copyright (C) 2023 Huawei (Yu Ankun)
11 Copyright (C) 2023 OpenCFD Ltd.
12-------------------------------------------------------------------------------
13License
14 This file is part of OpenFOAM.
15
16 OpenFOAM is free software: you can redistribute it and/or modify it
17 under the terms of the GNU General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24 for more details.
25
26 You should have received a copy of the GNU General Public License
27 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
31#include "GAMGSolver.H"
32#include "SubField.H"
33#include "PrecisionAdaptor.H"
34
35// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36
38(
39 scalarField& psi_s,
40 const scalarField& source,
41 const direction cmpt
42) const
43{
45 solveScalarField& psi = tpsi.ref();
46
48
49 // Setup class containing solver performance data
51
52 // Calculate A.psi used to calculate the initial residual
53 solveScalarField Apsi(psi.size());
54 matrix_.Amul(Apsi, psi, interfaceBouCoeffs_, interfaces_, cmpt);
55
56 // Create the storage for the finestCorrection which may be used as a
57 // temporary in normFactor
58 solveScalarField finestCorrection(psi.size());
59
60 // Calculate normalisation factor
61 solveScalar normFactor =
62 this->normFactor(psi, tsource(), Apsi, finestCorrection);
63
64 if ((log_ >= 2) || (debug >= 2))
65 {
66 Pout<< " Normalisation factor = " << normFactor << endl;
67 }
68
69 // Calculate initial finest-grid residual field
70 solveScalarField finestResidual(tsource() - Apsi);
71
73 (
76 true
77 );
78
79 // Calculate normalised residual for convergence test
80 solverPerf.initialResidual() = gSumMag
81 (
82 finestResidual,
83 matrix().mesh().comm()
84 )/normFactor;
85 solverPerf.finalResidual() = solverPerf.initialResidual();
86
87
88 // Check convergence, solve if not converged
89 if
90 (
91 minIter_ > 0
92 || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
93 )
94 {
95 // Create coarse grid correction fields
96 PtrList<solveScalarField> coarseCorrFields;
97
98 // Create coarse grid sources
99 PtrList<solveScalarField> coarseSources;
100
101 // Create the smoothers for all levels
103
104 // Scratch fields if processor-agglomerated coarse level meshes
105 // are bigger than original. Usually not needed
106 solveScalarField scratch1;
107 solveScalarField scratch2;
108
109 // Initialise the above data structures
110 initVcycle
111 (
112 coarseCorrFields,
113 coarseSources,
114 smoothers,
115 scratch1,
116 scratch2
117 );
118
119 do
120 {
121 Vcycle
122 (
123 smoothers,
124 psi,
125 source,
126 Apsi,
127 finestCorrection,
128 finestResidual,
129
130 (scratch1.size() ? scratch1 : Apsi),
131 (scratch2.size() ? scratch2 : finestCorrection),
132
133 coarseCorrFields,
134 coarseSources,
135 cmpt
136 );
137
138 // Calculate finest level residual field
139 matrix_.Amul(Apsi, psi, interfaceBouCoeffs_, interfaces_, cmpt);
140 finestResidual = tsource();
141 finestResidual -= Apsi;
142
143 solverPerf.finalResidual() = gSumMag
144 (
145 finestResidual,
146 matrix().mesh().comm()
147 )/normFactor;
148
149 if ((log_ >= 2) || (debug >= 2))
150 {
151 solverPerf.print(Info.masterStream(matrix().mesh().comm()));
152 }
153 } while
154 (
155 (
156 ++solverPerf.nIterations() < maxIter_
157 && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
158 )
159 || solverPerf.nIterations() < minIter_
160 );
161 }
162
164 (
167 false
168 );
169
170 return solverPerf;
171}
172
173
174void Foam::GAMGSolver::Vcycle
175(
176 const PtrList<lduMatrix::smoother>& smoothers,
178 const scalarField& source,
179 solveScalarField& Apsi,
180 solveScalarField& finestCorrection,
181 solveScalarField& finestResidual,
182
183 solveScalarField& scratch1,
184 solveScalarField& scratch2,
185
186 PtrList<solveScalarField>& coarseCorrFields,
187 PtrList<solveScalarField>& coarseSources,
188 const direction cmpt
189) const
190{
191 //debug = 2;
192
193 const label coarsestLevel = matrixLevels_.size() - 1;
194
195 // Restrict finest grid residual for the next level up.
196 agglomeration_.restrictField(coarseSources[0], finestResidual, 0, true);
197
198 if (nPreSweeps_ && ((log_ >= 2) || (debug >= 2)))
199 {
200 Pout<< "Pre-smoothing scaling factors: ";
201 }
202
203
204 // Residual restriction (going to coarser levels)
205 for (label leveli = 0; leveli < coarsestLevel; leveli++)
206 {
207 if (coarseSources.set(leveli + 1))
208 {
209 // If the optional pre-smoothing sweeps are selected
210 // smooth the coarse-grid field for the restricted source
211 if (nPreSweeps_)
212 {
213 coarseCorrFields[leveli] = 0.0;
214
215 smoothers[leveli + 1].scalarSmooth
216 (
217 coarseCorrFields[leveli],
218 coarseSources[leveli], //coarseSource,
219 cmpt,
220 min
221 (
222 nPreSweeps_ + preSweepsLevelMultiplier_*leveli,
223 maxPreSweeps_
224 )
225 );
226
227 // Scale coarse-grid correction field
228 // but not on the coarsest level because it evaluates to 1
229 if (scaleCorrection_ && leveli < coarsestLevel - 1)
230 {
232 (
233 scratch1,
234 coarseCorrFields[leveli].size()
235 );
236
237 scale
238 (
239 coarseCorrFields[leveli],
240 const_cast<solveScalarField&>
241 (
242 static_cast<const solveScalarField&>(ACf)
243 ),
244 matrixLevels_[leveli],
245 interfaceLevelsBouCoeffs_[leveli],
246 interfaceLevels_[leveli],
247 coarseSources[leveli],
248 cmpt
249 );
250 }
251
252 // Correct the residual with the new solution
253 // residual can be used by fusing Amul with b-Amul
254 matrixLevels_[leveli].residual
255 (
256 coarseSources[leveli],
257 coarseCorrFields[leveli],
259 (
260 coarseSources[leveli]
261 )(),
262 interfaceLevelsBouCoeffs_[leveli],
263 interfaceLevels_[leveli],
264 cmpt
265 );
266 }
267
268 // Residual is equal to source
269 agglomeration_.restrictField
270 (
271 coarseSources[leveli + 1],
272 coarseSources[leveli],
273 leveli + 1,
274 true
275 );
276 }
277 }
278
279 if (nPreSweeps_ && ((log_ >= 2) || (debug >= 2)))
280 {
281 Pout<< endl;
282 }
283
284
285 // Solve Coarsest level with either an iterative or direct solver
286 if (coarseCorrFields.set(coarsestLevel))
287 {
288 solveCoarsestLevel
289 (
290 coarseCorrFields[coarsestLevel],
291 coarseSources[coarsestLevel]
292 );
293 }
294
295 if ((log_ >= 2) || (debug >= 2))
296 {
297 Pout<< "Post-smoothing scaling factors: ";
298 }
299
300 // Smoothing and prolongation of the coarse correction fields
301 // (going to finer levels)
302
303 solveScalarField dummyField(0);
304
305 // Work storage for prolongation
306 solveScalarField work;
307
308 for (label leveli = coarsestLevel - 1; leveli >= 0; leveli--)
309 {
310 if (coarseCorrFields.set(leveli))
311 {
312 // Create a field for the pre-smoothed correction field
313 // as a sub-field of the finestCorrection which is not
314 // currently being used
315 solveScalarField::subField preSmoothedCoarseCorrField
316 (
317 scratch2,
318 coarseCorrFields[leveli].size()
319 );
320
321 // Only store the preSmoothedCoarseCorrField if pre-smoothing is
322 // used
323 if (nPreSweeps_)
324 {
325 preSmoothedCoarseCorrField = coarseCorrFields[leveli];
326 }
327
328
329 // Prolong correction to leveli
330 const auto& cf = agglomeration_.prolongField
331 (
332 coarseCorrFields[leveli], // current level
333 work,
334 (
335 coarseCorrFields.set(leveli + 1)
336 ? coarseCorrFields[leveli + 1]
337 : dummyField // dummy value
338 ),
339 leveli + 1
340 );
341
342
343 // Create A.psi for this coarse level as a sub-field of Apsi
345 (
346 scratch1,
347 coarseCorrFields[leveli].size()
348 );
349 auto& ACfRef = const_cast<solveScalarField&>
350 (
351 static_cast<const solveScalarField&>(ACf)
352 );
353
354 if (interpolateCorrection_)
355 {
356 // Normal operation : have both coarse level and fine
357 // level. No processor agglomeration
359 (
360 coarseCorrFields[leveli],
361 ACfRef,
362 matrixLevels_[leveli],
363 interfaceLevelsBouCoeffs_[leveli],
364 interfaceLevels_[leveli],
365 agglomeration_.restrictAddressing(leveli + 1),
366 cf,
367 cmpt
368 );
369 }
370
371 // Scale coarse-grid correction field
372 // but not on the coarsest level because it evaluates to 1
373 if
374 (
375 scaleCorrection_
376 && (interpolateCorrection_ || leveli < coarsestLevel - 1)
377 )
378 {
379 scale
380 (
381 coarseCorrFields[leveli],
382 ACfRef,
383 matrixLevels_[leveli],
384 interfaceLevelsBouCoeffs_[leveli],
385 interfaceLevels_[leveli],
386 coarseSources[leveli],
387 cmpt
388 );
389 }
390
391 // Only add the preSmoothedCoarseCorrField if pre-smoothing is
392 // used
393 if (nPreSweeps_)
394 {
395 coarseCorrFields[leveli] += preSmoothedCoarseCorrField;
396 }
397
398 smoothers[leveli + 1].scalarSmooth
399 (
400 coarseCorrFields[leveli],
401 coarseSources[leveli], //coarseSource,
402 cmpt,
403 min
404 (
405 nPostSweeps_ + postSweepsLevelMultiplier_*leveli,
406 maxPostSweeps_
407 )
408 );
409 }
410 }
411
412 // Prolong the finest level correction
413 agglomeration_.prolongField
414 (
415 finestCorrection,
416 coarseCorrFields[0],
417 0,
418 true
419 );
420
421 if (interpolateCorrection_)
422 {
424 (
425 finestCorrection,
426 Apsi,
427 matrix_,
428 interfaceBouCoeffs_,
429 interfaces_,
430 agglomeration_.restrictAddressing(0),
431 coarseCorrFields[0],
432 cmpt
433 );
434 }
435
436 if (scaleCorrection_)
437 {
438 // Scale the finest level correction
439 scale
440 (
441 finestCorrection,
442 Apsi,
443 matrix_,
444 interfaceBouCoeffs_,
445 interfaces_,
446 finestResidual,
447 cmpt
448 );
449 }
450
451 forAll(psi, i)
452 {
453 psi[i] += finestCorrection[i];
454 }
455
456 smoothers[0].smooth
457 (
458 psi,
459 source,
460 cmpt,
461 nFinestSweeps_
462 );
463}
464
465
466void Foam::GAMGSolver::initVcycle
467(
468 PtrList<solveScalarField>& coarseCorrFields,
469 PtrList<solveScalarField>& coarseSources,
471 solveScalarField& scratch1,
472 solveScalarField& scratch2
473) const
474{
475 label maxSize = matrix_.diag().size();
476
477 coarseCorrFields.setSize(matrixLevels_.size());
478 coarseSources.setSize(matrixLevels_.size());
479 smoothers.setSize(matrixLevels_.size() + 1);
480
481 // Create the smoother for the finest level
482 smoothers.set
483 (
484 0,
486 (
487 fieldName_,
488 matrix_,
489 interfaceBouCoeffs_,
490 interfaceIntCoeffs_,
491 interfaces_,
492 controlDict_
493 )
494 );
495
496 forAll(matrixLevels_, leveli)
497 {
498 if (agglomeration_.nCells(leveli) >= 0)
499 {
500 label nCoarseCells = agglomeration_.nCells(leveli);
501
502 coarseSources.set(leveli, new solveScalarField(nCoarseCells));
503 }
504
505 if (matrixLevels_.set(leveli))
506 {
507 const lduMatrix& mat = matrixLevels_[leveli];
508
509 label nCoarseCells = mat.diag().size();
510
511 maxSize = max(maxSize, nCoarseCells);
512
513 coarseCorrFields.set(leveli, new solveScalarField(nCoarseCells));
514
515 smoothers.set
516 (
517 leveli + 1,
519 (
520 fieldName_,
521 matrixLevels_[leveli],
522 interfaceLevelsBouCoeffs_[leveli],
523 interfaceLevelsIntCoeffs_[leveli],
524 interfaceLevels_[leveli],
525 controlDict_
526 )
527 );
528 }
529 }
530
531 if (maxSize > matrix_.diag().size())
532 {
533 // Allocate some scratch storage
534 scratch1.setSize(maxSize);
535 scratch2.setSize(maxSize);
536 }
537}
538
539
540Foam::dictionary Foam::GAMGSolver::PCGsolverDict
541(
542 const scalar tol,
543 const scalar relTol
544) const
545{
546 dictionary dict(IStringStream("solver PCG; preconditioner DIC;")());
547 dict.add("tolerance", tol);
548 dict.add("relTol", relTol);
549
550 return dict;
551}
552
553
554Foam::dictionary Foam::GAMGSolver::PBiCGStabSolverDict
555(
556 const scalar tol,
557 const scalar relTol
558) const
559{
560 dictionary dict(IStringStream("solver PBiCGStab; preconditioner DILU;")());
561 dict.add("tolerance", tol);
562 dict.add("relTol", relTol);
563
564 return dict;
565}
566
567
568void Foam::GAMGSolver::solveCoarsestLevel
569(
570 solveScalarField& coarsestCorrField,
571 const solveScalarField& coarsestSource
572) const
573{
574 const label coarsestLevel = matrixLevels_.size() - 1;
575
576 const label coarseComm = matrixLevels_[coarsestLevel].mesh().comm();
577
578 if (directSolveCoarsest_)
579 {
580 PrecisionAdaptor<scalar, solveScalar> tcorrField(coarsestCorrField);
581
582 coarsestLUMatrixPtr_->solve
583 (
584 tcorrField.ref(),
586 );
587 }
588 //else if
589 //(
590 // agglomeration_.processorAgglomerate()
591 // && procMatrixLevels_.set(coarsestLevel)
592 //)
593 //{
594 // //const labelList& agglomProcIDs = agglomeration_.agglomProcIDs
595 // //(
596 // // coarsestLevel
597 // //);
598 // //
599 // //scalarField allSource;
600 // //
601 // //globalIndex cellOffsets;
602 // //if (Pstream::myProcNo(coarseComm) == agglomProcIDs[0])
603 // //{
604 // // cellOffsets.offsets() =
605 // // agglomeration_.cellOffsets(coarsestLevel);
606 // //}
607 // //
608 // //cellOffsets.gather
609 // //(
610 // // coarseComm,
611 // // agglomProcIDs,
612 // // coarsestSource,
613 // // allSource
614 // //);
615 // //
616 // //scalarField allCorrField;
617 // //solverPerformance coarseSolverPerf;
618 //
619 // label solveComm = agglomeration_.procCommunicator(coarsestLevel);
620 //
621 // coarsestCorrField = 0;
622 // solverPerformance coarseSolverPerf;
623 //
624 // if (Pstream::myProcNo(solveComm) != -1)
625 // {
626 // const lduMatrix& allMatrix = procMatrixLevels_[coarsestLevel];
627 //
628 // {
629 // Pout<< "** Master:Solving on comm:" << solveComm
630 // << " with procs:" << UPstream::procID(solveComm) << endl;
631 //
632 // if (allMatrix.asymmetric())
633 // {
634 // coarseSolverPerf = PBiCGStab
635 // (
636 // "coarsestLevelCorr",
637 // allMatrix,
638 // procInterfaceLevelsBouCoeffs_[coarsestLevel],
639 // procInterfaceLevelsIntCoeffs_[coarsestLevel],
640 // procInterfaceLevels_[coarsestLevel],
641 // PBiCGStabSolverDict(tolerance_, relTol_)
642 // ).solve
643 // (
644 // coarsestCorrField,
645 // coarsestSource
646 // );
647 // }
648 // else
649 // {
650 // coarseSolverPerf = PCG
651 // (
652 // "coarsestLevelCorr",
653 // allMatrix,
654 // procInterfaceLevelsBouCoeffs_[coarsestLevel],
655 // procInterfaceLevelsIntCoeffs_[coarsestLevel],
656 // procInterfaceLevels_[coarsestLevel],
657 // PCGsolverDict(tolerance_, relTol_)
658 // ).solve
659 // (
660 // coarsestCorrField,
661 // coarsestSource
662 // );
663 // }
664 // }
665 // }
666 //
667 // Pout<< "done master solve." << endl;
668 //
669 // //// Scatter to all processors
670 // //coarsestCorrField.setSize(coarsestSource.size());
671 // //cellOffsets.scatter
672 // //(
673 // // coarseComm,
674 // // agglomProcIDs,
675 // // allCorrField,
676 // // coarsestCorrField
677 // //);
678 //
679 // if (debug >= 2)
680 // {
681 // coarseSolverPerf.print(Info.masterStream(coarseComm));
682 // }
683 //
684 // Pout<< "procAgglom: coarsestSource :" << coarsestSource << endl;
685 // Pout<< "procAgglom: coarsestCorrField:" << coarsestCorrField << endl;
686 //}
687 else
688 {
689 coarsestCorrField = 0;
690 const solverPerformance coarseSolverPerf
691 (
692 coarsestSolverPtr_->scalarSolve
693 (
694 coarsestCorrField,
695 coarsestSource
696 )
697 );
698
699 if ((log_ >= 2) || debug)
700 {
701 coarseSolverPerf.print(Info.masterStream(coarseComm));
702 }
703 }
704}
705
706
707// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
SubField< solveScalar > subField
Definition Field.H:183
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....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
const Type & finalResidual() const noexcept
Return final residual.
const labelType & nIterations() const noexcept
Return number of iterations.
void print(Ostream &os) const
Print summary of solver performance to the given stream.
const Type & initialResidual() const noexcept
Return initial residual.
bool checkConvergence(const Type &tolerance, const Type &relTolerance, const int logLevel=0)
Check, store and return convergence.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
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_
Definition lduMatrix.H:158
label maxIter_
Maximum number of iterations in the solver.
Definition lduMatrix.H:181
scalar tolerance_
Final convergence tolerance.
Definition lduMatrix.H:191
lduInterfaceFieldPtrsList interfaces_
Definition lduMatrix.H:161
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition lduMatrix.H:159
label minIter_
Minimum number of iterations in the solver.
Definition lduMatrix.H:176
const lduMatrix & matrix() const noexcept
Definition lduMatrix.H:328
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.
Definition lduMatrix.H:171
scalar relTol_
Convergence tolerance relative to the initial.
Definition lduMatrix.H:196
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
const scalarField & diag() const
Definition lduMatrix.C:195
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.
Definition lduMatrix.C:457
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
const volScalarField & psi
dynamicFvMesh & mesh
Namespace for handling debugging switches.
Definition debug.C:45
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
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.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299