lduMatrix.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::lduMatrix
26 
27 Description
28  lduMatrix is a general matrix class in which the coefficients are
29  stored as three arrays, one for the upper triangle, one for the
30  lower triangle and a third for the diagonal.
31 
32  Addressing arrays must be supplied for the upper and lower triangles.
33 
34  It might be better if this class were organised as a hierachy starting
35  from an empty matrix, then deriving diagonal, symmetric and asymmetric
36  matrices.
37 
38 SourceFiles
39  lduMatrixATmul.C
40  lduMatrix.C
41  lduMatrixTemplates.C
42  lduMatrixOperations.C
43  lduMatrixSolver.C
44  lduMatrixPreconditioner.C
45  lduMatrixTests.C
46  lduMatrixUpdateMatrixInterfaces.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef lduMatrix_H
51 #define lduMatrix_H
52 
53 #include "lduMesh.H"
54 #include "primitiveFieldsFwd.H"
55 #include "FieldField.H"
57 #include "typeInfo.H"
58 #include "autoPtr.H"
59 #include "runTimeSelectionTables.H"
60 #include "solverPerformance.H"
61 #include "InfoProxy.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 
68 // Forward declaration of friend functions and operators
69 
70 class lduMatrix;
71 Ostream& operator<<(Ostream&, const lduMatrix&);
72 
73 
74 /*---------------------------------------------------------------------------*\
75  Class lduMatrix Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class lduMatrix
79 {
80  // private data
81 
82  //- LDU mesh reference
83  const lduMesh& lduMesh_;
84 
85  //- Coefficients (not including interfaces)
87 
88 
89 public:
90 
91  //- Abstract base-class for lduMatrix solvers
92  class solver
93  {
94  protected:
95 
96  // Protected data
97 
103 
104  //- Dictionary of controls
106 
107  //- Maximum number of iterations in the solver
108  label maxIter_;
109 
110  //- Minimum number of iterations in the solver
111  label minIter_;
112 
113  //- Final convergence tolerance
114  scalar tolerance_;
115 
116  //- Convergence tolerance relative to the initial
117  scalar relTol_;
118 
119 
120  // Protected Member Functions
121 
122  //- Read the control parameters from the controlDict_
123  virtual void readControls();
124 
125 
126  public:
127 
128  //- Runtime type information
129  virtual const word& type() const = 0;
130 
131 
132  // Declare run-time constructor selection tables
133 
135  (
136  autoPtr,
137  solver,
138  symMatrix,
139  (
140  const word& fieldName,
141  const lduMatrix& matrix,
145  const dictionary& solverControls
146  ),
147  (
148  fieldName,
149  matrix,
152  interfaces,
153  solverControls
154  )
155  );
156 
158  (
159  autoPtr,
160  solver,
161  asymMatrix,
162  (
163  const word& fieldName,
164  const lduMatrix& matrix,
168  const dictionary& solverControls
169  ),
170  (
171  fieldName,
172  matrix,
175  interfaces,
176  solverControls
177  )
178  );
179 
180 
181  // Constructors
182 
183  solver
184  (
185  const word& fieldName,
186  const lduMatrix& matrix,
190  const dictionary& solverControls
191  );
192 
193  // Selectors
194 
195  //- Return a new solver
196  static autoPtr<solver> New
197  (
198  const word& fieldName,
199  const lduMatrix& matrix,
203  const dictionary& solverControls
204  );
205 
206 
207 
208  //- Destructor
209  virtual ~solver()
210  {}
211 
212 
213  // Member functions
214 
215  // Access
216 
217  const word& fieldName() const
218  {
219  return fieldName_;
220  }
221 
222  const lduMatrix& matrix() const
223  {
224  return matrix_;
225  }
226 
228  {
229  return interfaceBouCoeffs_;
230  }
231 
233  {
234  return interfaceIntCoeffs_;
235  }
236 
238  {
239  return interfaces_;
240  }
241 
242 
243  //- Read and reset the solver parameters from the given stream
244  virtual void read(const dictionary&);
245 
246  virtual solverPerformance solve
247  (
248  scalarField& psi,
249  const scalarField& source,
250  const direction cmpt=0
251  ) const = 0;
252 
253  //- Return the matrix norm used to normalise the residual for the
254  // stopping criterion
255  scalar normFactor
256  (
257  const scalarField& psi,
258  const scalarField& source,
259  const scalarField& Apsi,
260  scalarField& tmpField
261  ) const;
262  };
263 
264 
265  //- Abstract base-class for lduMatrix smoothers
266  class smoother
267  {
268  protected:
269 
270  // Protected data
271 
277 
278 
279  public:
280 
281  //- Find the smoother name (directly or from a sub-dictionary)
282  static word getName(const dictionary&);
283 
284  //- Runtime type information
285  virtual const word& type() const = 0;
286 
287 
288  // Declare run-time constructor selection tables
289 
291  (
292  autoPtr,
293  smoother,
294  symMatrix,
295  (
296  const word& fieldName,
297  const lduMatrix& matrix,
301  ),
302  (
303  fieldName,
304  matrix,
307  interfaces
308  )
309  );
310 
312  (
313  autoPtr,
314  smoother,
315  asymMatrix,
316  (
317  const word& fieldName,
318  const lduMatrix& matrix,
322  ),
323  (
324  fieldName,
325  matrix,
328  interfaces
329  )
330  );
331 
332 
333  // Constructors
334 
335  smoother
336  (
337  const word& fieldName,
338  const lduMatrix& matrix,
342  );
343 
344 
345  // Selectors
346 
347  //- Return a new smoother
348  static autoPtr<smoother> New
349  (
350  const word& fieldName,
351  const lduMatrix& matrix,
355  const dictionary& solverControls
356  );
357 
358 
359  //- Destructor
360  virtual ~smoother()
361  {}
362 
363 
364  // Member functions
365 
366  // Access
367 
368  const word& fieldName() const
369  {
370  return fieldName_;
371  }
372 
373  const lduMatrix& matrix() const
374  {
375  return matrix_;
376  }
377 
379  {
380  return interfaceBouCoeffs_;
381  }
382 
384  {
385  return interfaceIntCoeffs_;
386  }
387 
389  {
390  return interfaces_;
391  }
392 
393 
394  //- Smooth the solution for a given number of sweeps
395  virtual void smooth
396  (
397  scalarField& psi,
398  const scalarField& source,
399  const direction cmpt,
400  const label nSweeps
401  ) const = 0;
402  };
403 
404 
405  //- Abstract base-class for lduMatrix preconditioners
406  class preconditioner
407  {
408  protected:
409 
410  // Protected data
411 
412  //- Reference to the base-solver this preconditioner is used with
413  const solver& solver_;
414 
415 
416  public:
417 
418  //- Find the preconditioner name (directly or from a sub-dictionary)
419  static word getName(const dictionary&);
420 
421  //- Runtime type information
422  virtual const word& type() const = 0;
423 
424 
425  // Declare run-time constructor selection tables
426 
428  (
429  autoPtr,
431  symMatrix,
432  (
433  const solver& sol,
434  const dictionary& solverControls
435  ),
436  (sol, solverControls)
437  );
438 
440  (
441  autoPtr,
443  asymMatrix,
444  (
445  const solver& sol,
446  const dictionary& solverControls
447  ),
448  (sol, solverControls)
449  );
450 
451 
452  // Constructors
453 
455  (
456  const solver& sol
457  )
458  :
459  solver_(sol)
460  {}
461 
462 
463  // Selectors
464 
465  //- Return a new preconditioner
467  (
468  const solver& sol,
469  const dictionary& solverControls
470  );
471 
472 
473  //- Destructor
474  virtual ~preconditioner()
475  {}
476 
477 
478  // Member functions
479 
480  //- Read and reset the preconditioner parameters
481  // from the given stream
482  virtual void read(const dictionary&)
483  {}
484 
485  //- Return wA the preconditioned form of residual rA
486  virtual void precondition
487  (
488  scalarField& wA,
489  const scalarField& rA,
490  const direction cmpt=0
491  ) const = 0;
492 
493  //- Return wT the transpose-matrix preconditioned form of
494  // residual rT.
495  // This is only required for preconditioning asymmetric matrices.
496  virtual void preconditionT
497  (
498  scalarField& wT,
499  const scalarField& rT,
500  const direction cmpt=0
501  ) const
502  {
504  }
505  };
506 
507 
508  // Static data
509 
510  // Declare name of the class and its debug switch
511  ClassName("lduMatrix");
512 
513 
514  // Constructors
515 
516  //- Construct given an LDU addressed mesh.
517  // The coefficients are initially empty for subsequent setting.
518  lduMatrix(const lduMesh&);
519 
520  //- Construct as copy
521  lduMatrix(const lduMatrix&);
522 
523  //- Construct as copy or re-use as specified.
524  lduMatrix(lduMatrix&, bool reUse);
525 
526  //- Construct given an LDU addressed mesh and an Istream
527  // from which the coefficients are read
528  lduMatrix(const lduMesh&, Istream&);
529 
530 
531  //- Destructor
532  ~lduMatrix();
533 
534 
535  // Member functions
536 
537  // Access to addressing
538 
539  //- Return the LDU mesh from which the addressing is obtained
540  const lduMesh& mesh() const
541  {
542  return lduMesh_;
543  }
544 
545  //- Return the LDU addressing
546  const lduAddressing& lduAddr() const
547  {
548  return lduMesh_.lduAddr();
549  }
550 
551  //- Return the patch evaluation schedule
552  const lduSchedule& patchSchedule() const
553  {
554  return lduAddr().patchSchedule();
555  }
556 
557 
558  // Access to coefficients
559 
560  scalarField& lower();
561  scalarField& diag();
562  scalarField& upper();
563 
564  // Size with externally provided sizes (for constructing with 'fake'
565  // mesh in GAMG)
566 
567  scalarField& lower(const label size);
568  scalarField& diag(const label nCoeffs);
569  scalarField& upper(const label nCoeffs);
570 
571 
572  const scalarField& lower() const;
573  const scalarField& diag() const;
574  const scalarField& upper() const;
575 
576  bool hasDiag() const
577  {
578  return (diagPtr_);
579  }
580 
581  bool hasUpper() const
582  {
583  return (upperPtr_);
584  }
585 
586  bool hasLower() const
587  {
588  return (lowerPtr_);
589  }
590 
591  bool diagonal() const
592  {
593  return (diagPtr_ && !lowerPtr_ && !upperPtr_);
594  }
595 
596  bool symmetric() const
597  {
598  return (diagPtr_ && (!lowerPtr_ && upperPtr_));
599  }
600 
601  bool asymmetric() const
602  {
603  return (diagPtr_ && lowerPtr_ && upperPtr_);
604  }
605 
606 
607  // operations
608 
609  void sumDiag();
610  void negSumDiag();
611 
612  void sumMagOffDiag(scalarField& sumOff) const;
613 
614  //- Matrix multiplication with updated interfaces.
615  void Amul
616  (
617  scalarField&,
618  const tmp<scalarField>&,
621  const direction cmpt
622  ) const;
623 
624  //- Matrix transpose multiplication with updated interfaces.
625  void Tmul
626  (
627  scalarField&,
628  const tmp<scalarField>&,
631  const direction cmpt
632  ) const;
633 
634 
635  //- Sum the coefficients on each row of the matrix
636  void sumA
637  (
638  scalarField&,
641  ) const;
642 
643 
644  void residual
645  (
646  scalarField& rA,
647  const scalarField& psi,
648  const scalarField& source,
649  const FieldField<Field, scalar>& interfaceBouCoeffs,
650  const lduInterfaceFieldPtrsList& interfaces,
651  const direction cmpt
652  ) const;
653 
655  (
656  const scalarField& psi,
657  const scalarField& source,
658  const FieldField<Field, scalar>& interfaceBouCoeffs,
659  const lduInterfaceFieldPtrsList& interfaces,
660  const direction cmpt
661  ) const;
662 
663 
664  //- Initialise the update of interfaced interfaces
665  // for matrix operations
667  (
668  const FieldField<Field, scalar>& interfaceCoeffs,
669  const lduInterfaceFieldPtrsList& interfaces,
670  const scalarField& psiif,
671  scalarField& result,
672  const direction cmpt
673  ) const;
674 
675  //- Update interfaced interfaces for matrix operations
677  (
678  const FieldField<Field, scalar>& interfaceCoeffs,
679  const lduInterfaceFieldPtrsList& interfaces,
680  const scalarField& psiif,
681  scalarField& result,
682  const direction cmpt
683  ) const;
684 
685 
686  template<class Type>
687  tmp<Field<Type> > H(const Field<Type>&) const;
688 
689  template<class Type>
690  tmp<Field<Type> > H(const tmp<Field<Type> >&) const;
691 
692  tmp<scalarField> H1() const;
693 
694  template<class Type>
695  tmp<Field<Type> > faceH(const Field<Type>&) const;
696 
697  template<class Type>
698  tmp<Field<Type> > faceH(const tmp<Field<Type> >&) const;
699 
700 
701  // Info
702 
703  //- Return info proxy.
704  // Used to print matrix information to a stream
705  InfoProxy<lduMatrix> info() const
706  {
707  return *this;
708  }
709 
710 
711  // Member operators
712 
713  void operator=(const lduMatrix&);
714 
715  void negate();
716 
717  void operator+=(const lduMatrix&);
718  void operator-=(const lduMatrix&);
719 
720  void operator*=(const scalarField&);
721  void operator*=(scalar);
722 
723 
724  // Ostream operator
725 
726  friend Ostream& operator<<(Ostream&, const lduMatrix&);
728 };
729 
730 
731 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
732 
733 } // End namespace Foam
734 
735 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
736 
737 #ifdef NoRepository
738 # include "lduMatrixTemplates.C"
739 #endif
740 
741 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
742 
743 #endif
744 
745 // ************************************************************************* //
Foam::lduMatrix::solver::solve
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const =0
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:111
Foam::lduMatrix::preconditioner::getName
static word getName(const dictionary &)
Find the preconditioner name (directly or from a sub-dictionary)
Definition: lduMatrixPreconditioner.C:40
Foam::lduMatrix::operator-=
void operator-=(const lduMatrix &)
Definition: lduMatrixOperations.C:223
Foam::lduMatrix::lduAddr
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:545
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
primitiveFieldsFwd.H
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
Foam::lduMatrix::solver::normFactor
scalar normFactor(const scalarField &psi, const scalarField &source, const scalarField &Apsi, scalarField &tmpField) const
Return the matrix norm used to normalise the residual for the.
Definition: lduMatrixSolver.C:175
Foam::lduMatrix::solver::interfaces
const lduInterfaceFieldPtrsList & interfaces() const
Definition: lduMatrix.H:236
Foam::lduMatrix::info
InfoProxy< lduMatrix > info() const
Return info proxy.
Definition: lduMatrix.H:704
Foam::lduMatrix::preconditioner::~preconditioner
virtual ~preconditioner()
Destructor.
Definition: lduMatrix.H:473
Foam::lduMatrix::hasLower
bool hasLower() const
Definition: lduMatrix.H:585
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::InfoProxy
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:45
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::lduMatrix::Amul
void Amul(scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix multiplication with updated interfaces.
Definition: lduMatrixATmul.C:35
typeInfo.H
FieldField.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::lduMatrix::smoother::declareRunTimeSelectionTable
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))
Foam::lduMatrix::operator+=
void operator+=(const lduMatrix &)
Definition: lduMatrixOperations.C:144
Foam::lduMatrix::~lduMatrix
~lduMatrix()
Destructor.
Definition: lduMatrix.C:146
Foam::lduMatrix::solver::~solver
virtual ~solver()
Destructor.
Definition: lduMatrix.H:208
Foam::lduMatrix::preconditioner::read
virtual void read(const dictionary &)
Read and reset the preconditioner parameters.
Definition: lduMatrix.H:481
Foam::lduMatrix::ClassName
ClassName("lduMatrix")
InfoProxy.H
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:77
Foam::lduMatrix::sumMagOffDiag
void sumMagOffDiag(scalarField &sumOff) const
Definition: lduMatrixOperations.C:68
Foam::lduMatrix::hasDiag
bool hasDiag() const
Definition: lduMatrix.H:575
Foam::lduMatrix::smoother::interfaceBouCoeffs
const FieldField< Field, scalar > & interfaceBouCoeffs() const
Definition: lduMatrix.H:377
Foam::lduMatrix::symmetric
bool symmetric() const
Definition: lduMatrix.H:595
Foam::lduMatrix::preconditioner::New
static autoPtr< preconditioner > New(const solver &sol, const dictionary &solverControls)
Return a new preconditioner.
Definition: lduMatrixPreconditioner.C:63
Foam::lduMatrix::solver::matrix
const lduMatrix & matrix() const
Definition: lduMatrix.H:221
Foam::lduMatrix::solver
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:91
lduMatrixTemplates.C
lduMatrix member H operations.
Foam::lduMatrix::smoother::fieldName
const word & fieldName() const
Definition: lduMatrix.H:367
Foam::lduMatrix::updateMatrixInterfaces
void updateMatrixInterfaces(const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const
Update interfaced interfaces for matrix operations.
Definition: lduMatrixUpdateMatrixInterfaces.C:97
Foam::lduMatrix::smoother::matrix_
const lduMatrix & matrix_
Definition: lduMatrix.H:272
lduInterfaceFieldPtrsList.H
Foam::lduMatrix::upper
scalarField & upper()
Definition: lduMatrix.C:194
Foam::lduMatrix::smoother::type
virtual const word & type() const =0
Runtime type information.
Foam::lduMatrix::lowerPtr_
scalarField * lowerPtr_
Coefficients (not including interfaces)
Definition: lduMatrix.H:85
Foam::lduMatrix::lduMatrix
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:40
Foam::lduMatrix::smoother::interfaceIntCoeffs
const FieldField< Field, scalar > & interfaceIntCoeffs() const
Definition: lduMatrix.H:382
Foam::lduMatrix::solver::interfaceBouCoeffs_
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:99
Foam::lduMatrix::negSumDiag
void negSumDiag()
Definition: lduMatrixOperations.C:50
Foam::lduMatrix::solver::minIter_
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:110
Foam::lduMatrix::lower
scalarField & lower()
Definition: lduMatrix.C:165
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:365
Foam::lduMatrix::smoother::interfaceBouCoeffs_
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:273
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::lduMatrix::solver::read
virtual void read(const dictionary &)
Read and reset the solver parameters from the given stream.
Definition: lduMatrixSolver.C:167
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::lduMatrix::smoother::matrix
const lduMatrix & matrix() const
Definition: lduMatrix.H:372
Foam::lduMatrix::initMatrixInterfaces
void initMatrixInterfaces(const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const
Initialise the update of interfaced interfaces.
Definition: lduMatrixUpdateMatrixInterfaces.C:31
Foam::lduMatrix::solver::maxIter_
label maxIter_
Maximum number of iterations in the solver.
Definition: lduMatrix.H:107
Foam::UPtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:53
solverPerformance.H
Foam::lduAddressing::patchSchedule
virtual const lduSchedule & patchSchedule() const =0
Foam::lduMatrix::diagonal
bool diagonal() const
Definition: lduMatrix.H:590
Foam::operator<<
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:130
Foam::lduMatrix::solver::fieldName_
word fieldName_
Definition: lduMatrix.H:97
Foam::lduMatrix::H1
tmp< scalarField > H1() const
Definition: lduMatrixATmul.C:298
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::lduMatrix::solver::readControls
virtual void readControls()
Read the control parameters from the controlDict_.
Definition: lduMatrixSolver.C:158
Foam::lduMatrix::smoother
Abstract base-class for lduMatrix smoothers.
Definition: lduMatrix.H:265
Foam::lduMatrix::preconditioner::preconditionT
virtual void preconditionT(scalarField &wT, const scalarField &rT, const direction cmpt=0) const
Return wT the transpose-matrix preconditioned form of.
Definition: lduMatrix.H:496
Foam::lduMatrix::smoother::getName
static word getName(const dictionary &)
Find the smoother name (directly or from a sub-dictionary)
Definition: lduMatrixSmoother.C:40
Foam::lduMatrix::solver::New
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)
Return a new solver.
Definition: lduMatrixSolver.C:41
Foam::lduMatrix::negate
void negate()
Definition: lduMatrixOperations.C:125
Foam::lduMatrix::diagPtr_
scalarField * diagPtr_
Definition: lduMatrix.H:85
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::lduMatrix::H
tmp< Field< Type > > H(const Field< Type > &) const
Foam::lduMatrix::solver::relTol_
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:116
Foam::lduMatrix::smoother::New
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.
Definition: lduMatrixSmoother.C:62
Foam::lduMatrix::solver::matrix_
const lduMatrix & matrix_
Definition: lduMatrix.H:98
Foam::lduMatrix::sumA
void sumA(scalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
Sum the coefficients on each row of the matrix.
Definition: lduMatrixATmul.C:155
Foam::lduMatrix::solver::interfaceIntCoeffs_
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:100
Foam::lduMatrix::diag
scalarField & diag()
Definition: lduMatrix.C:183
Foam::lduMatrix::residual
void residual(scalarField &rA, const scalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
Definition: lduMatrixATmul.C:204
Foam::lduMatrix::patchSchedule
const lduSchedule & patchSchedule() const
Return the patch evaluation schedule.
Definition: lduMatrix.H:551
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::lduMatrix::solver::controlDict_
dictionary controlDict_
Dictionary of controls.
Definition: lduMatrix.H:104
Foam::lduMatrix::solver::tolerance_
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:113
Foam::lduMatrix::solver::interfaces_
lduInterfaceFieldPtrsList interfaces_
Definition: lduMatrix.H:101
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::lduMatrix::operator=
void operator=(const lduMatrix &)
Definition: lduMatrixOperations.C:88
Foam::lduMatrix::sumDiag
void sumDiag()
Definition: lduMatrixOperations.C:33
Foam::lduMatrix::smoother::interfaces_
const lduInterfaceFieldPtrsList & interfaces_
Definition: lduMatrix.H:275
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::lduMatrix::smoother::interfaces
const lduInterfaceFieldPtrsList & interfaces() const
Definition: lduMatrix.H:387
Foam::lduMatrix::smoother::smoother
smoother(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces)
Definition: lduMatrixSmoother.C:156
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::lduMatrix::solver::declareRunTimeSelectionTable
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))
Foam::lduMatrix::solver::interfaceBouCoeffs
const FieldField< Field, scalar > & interfaceBouCoeffs() const
Definition: lduMatrix.H:226
Foam::lduMatrix::operator*=
void operator*=(const scalarField &)
Definition: lduMatrixOperations.C:302
Foam::lduMatrix::faceH
tmp< Field< Type > > faceH(const Field< Type > &) const
Foam::lduMatrix::operator<<
friend Ostream & operator<<(Ostream &, const lduMatrix &)
Foam::lduMatrix::mesh
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:539
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::lduMatrix::preconditioner::preconditioner
preconditioner(const solver &sol)
Definition: lduMatrix.H:454
Foam::lduMatrix::solver::fieldName
const word & fieldName() const
Definition: lduMatrix.H:216
Foam::lduMatrix::smoother::fieldName_
word fieldName_
Definition: lduMatrix.H:271
Foam::lduMatrix::asymmetric
bool asymmetric() const
Definition: lduMatrix.H:600
Foam::lduMatrix::hasUpper
bool hasUpper() const
Definition: lduMatrix.H:580
Foam::lduMatrix::preconditioner::type
virtual const word & type() const =0
Runtime type information.
Foam::lduMatrix::lduMesh_
const lduMesh & lduMesh_
LDU mesh reference.
Definition: lduMatrix.H:82
Foam::lduMatrix::smoother::~smoother
virtual ~smoother()
Destructor.
Definition: lduMatrix.H:359
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::lduMatrix::preconditioner
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:405
Foam::lduMatrix::preconditioner::solver_
const solver & solver_
Reference to the base-solver this preconditioner is used with.
Definition: lduMatrix.H:412
lduMesh.H
Foam::lduMatrix::solver::interfaceIntCoeffs
const FieldField< Field, scalar > & interfaceIntCoeffs() const
Definition: lduMatrix.H:231
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:49
Foam::lduMatrix::solver::type
virtual const word & type() const =0
Runtime type information.
Foam::lduMatrix::smoother::smooth
virtual void smooth(scalarField &psi, const scalarField &source, const direction cmpt, const label nSweeps) const =0
Smooth the solution for a given number of sweeps.
Foam::lduMesh::lduAddr
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
Foam::lduMatrix::upperPtr_
scalarField * upperPtr_
Definition: lduMatrix.H:85
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:51
Foam::lduMatrix::solver::solver
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Definition: lduMatrixSolver.C:136
Foam::lduMatrix::Tmul
void Tmul(scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix transpose multiplication with updated interfaces.
Definition: lduMatrixATmul.C:96
Foam::lduMatrix::preconditioner::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, preconditioner, symMatrix,(const solver &sol, const dictionary &solverControls),(sol, solverControls))
Foam::lduMatrix::smoother::interfaceIntCoeffs_
const FieldField< Field, scalar > & interfaceIntCoeffs_
Definition: lduMatrix.H:274
autoPtr.H
Foam::lduMatrix::preconditioner::precondition
virtual void precondition(scalarField &wA, const scalarField &rA, const direction cmpt=0) const =0
Return wA the preconditioned form of residual rA.