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
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/MatCreateLMVMDiagBroyden.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>MatCreateLMVMDiagBroyden</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.14.5 2021-03-03</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.14.5 v3.14.5 docs/manualpages/KSP/MatCreateLMVMDiagBroyden.html "><small>Report Typos and Errors</small></a></div>
<A NAME="MatCreateLMVMDiagBroyden"><H1>MatCreateLMVMDiagBroyden</H1></A>
DiagBrdn creates a symmetric Broyden-type diagonal matrix used for approximating Hessians. It consists of a convex combination of DFP and BFGS diagonal approximation schemes, such that DiagBrdn = (1-theta)*BFGS + theta*DFP. To preserve symmetric positive-definiteness, we restrict theta to be in [0, 1]. We also ensure positive definiteness by taking the <A HREF="../Vec/VecAbs.html#VecAbs">VecAbs</A>() of the final vector.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscksp.h"
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../KSP/MatCreateLMVMDiagBroyden.html#MatCreateLMVMDiagBroyden">MatCreateLMVMDiagBroyden</A>(<A HREF="../Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A> comm, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> n, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> N, <A HREF="../Mat/Mat.html#Mat">Mat</A> *B)
</PRE>
There are two ways of approximating the diagonal: using the forward (B) update
schemes for BFGS and DFP and then taking the inverse, or directly working with
the inverse (H) update schemes for the BFGS and DFP updates, derived using the
Sherman-Morrison-Woodbury formula. We have implemented both, controlled by a
parameter below.
<P>
In order to use the DiagBrdn matrix with other vector types, i.e. doing MatMults
and <A HREF="../Mat/MatSolves.html#MatSolves">MatSolves</A>, the matrix must first be created using <A HREF="../Mat/MatCreate.html#MatCreate">MatCreate</A>() and <A HREF="../Mat/MatSetType.html#MatSetType">MatSetType</A>(),
followed by <A HREF="../KSP/MatLMVMAllocate.html#MatLMVMAllocate">MatLMVMAllocate</A>(). Then it will be available for updating
(via <A HREF="../KSP/MatLMVMUpdate.html#MatLMVMUpdate">MatLMVMUpdate</A>) in one's favored solver implementation.
This allows for MPI compatibility.
<P>
Collective
<P>
<H3><FONT COLOR="#CC3333">Input Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>comm </B></TD><TD>- MPI communicator, set to <A HREF="../Sys/PETSC_COMM_SELF.html#PETSC_COMM_SELF">PETSC_COMM_SELF</A>
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>n </B></TD><TD>- number of local rows for storage vectors
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>N </B></TD><TD>- global size of the storage vectors
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>B </B></TD><TD>- the matrix
</TD></TR></TABLE>
<P>
It is recommended that one use the <A HREF="../Mat/MatCreate.html#MatCreate">MatCreate</A>(), <A HREF="../Mat/MatSetType.html#MatSetType">MatSetType</A>() and/or <A HREF="../Mat/MatSetFromOptions.html#MatSetFromOptions">MatSetFromOptions</A>()
paradigm instead of this routine directly.
<P>
<H3><FONT COLOR="#CC3333">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-mat_lmvm_theta </B></TD><TD>- (developer) convex ratio between BFGS and DFP components of the diagonal J0 scaling
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-mat_lmvm_rho </B></TD><TD>- (developer) update limiter for the J0 scaling
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-mat_lmvm_alpha </B></TD><TD>- (developer) coefficient factor for the quadratic subproblem in J0 scaling
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-mat_lmvm_beta </B></TD><TD>- (developer) exponential factor for the diagonal J0 scaling
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-mat_lmvm_sigma_hist </B></TD><TD>- (developer) number of past updates to use in J0 scaling.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-mat_lmvm_tol </B></TD><TD>- (developer) tolerance for bounding the denominator of the rescaling away from 0.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-mat_lmvm_forward </B></TD><TD>- (developer) whether or not to use the forward or backward Broyden update to the diagonal
</TD></TR></TABLE>
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../Mat/MatCreate.html#MatCreate">MatCreate</A>(), MATLMVM, MATLMVMDIAGBRDN, <A HREF="../KSP/MatCreateLMVMDFP.html#MatCreateLMVMDFP">MatCreateLMVMDFP</A>(), <A HREF="../KSP/MatCreateLMVMSR1.html#MatCreateLMVMSR1">MatCreateLMVMSR1</A>(),
<BR><A HREF="../KSP/MatCreateLMVMBFGS.html#MatCreateLMVMBFGS">MatCreateLMVMBFGS</A>(), MatCreateLMVMBrdn(), MatCreateLMVMSymBrdn()
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/ksp/utils/lmvm/diagbrdn/diagbrdn.c.html#MatCreateLMVMDiagBroyden">src/ksp/ksp/utils/lmvm/diagbrdn/diagbrdn.c</A>
<BR><A HREF="./index.html">Index of all KSP routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|