MatCreateLMVMSymBadBroyden#
Creates a limited-memory Symmetric “Bad” Broyden-type matrix used for approximating Jacobians.
Synopsis#
#include "petscksp.h"
PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm comm, PetscInt n, PetscInt N, Mat *B)
Collective
Input Parameters#
comm - MPI communicator
n - number of local rows for storage vectors
N - global size of the storage vectors
Output Parameter#
B - the matrix
Options Database Keys#
-mat_lmvm_hist_size - the number of history vectors to keep
-mat_lmvm_psi - convex ratio between BFGS and DFP components of the update
-mat_lmvm_scale_type - type of scaling applied to J0 (none, scalar, diagonal)
-mat_lmvm_mult_algorithm - the algorithm to use for multiplication (recursive, dense, compact_dense)
-mat_lmvm_cache_J0_products - whether products between the base Jacobian J0 and history vectors should be cached or recomputed
-mat_lmvm_eps - (developer) numerical zero tolerance for testing when an update should be skipped
-mat_lmvm_debug - (developer) perform internal debugging checks
-mat_lmvm_theta - (developer) convex ratio between BFGS and DFP components of the diagonal J0 scaling
-mat_lmvm_rho - (developer) update limiter for the J0 scaling
-mat_lmvm_alpha - (developer) coefficient factor for the quadratic subproblem in J0 scaling
-mat_lmvm_beta - (developer) exponential factor for the diagonal J0 scaling
-mat_lmvm_sigma_hist - (developer) number of past updates to use in J0 scaling
Notes#
L-SymBadBrdn is a convex combination of L-DFP and L-BFGS such that \(B^{-1} = (1 - \psi)*B_{\text{DFP}}^{-1} +
\psi*B_{\text{BFGS}}^{-1}\). The combination factor \(\psi\) is restricted to the range \([0, 1]\), where the L-SymBadBrdn matrix
is guaranteed to be symmetric positive-definite. Note that this combination is on the inverses and not on the
forwards. For forward convex combinations, use the L-SymBrdn matrix (MATLMVMSYMBROYDEN
).
To use the L-SymBrdn matrix with other vector types, the matrix must be created using MatCreate()
and
MatSetType()
, followed by MatLMVMAllocate()
. This ensures that the internal storage and work vectors are
duplicated from the correct type of vector.
It is recommended that one use the MatCreate()
, MatSetType()
and/or MatSetFromOptions()
paradigm instead of this
routine directly.
See Also#
KSP: Linear System Solvers, LMVM Matrices, MatCreate()
, MATLMVM
, MATLMVMSYMBROYDEN
, MatCreateLMVMDFP()
, MatCreateLMVMSR1()
,
MatCreateLMVMBFGS()
, MatCreateLMVMBroyden()
, MatCreateLMVMBadBroyden()
Level#
intermediate
Location#
src/ksp/ksp/utils/lmvm/symbrdn/symbadbrdn.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages