/usr/include/freefoam/OpenFOAM/lduMatrix.H is in libfreefoam-dev 0.1.0+dfsg-1build1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 | /*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
\\/ 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::lduMatrix
Description
lduMatrix is a general matrix class in which the coefficients are
stored as three arrays, one for the upper triangle, one for the
lower triangle and a third for the diagonal.
Addressing arrays must be supplied for the upper and lower triangles.
It might be better if this class were organised as a hierachy starting
from an empty matrix, then deriving diagonal, symmetric and asymmetric
matrices.
SourceFiles
lduMatrixATmul.C
lduMatrix.C
lduMatrixTemplates.C
lduMatrixOperations.C
lduMatrixSolver.C
lduMatrixPreconditioner.C
lduMatrixTests.C
lduMatrixUpdateMatrixInterfaces.C
\*---------------------------------------------------------------------------*/
#ifndef lduMatrix_H
#define lduMatrix_H
#include <OpenFOAM/lduMesh.H>
#include <OpenFOAM/primitiveFieldsFwd.H>
#include <OpenFOAM/FieldField.H>
#include <OpenFOAM/lduInterfaceFieldPtrsList.H>
#include <OpenFOAM/typeInfo.H>
#include <OpenFOAM/autoPtr.H>
#include <OpenFOAM/runTimeSelectionTables.H>
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// Forward declaration of friend functions and operators
class lduMatrix;
Ostream& operator<<(Ostream&, const lduMatrix&);
/*---------------------------------------------------------------------------*\
Class lduMatrix Declaration
\*---------------------------------------------------------------------------*/
class lduMatrix
{
// private data
//- LDU mesh reference
const lduMesh& lduMesh_;
//- Coefficients (not including interfaces)
scalarField *lowerPtr_, *diagPtr_, *upperPtr_;
public:
//- Class returned by the solver, containing performance statistics
class solverPerformance
{
word solverName_;
word fieldName_;
scalar initialResidual_;
scalar finalResidual_;
label noIterations_;
bool converged_;
bool singular_;
public:
// Constructors
solverPerformance()
:
initialResidual_(0),
finalResidual_(0),
noIterations_(0),
converged_(false),
singular_(false)
{}
solverPerformance
(
const word& solverName,
const word& fieldName,
const scalar iRes = 0,
const scalar fRes = 0,
const label nIter = 0,
const bool converged = false,
const bool singular = false
)
:
solverName_(solverName),
fieldName_(fieldName),
initialResidual_(iRes),
finalResidual_(fRes),
noIterations_(nIter),
converged_(converged),
singular_(singular)
{}
// Member functions
//- Return solver name
const word& solverName() const
{
return solverName_;
}
//- Return initial residual
scalar initialResidual() const
{
return initialResidual_;
}
//- Return initial residual
scalar& initialResidual()
{
return initialResidual_;
}
//- Return final residual
scalar finalResidual() const
{
return finalResidual_;
}
//- Return final residual
scalar& finalResidual()
{
return finalResidual_;
}
//- Return number of iterations
label nIterations() const
{
return noIterations_;
}
//- Return number of iterations
label& nIterations()
{
return noIterations_;
}
//- Has the solver converged?
bool converged() const
{
return converged_;
}
//- Is the matrix singular?
bool singular() const
{
return singular_;
}
//- Convergence test
bool checkConvergence
(
const scalar tolerance,
const scalar relTolerance
);
//- Singularity test
bool checkSingularity(const scalar residual);
//- Print summary of solver performance
void print() const;
};
//- Abstract base-class for lduMatrix solvers
class solver
{
protected:
// Protected data
word fieldName_;
const lduMatrix& matrix_;
const FieldField<Field, scalar>& interfaceBouCoeffs_;
const FieldField<Field, scalar>& interfaceIntCoeffs_;
lduInterfaceFieldPtrsList interfaces_;
//- dictionary of controls
dictionary controlDict_;
//- Maximum number of iterations in the solver
label maxIter_;
//- Final convergence tolerance
scalar tolerance_;
//- Convergence tolerance relative to the initial
scalar relTol_;
// Protected Member Functions
//- Read the control parameters from the controlDict_
virtual void readControls();
public:
//- Runtime type information
virtual const word& type() const = 0;
// Declare run-time constructor selection tables
declareRunTimeSelectionTable
(
autoPtr,
solver,
symMatrix,
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const FieldField<Field, scalar>& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const dictionary& solverControls
),
(
fieldName,
matrix,
interfaceBouCoeffs,
interfaceIntCoeffs,
interfaces,
solverControls
)
);
declareRunTimeSelectionTable
(
autoPtr,
solver,
asymMatrix,
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const FieldField<Field, scalar>& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const dictionary& solverControls
),
(
fieldName,
matrix,
interfaceBouCoeffs,
interfaceIntCoeffs,
interfaces,
solverControls
)
);
// Constructors
solver
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const FieldField<Field, scalar>& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const dictionary& solverControls
);
// Selectors
//- Return a new solver
static autoPtr<solver> New
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const FieldField<Field, scalar>& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const dictionary& solverControls
);
// Destructor
virtual ~solver()
{}
// Member functions
// Access
const word& fieldName() const
{
return fieldName_;
}
const lduMatrix& matrix() const
{
return matrix_;
}
const FieldField<Field, scalar>& interfaceBouCoeffs() const
{
return interfaceBouCoeffs_;
}
const FieldField<Field, scalar>& interfaceIntCoeffs() const
{
return interfaceIntCoeffs_;
}
const lduInterfaceFieldPtrsList& interfaces() const
{
return interfaces_;
}
//- Read and reset the solver parameters from the given stream
virtual void read(const dictionary&);
virtual solverPerformance solve
(
scalarField& psi,
const scalarField& source,
const direction cmpt=0
) const = 0;
//- Return the matrix norm used to normalise the residual for the
// stopping criterion
scalar normFactor
(
const scalarField& psi,
const scalarField& source,
const scalarField& Apsi,
scalarField& tmpField
) const;
};
//- Abstract base-class for lduMatrix smoothers
class smoother
{
protected:
// Protected data
word fieldName_;
const lduMatrix& matrix_;
const FieldField<Field, scalar>& interfaceBouCoeffs_;
const FieldField<Field, scalar>& interfaceIntCoeffs_;
const lduInterfaceFieldPtrsList& interfaces_;
public:
//- Find the smoother name (directly or from a sub-dictionary)
static word getName(const dictionary&);
//- Runtime type information
virtual const word& type() const = 0;
// Declare run-time constructor selection tables
declareRunTimeSelectionTable
(
autoPtr,
smoother,
symMatrix,
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const FieldField<Field, scalar>& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces
),
(
fieldName,
matrix,
interfaceBouCoeffs,
interfaceIntCoeffs,
interfaces
)
);
declareRunTimeSelectionTable
(
autoPtr,
smoother,
asymMatrix,
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const FieldField<Field, scalar>& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces
),
(
fieldName,
matrix,
interfaceBouCoeffs,
interfaceIntCoeffs,
interfaces
)
);
// Constructors
smoother
(
const word& fieldName,
const lduMatrix& matrix,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const FieldField<Field, scalar>& interfaceIntCoeffs,
const lduInterfaceFieldPtrsList& interfaces
);
// Selectors
//- Return a new smoother
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
);
// Destructor
virtual ~smoother()
{}
// Member functions
// Access
const word& fieldName() const
{
return fieldName_;
}
const lduMatrix& matrix() const
{
return matrix_;
}
const FieldField<Field, scalar>& interfaceBouCoeffs() const
{
return interfaceBouCoeffs_;
}
const FieldField<Field, scalar>& interfaceIntCoeffs() const
{
return interfaceIntCoeffs_;
}
const lduInterfaceFieldPtrsList& interfaces() const
{
return interfaces_;
}
//- Smooth the solution for a given number of sweeps
virtual void smooth
(
scalarField& psi,
const scalarField& source,
const direction cmpt,
const label nSweeps
) const = 0;
};
//- Abstract base-class for lduMatrix preconditioners
class preconditioner
{
protected:
// Protected data
//- Reference to the base-solver this preconditioner is used with
const solver& solver_;
public:
//- Find the preconditioner name (directly or from a sub-dictionary)
static word getName(const dictionary&);
//- Runtime type information
virtual const word& type() const = 0;
// Declare run-time constructor selection tables
declareRunTimeSelectionTable
(
autoPtr,
preconditioner,
symMatrix,
(
const solver& sol,
const dictionary& solverControls
),
(sol, solverControls)
);
declareRunTimeSelectionTable
(
autoPtr,
preconditioner,
asymMatrix,
(
const solver& sol,
const dictionary& solverControls
),
(sol, solverControls)
);
// Constructors
preconditioner
(
const solver& sol
)
:
solver_(sol)
{}
// Selectors
//- Return a new preconditioner
static autoPtr<preconditioner> New
(
const solver& sol,
const dictionary& solverControls
);
// Destructor
virtual ~preconditioner()
{}
// Member functions
//- Read and reset the preconditioner parameters
// from the given stream
virtual void read(const dictionary&)
{}
//- Return wA the preconditioned form of residual rA
virtual void precondition
(
scalarField& wA,
const scalarField& rA,
const direction cmpt=0
) const = 0;
//- Return wT the transpose-matrix preconditioned form of
// residual rT.
// This is only required for preconditioning asymmetric matrices.
virtual void preconditionT
(
scalarField& wT,
const scalarField& rT,
const direction cmpt=0
) const
{
notImplemented
(
type() +"::preconditionT"
"(scalarField& wT, const scalarField& rT, "
"const direction cmpt)"
);
}
};
// Static data
// Declare name of the class and its debug switch
ClassName("lduMatrix");
//- Large scalar for the use in solvers
static const scalar great_;
//- Small scalar for the use in solvers
static const scalar small_;
// Constructors
//- Construct given an LDU addressed mesh.
// The coefficients are initially empty for subsequent setting.
lduMatrix(const lduMesh&);
//- Construct as copy
lduMatrix(const lduMatrix&);
//- Construct as copy or re-use as specified.
lduMatrix(lduMatrix&, bool reUse);
//- Construct given an LDU addressed mesh and an Istream
// from which the coefficients are read
lduMatrix(const lduMesh&, Istream&);
// Destructor
~lduMatrix();
// Member functions
// Access to addressing
//- Return the LDU mesh from which the addressing is obtained
const lduMesh& mesh() const
{
return lduMesh_;
}
//- Return the LDU addressing
const lduAddressing& lduAddr() const
{
return lduMesh_.lduAddr();
}
//- Return the patch evaluation schedule
const lduSchedule& patchSchedule() const
{
return lduAddr().patchSchedule();
}
// Access to coefficients
scalarField& lower();
scalarField& diag();
scalarField& upper();
const scalarField& lower() const;
const scalarField& diag() const;
const scalarField& upper() const;
bool hasDiag() const
{
return (diagPtr_);
}
bool hasUpper() const
{
return (upperPtr_);
}
bool hasLower() const
{
return (lowerPtr_);
}
bool diagonal() const
{
return (diagPtr_ && !lowerPtr_ && !upperPtr_);
}
bool symmetric() const
{
return (diagPtr_ && (!lowerPtr_ && upperPtr_));
}
bool asymmetric() const
{
return (diagPtr_ && lowerPtr_ && upperPtr_);
}
// operations
void sumDiag();
void negSumDiag();
void sumMagOffDiag(scalarField& sumOff) const;
//- Matrix multiplication with updated interfaces.
void Amul
(
scalarField&,
const tmp<scalarField>&,
const FieldField<Field, scalar>&,
const lduInterfaceFieldPtrsList&,
const direction cmpt
) const;
//- Matrix transpose multiplication with updated interfaces.
void Tmul
(
scalarField&,
const tmp<scalarField>&,
const FieldField<Field, scalar>&,
const lduInterfaceFieldPtrsList&,
const direction cmpt
) const;
//- Sum the coefficients on each row of the matrix
void sumA
(
scalarField&,
const FieldField<Field, scalar>&,
const lduInterfaceFieldPtrsList&
) const;
void residual
(
scalarField& rA,
const scalarField& psi,
const scalarField& source,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const direction cmpt
) const;
tmp<scalarField> residual
(
const scalarField& psi,
const scalarField& source,
const FieldField<Field, scalar>& interfaceBouCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const direction cmpt
) const;
//- Initialise the update of interfaced interfaces
// for matrix operations
void initMatrixInterfaces
(
const FieldField<Field, scalar>& interfaceCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const scalarField& psiif,
scalarField& result,
const direction cmpt
) const;
//- Update interfaced interfaces for matrix operations
void updateMatrixInterfaces
(
const FieldField<Field, scalar>& interfaceCoeffs,
const lduInterfaceFieldPtrsList& interfaces,
const scalarField& psiif,
scalarField& result,
const direction cmpt
) const;
template<class Type>
tmp<Field<Type> > H(const Field<Type>&) const;
template<class Type>
tmp<Field<Type> > H(const tmp<Field<Type> >&) const;
tmp<scalarField> H1() const;
template<class Type>
tmp<Field<Type> > faceH(const Field<Type>&) const;
template<class Type>
tmp<Field<Type> > faceH(const tmp<Field<Type> >&) const;
// Member operators
void operator=(const lduMatrix&);
void negate();
void operator+=(const lduMatrix&);
void operator-=(const lduMatrix&);
void operator*=(const scalarField&);
void operator*=(scalar);
// Ostream operator
friend Ostream& operator<<(Ostream&, const lduMatrix&);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
# include <OpenFOAM/lduMatrixTemplates.C>
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************ vim: set sw=4 sts=4 et: ************************ //
|