Home | Trees | Indices | Help |
---|
|
Author: Lutz Gross, l.gross@uq.edu.au
Copyright: Copyright (c) 2003-2009 by University of Queensland Earth Systems Science Computational Center (ESSCC) http://www.uq.edu.au/esscc Primary Business: Queensland, Australia
License: Licensed under the Open Software License version 3.0 http://www.opensource.org/licenses/osl-3.0.php
|
|||
ArithmeticTuple Tuple supporting inplace update x+=y and scaling x=a*y where x,y is an ArithmeticTuple and a is a
float.
|
|||
CorrectionFailed Exception thrown if no convergence has been achieved in the solution correction scheme. |
|||
Defect Defines a non-linear defect F(x) of a variable x. |
|||
HomogeneousSaddlePointProblem This class provides a framework for solving linear homogeneous saddle point problems of the form: |
|||
IndefinitePreconditioner Exception thrown if the preconditioner is not positive definite. |
|||
IterationBreakDown Exception thrown if the iteration scheme encountered an incurable breakdown. |
|||
Locator Locator provides access to the values of data objects at a given spatial coordinate x. |
|||
MaxIterReached Exception thrown if the maximum number of iteration steps is reached. |
|||
NegativeNorm Exception thrown if a norm calculation returns a negative norm. |
|||
NoPDE Solves the following problem for u: |
|||
Projector The Projector is a factory which projects a discontinuous function onto a continuous function on a given domain. |
|||
SolverSchemeException This is a generic exception thrown by solvers. |
|||
TimeIntegrationManager A simple mechanism to manage time dependend values. |
|
|||
tuple
|
|
||
tuple
|
|
||
escript.Data of rank 0 |
|
||
escript.Data of rank 0 |
|
||
same type as the initial guess |
|
||
tuple
|
|
||
tuple
|
|
||
|
|||
|
|||
|
|
|||
__url__ =
|
|
Solver for Ax=b with a general operator A (more details required!). It uses the generalized minimum residual method (GMRES). The iteration is terminated if |r| <= atol+rtol*|r0| where r0 is the initial residual and |.| is the energy norm. In fact |r| = sqrt( bilinearform(r,r))
Warning:
|
Solver for Ax=b with a symmetric and positive definite operator A (more details required!). It uses the minimum residual method (MINRES) with preconditioner M providing an approximation of A. The iteration is terminated if |r| <= atol+rtol*|r0| where r0 is the initial residual and |.| is the energy norm. In fact |r| = sqrt( bilinearform(Msolve(r),r)) For details on the preconditioned conjugate gradient method see the book: Templates for the Solution of Linear Systems by R. Barrett, M. Berry, T.F. Chan, J. Demmel, J. Donato, J. Dongarra, V. Eijkhout, R. Pozo, C. Romine, and H. van der Vorst.
Warning:
|
Creates a mask on the Solution(domain) function space where the value is one for samples that touch the boundary tagged by tags. Usage: m=MaskFromBoundaryTag(domain, "left", "right")
|
Creates a mask on the Solution(domain) function space where the value is one for samples that touch regions tagged by tags. Usage: m=MaskFromTag(domain, "ham")
|
Solves a non-linear problem F(x)=0 for unknown x using the stopping criterion: norm(F(x) <= atol + rtol * norm(F(x0) where x0 is the initial guess.
|
Solver for Ax=b with a symmetric and positive definite operator A (more details required!). It uses the conjugate gradient method with preconditioner M providing an approximation of A. The iteration is terminated if |r| <= atol+rtol*|r0| where r0 is the initial residual and |.| is the energy norm. In fact |r| = sqrt( bilinearform(Msolve(r),r)) For details on the preconditioned conjugate gradient method see the book: Templates for the Solution of Linear Systems by R. Barrett, M. Berry, T.F. Chan, J. Demmel, J. Donato, J. Dongarra, V. Eijkhout, R. Pozo, C. Romine, and H. van der Vorst.
Warning:
|
Solver for Ax=b with a general operator A (more details required!). It uses the Transpose-Free Quasi-Minimal Residual method (TFQMR). The iteration is terminated if |r| <= atol+rtol*|r0| where r0 is the initial residual and |.| is the energy norm. In fact |r| = sqrt( bilinearform(r,r))
Warning:
|
Applies a sequence of Givens rotations (c,s) recursively to the vector
Warning:
|
Home | Trees | Indices | Help |
---|
Generated by Epydoc 3.0.1 on Mon Aug 10 10:04:37 2009 | http://epydoc.sourceforge.net |