petsc-3.10.5 2019-03-28
MatCreateLMVMBadBrdn
Creates a limited-memory modified (aka "bad") Broyden-type approximation matrix used for a Jacobian. L-BadBrdn is not guaranteed to be symmetric or positive-definite.
Synopsis
#include "petscksp.h"
PetscErrorCode MatCreateLMVMBadBrdn(MPI_Comm comm, PetscInt n, PetscInt N, Mat *B)
The provided local and global sizes must match the solution and function vectors
used with MatLMVMUpdate() and MatSolve(). The resulting L-BadBrdn matrix will have
storage vectors allocated with VecCreateSeq() in serial and VecCreateMPI() in
parallel. To use the L-BadBrdn 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.
Collective on MPI_Comm
Input Parameters
| comm | - MPI communicator, set to PETSC_COMM_SELF
|
| n | - number of local rows for storage vectors
|
| N | - global size of the storage vectors
|
Output Parameter
B -the matrix
It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions()
paradigm instead of this routine directly.
Options Database Keys
-mat_lmvm_num_vecs -maximum number of correction vectors (i.e.: updates) stored
-mat_lmvm_scale_type -(developer) type of scaling applied to J0 (none, scalar, diagonal)
-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
See Also
MatCreate(), MATLMVM, MATLMVMBADBRDN, MatCreateLMVMDFP(), MatCreateLMVMSR1(),
MatCreateLMVMBFGS(), MatCreateLMVMBrdn(), MatCreateLMVMSymBrdn()
Level
intermediate
Location
src/ksp/ksp/utils/lmvm/badbrdn/badbrdn.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages