<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd"> <html xmlns="http://www.w3.org/1999/xhtml"> <head> <meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/> <title>DSDP: src/solver/dsdpcg.c Source File</title> <link href="tabs.css" rel="stylesheet" type="text/css"/> <link href="doxygen.css" rel="stylesheet" type="text/css"/> </head> <body> <!-- Generated by Doxygen 1.7.4 --> <div id="top"> <div id="titlearea"> <table cellspacing="0" cellpadding="0"> <tbody> <tr style="height: 56px;"> <td style="padding-left: 0.5em;"> <div id="projectname">DSDP</div> </td> </tr> </tbody> </table> </div> <div id="navrow1" class="tabs"> <ul class="tablist"> <li><a href="index.html"><span>Main Page</span></a></li> <li><a href="pages.html"><span>Related Pages</span></a></li> <li><a href="modules.html"><span>Modules</span></a></li> <li><a href="annotated.html"><span>Data Structures</span></a></li> <li class="current"><a href="files.html"><span>Files</span></a></li> <li><a href="dirs.html"><span>Directories</span></a></li> </ul> </div> <div id="navrow2" class="tabs2"> <ul class="tablist"> <li><a href="files.html"><span>File List</span></a></li> <li><a href="globals.html"><span>Globals</span></a></li> </ul> </div> <div id="nav-path" class="navpath"> <ul> <li class="navelem"><a class="el" href="dir_23046874d7fed141927c769a66d8e3a5.html">src</a> </li> <li class="navelem"><a class="el" href="dir_03c647d51c93e018646ff83aa2eeb169.html">solver</a> </li> </ul> </div> </div> <div class="header"> <div class="headertitle"> <div class="title">dsdpcg.c</div> </div> </div> <div class="contents"> <a href="dsdpcg_8c.html">Go to the documentation of this file.</a><div class="fragment"><pre class="fragment"><a name="l00001"></a>00001 <span class="preprocessor">#include "<a class="code" href="dsdpcg_8h.html" title="Internal data structure for CG method.">dsdpcg.h</a>"</span> <a name="l00002"></a>00002 <span class="preprocessor">#include "<a class="code" href="dsdpschurmat_8h.html" title="Methods of a Schur Matrix.">dsdpschurmat.h</a>"</span> <a name="l00003"></a>00003 <span class="preprocessor">#include "<a class="code" href="dsdpvec_8h.html" title="Vector operations used by the solver.">dsdpvec.h</a>"</span> <a name="l00004"></a>00004 <span class="preprocessor">#include "<a class="code" href="dsdpsys_8h.html" title="Error handling, printing, and profiling.">dsdpsys.h</a>"</span> <a name="l00005"></a>00005 <span class="preprocessor">#include "<a class="code" href="dsdp_8h.html" title="Internal data structure for the DSDP solver.">dsdp.h</a>"</span> <a name="l00012"></a>00012 <span class="keyword">typedef</span> <span class="keyword">enum</span> { DSDPNoMatrix=1, DSDPUnfactoredMatrix=2, DSDPFactoredMatrix=3} DSDPCGType; <a name="l00013"></a>00013 <a name="l00014"></a>00014 <span class="keyword">typedef</span> <span class="keyword">struct</span>{ <a name="l00015"></a>00015 DSDPCGType type; <a name="l00016"></a>00016 <a class="code" href="structDSDPSchurMat__C.html" title="Schur complement matrix whose solution is the Newton direction.">DSDPSchurMat</a> M; <a name="l00017"></a>00017 <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> Diag; <a name="l00018"></a>00018 <a class="code" href="structDSDP__C.html" title="Internal structures for the DSDP solver.">DSDP</a> dsdp; <a name="l00019"></a>00019 } DSDPCGMat; <a name="l00020"></a>00020 <a name="l00021"></a>00021 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00022"></a>00022 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGMatMult"</span> <a name="l00023"></a>00023 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPCGMatMult(DSDPCGMat M, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> Y){ <a name="l00024"></a>00024 <span class="keywordtype">int</span> info; <a name="l00025"></a>00025 DSDPFunctionBegin; <a name="l00026"></a>00026 info=DSDPVecZero(Y); DSDPCHKERR(info); <a name="l00027"></a>00027 <span class="keywordflow">if</span> (M.type==DSDPUnfactoredMatrix){ <a name="l00028"></a>00028 info=<a class="code" href="dsdpschurmat_8c.html#a3ac77a7f5f31d116197d0dc07bef782f" title="Multiply M by a vector. y = M x.">DSDPSchurMatMultiply</a>(M.M,X,Y); DSDPCHKERR(info); <a name="l00029"></a>00029 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPFactoredMatrix){ <a name="l00030"></a>00030 info=DSDPSchurMatMultR(M.M,X,Y); DSDPCHKERR(info); <a name="l00031"></a>00031 info=DSDPVecAXPY(-0*M.dsdp->Mshift,X,Y); DSDPCHKERR(info); <a name="l00032"></a>00032 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPNoMatrix){ <a name="l00033"></a>00033 info=<a class="code" href="dsdp_8h.html#a9174f9ff18c7b5ccaa0550f12b22ab5a" title="Add the product of Schur matrix with v.">DSDPHessianMultiplyAdd</a>(M.dsdp,X,Y);DSDPCHKERR(info); <a name="l00034"></a>00034 } <a name="l00035"></a>00035 DSDPFunctionReturn(0); <a name="l00036"></a>00036 } <a name="l00037"></a>00037 <a name="l00038"></a>00038 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00039"></a>00039 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGMatPre"</span> <a name="l00040"></a>00040 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPCGMatPre(DSDPCGMat M, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> Y){ <a name="l00041"></a>00041 <span class="keywordtype">int</span> info; <a name="l00042"></a>00042 DSDPFunctionBegin; <a name="l00043"></a>00043 info=DSDPVecZero(Y); DSDPCHKERR(info); <a name="l00044"></a>00044 <span class="keywordflow">if</span> (M.type==DSDPUnfactoredMatrix){ <a name="l00045"></a>00045 info=DSDPVecPointwiseMult(X,M.Diag,Y);DSDPCHKERR(info); <a name="l00046"></a>00046 info=DSDPVecPointwiseMult(Y,M.Diag,Y);DSDPCHKERR(info); <a name="l00047"></a>00047 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPFactoredMatrix){ <a name="l00048"></a>00048 info=<a class="code" href="dsdpschurmat_8c.html#a4d6f4dcdea100c12545aac3890eb0b9b" title="Solve the linear system.">DSDPSchurMatSolve</a>(M.M,X,Y);DSDPCHKERR(info); <a name="l00049"></a>00049 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPNoMatrix){ <a name="l00050"></a>00050 info=DSDPVecCopy(X,Y);DSDPCHKERR(info); <a name="l00051"></a>00051 } <a name="l00052"></a>00052 DSDPFunctionReturn(0); <a name="l00053"></a>00053 } <a name="l00054"></a>00054 <a name="l00055"></a>00055 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00056"></a>00056 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGMatPreLeft"</span> <a name="l00057"></a>00057 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPCGMatPreLeft(DSDPCGMat M, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> Y){ <a name="l00058"></a>00058 <span class="keywordtype">int</span> info; <a name="l00059"></a>00059 DSDPFunctionBegin; <a name="l00060"></a>00060 info=DSDPVecZero(Y); DSDPCHKERR(info); <a name="l00061"></a>00061 <span class="keywordflow">if</span> (M.type==DSDPUnfactoredMatrix){ <a name="l00062"></a>00062 info=DSDPVecPointwiseMult(X,M.Diag,Y);DSDPCHKERR(info); <a name="l00063"></a>00063 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPFactoredMatrix){ <a name="l00064"></a>00064 info=<a class="code" href="dsdpschurmat_8c.html#a4d6f4dcdea100c12545aac3890eb0b9b" title="Solve the linear system.">DSDPSchurMatSolve</a>(M.M,X,Y);DSDPCHKERR(info); <a name="l00065"></a>00065 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPNoMatrix){ <a name="l00066"></a>00066 info=DSDPVecCopy(X,Y);DSDPCHKERR(info); <a name="l00067"></a>00067 } <a name="l00068"></a>00068 DSDPFunctionReturn(0); <a name="l00069"></a>00069 } <a name="l00070"></a>00070 <a name="l00071"></a>00071 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00072"></a>00072 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGMatPreRight"</span> <a name="l00073"></a>00073 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPCGMatPreRight(DSDPCGMat M, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> Y){ <a name="l00074"></a>00074 <span class="keywordtype">int</span> info; <a name="l00075"></a>00075 DSDPFunctionBegin; <a name="l00076"></a>00076 info=DSDPVecZero(Y); DSDPCHKERR(info); <a name="l00077"></a>00077 <span class="keywordflow">if</span> (M.type==DSDPNoMatrix){ <a name="l00078"></a>00078 info=DSDPVecPointwiseMult(X,M.Diag,Y);DSDPCHKERR(info); <a name="l00079"></a>00079 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPFactoredMatrix){ <a name="l00080"></a>00080 info=DSDPVecCopy(X,Y);DSDPCHKERR(info); <a name="l00081"></a>00081 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (M.type==DSDPUnfactoredMatrix){ <a name="l00082"></a>00082 info=DSDPVecCopy(X,Y);DSDPCHKERR(info); <a name="l00083"></a>00083 } <a name="l00084"></a>00084 DSDPFunctionReturn(0); <a name="l00085"></a>00085 } <a name="l00086"></a>00086 <a name="l00087"></a>00087 <a name="l00088"></a>00088 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00089"></a>00089 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPConjugateResidual"</span> <a name="l00090"></a>00090 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPConjugateResidual(DSDPCGMat B, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> D, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> R, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> BR, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> P, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> BP, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> TT3, <span class="keywordtype">int</span> maxits, <span class="keywordtype">int</span> *iter){ <a name="l00091"></a>00091 <a name="l00092"></a>00092 <span class="keywordtype">int</span> i,n,info; <a name="l00093"></a>00093 <span class="keywordtype">double</span> zero=0.0,minus_one=-1.0; <a name="l00094"></a>00094 <span class="keywordtype">double</span> alpha,beta,bpbp,rBr,rBrOld; <a name="l00095"></a>00095 <span class="keywordtype">double</span> r0,rerr=1.0e20; <a name="l00096"></a>00096 <a name="l00097"></a>00097 DSDPFunctionBegin; <a name="l00098"></a>00098 info=DSDPVecNorm2(X,&rBr); DSDPCHKERR(info); <a name="l00099"></a>00099 <span class="keywordflow">if</span> (rBr>0){ <a name="l00100"></a>00100 info=DSDPVecCopy(X,P); DSDPCHKERR(info); <a name="l00101"></a>00101 info=DSDPCGMatPreRight(B,P,X);DSDPCHKERR(info); <a name="l00102"></a>00102 info=DSDPCGMatMult(B,X,R); DSDPCHKERR(info); <a name="l00103"></a>00103 } <span class="keywordflow">else</span> { <a name="l00104"></a>00104 info=DSDPVecSet(zero,R); DSDPCHKERR(info); <a name="l00105"></a>00105 } <a name="l00106"></a>00106 info=DSDPVecAYPX(minus_one,D,R); DSDPCHKERR(info); <a name="l00107"></a>00107 <a name="l00108"></a>00108 info=DSDPCGMatPreLeft(B,D,R);DSDPCHKERR(info); <a name="l00109"></a>00109 info=DSDPVecCopy(R,P); DSDPCHKERR(info); <a name="l00110"></a>00110 <a name="l00111"></a>00111 info=DSDPCGMatPreRight(B,R,BR);DSDPCHKERR(info); <a name="l00112"></a>00112 info=DSDPCGMatMult(B,BR,TT3); DSDPCHKERR(info); <a name="l00113"></a>00113 info=DSDPCGMatPreRight(B,TT3,BR);DSDPCHKERR(info); <a name="l00114"></a>00114 <a name="l00115"></a>00115 info=DSDPVecCopy(BR,BP); DSDPCHKERR(info); <a name="l00116"></a>00116 info=DSDPVecDot(BR,R,&rBr); DSDPCHKERR(info); <a name="l00117"></a>00117 info=DSDPVecGetSize(X,&n); DSDPCHKERR(info); <a name="l00118"></a>00118 r0=rBr; <a name="l00119"></a>00119 <a name="l00120"></a>00120 <span class="keywordflow">for</span> (i=0;i<maxits;i++){ <a name="l00121"></a>00121 <a name="l00122"></a>00122 <span class="keywordflow">if</span> (rerr/n < 1.0e-30 || rBr/n <= 1.0e-30 || rBr< 1.0e-12 * r0 ) <span class="keywordflow">break</span>; <a name="l00123"></a>00123 <a name="l00124"></a>00124 info=DSDPVecDot(BP,BP,&bpbp); DSDPCHKERR(info); <a name="l00125"></a>00125 alpha = rBr / bpbp; <a name="l00126"></a>00126 info= DSDPVecAXPY(alpha,P,X); DSDPCHKERR(info); <a name="l00127"></a>00127 alpha = -alpha; <a name="l00128"></a>00128 info=DSDPVecAXPY(alpha,BP,R); DSDPCHKERR(info); <a name="l00129"></a>00129 <a name="l00130"></a>00130 info=DSDPCGMatPreRight(B,R,BR);DSDPCHKERR(info); <a name="l00131"></a>00131 info=DSDPCGMatMult(B,BR,TT3); DSDPCHKERR(info); <a name="l00132"></a>00132 info=DSDPCGMatPreLeft(B,TT3,BR);DSDPCHKERR(info); <a name="l00133"></a>00133 <a name="l00134"></a>00134 rBrOld=rBr; <a name="l00135"></a>00135 info=DSDPVecNorm2(R,&rerr); DSDPCHKERR(info); <a name="l00136"></a>00136 info=DSDPVecDot(BR,R,&rBr); DSDPCHKERR(info); <a name="l00137"></a>00137 <a name="l00138"></a>00138 DSDPLogInfo(0,11,<span class="stringliteral">"CG: rerr: %4.4e, rBr: %4.4e \n"</span>,rerr,rBr); <a name="l00139"></a>00139 <a name="l00140"></a>00140 beta = rBr/rBrOld; <a name="l00141"></a>00141 info=DSDPVecAYPX(beta,R,P); DSDPCHKERR(info); <a name="l00142"></a>00142 info=DSDPVecAYPX(beta,BR,BP); DSDPCHKERR(info); <a name="l00143"></a>00143 <a name="l00144"></a>00144 } <a name="l00145"></a>00145 info=DSDPVecCopy(X,BR);DSDPCHKERR(info); <a name="l00146"></a>00146 info=DSDPCGMatPreRight(B,BR,X);DSDPCHKERR(info); <a name="l00147"></a>00147 <a name="l00148"></a>00148 DSDPLogInfo(0,2,<span class="stringliteral">"Conjugate Residual, Initial rMr, %4.2e, Final rMr: %4.2e, Iterates: %d\n"</span>,r0,rBr,i); <a name="l00149"></a>00149 <a name="l00150"></a>00150 *iter=i; <a name="l00151"></a>00151 <a name="l00152"></a>00152 DSDPFunctionReturn(0); <a name="l00153"></a>00153 } <a name="l00154"></a>00154 <a name="l00155"></a>00155 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00156"></a>00156 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPConjugateGradient"</span> <a name="l00157"></a>00157 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPConjugateGradient(DSDPCGMat B, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> D, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> R, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> BR, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> P, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> BP, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> Z, <span class="keywordtype">double</span> cgtol, <span class="keywordtype">int</span> maxits, <span class="keywordtype">int</span> *iter){ <a name="l00158"></a>00158 <a name="l00159"></a>00159 <span class="keywordtype">int</span> i,n,info; <a name="l00160"></a>00160 <span class="keywordtype">double</span> alpha,beta=0,bpbp; <a name="l00161"></a>00161 <span class="keywordtype">double</span> r0,rerr=1.0e20; <a name="l00162"></a>00162 <span class="keywordtype">double</span> rz,rzold,rz0; <a name="l00163"></a>00163 <a name="l00164"></a>00164 DSDPFunctionBegin; <a name="l00165"></a>00165 <span class="keywordflow">if</span> (maxits<=0){ <a name="l00166"></a>00166 *iter=0; <a name="l00167"></a>00167 DSDPFunctionReturn(0); <a name="l00168"></a>00168 } <a name="l00169"></a>00169 info=DSDPVecGetSize(X,&n); DSDPCHKERR(info); <a name="l00170"></a>00170 info=DSDPVecNorm2(X,&rerr); DSDPCHKERR(info); <a name="l00171"></a>00171 info=DSDPVecZero(R); DSDPCHKERR(info); <a name="l00172"></a>00172 <span class="keywordflow">if</span> (rerr>0){ <a name="l00173"></a>00173 info=DSDPCGMatMult(B,X,R);DSDPCHKERR(info); <a name="l00174"></a>00174 } <a name="l00175"></a>00175 alpha=-1.0; <a name="l00176"></a>00176 info=DSDPVecAYPX(alpha,D,R); DSDPCHKERR(info); <a name="l00177"></a>00177 info=DSDPVecNorm2(R,&rerr); DSDPCHKERR(info); <a name="l00178"></a>00178 <span class="keywordflow">if</span> (rerr*sqrt(1.0*n)<1e-11 +0*cgtol*cgtol){ <a name="l00179"></a>00179 *iter=1; <a name="l00180"></a>00180 DSDPFunctionReturn(0); <a name="l00181"></a>00181 } <a name="l00182"></a>00182 <a name="l00183"></a>00183 <span class="keywordflow">if</span> (rerr>0){ <a name="l00184"></a>00184 info=DSDPVecCopy(R,P); DSDPCHKERR(info); <a name="l00185"></a>00185 info=DSDPVecSetR(P,0);DSDPCHKERR(info); <a name="l00186"></a>00186 info=DSDPVecNorm2(P,&rerr); DSDPCHKERR(info); <a name="l00187"></a>00187 info=DSDPCGMatPre(B,R,Z);DSDPCHKERR(info); <a name="l00188"></a>00188 } <a name="l00189"></a>00189 <a name="l00190"></a>00190 info=DSDPVecCopy(Z,P); DSDPCHKERR(info); <a name="l00191"></a>00191 info=DSDPVecDot(R,Z,&rz); DSDPCHKERR(info); <a name="l00192"></a>00192 r0=rerr;rz0=rz; <a name="l00193"></a>00193 <a name="l00194"></a>00194 <span class="keywordflow">for</span> (i=0;i<maxits;i++){ <a name="l00195"></a>00195 <a name="l00196"></a>00196 info=DSDPCGMatMult(B,P,BP);DSDPCHKERR(info); <a name="l00197"></a>00197 info=DSDPVecDot(BP,P,&bpbp); DSDPCHKERR(info); <a name="l00198"></a>00198 <span class="keywordflow">if</span> (bpbp==0) <span class="keywordflow">break</span>; <a name="l00199"></a>00199 alpha = rz/bpbp; <a name="l00200"></a>00200 info=DSDPVecAXPY(alpha,P,X); DSDPCHKERR(info); <a name="l00201"></a>00201 info=DSDPVecAXPY(-alpha,BP,R); DSDPCHKERR(info); <a name="l00202"></a>00202 info=DSDPVecNorm2(R,&rerr); DSDPCHKERR(info); <a name="l00203"></a>00203 <a name="l00204"></a>00204 DSDPLogInfo(0,15,<span class="stringliteral">"CG: iter: %d, rerr: %4.4e, alpha: %4.4e, beta=%4.4e, rz: %4.4e \n"</span>,i,rerr,alpha,beta,rz); <a name="l00205"></a>00205 <span class="keywordflow">if</span> (i>1){ <a name="l00206"></a>00206 <span class="keywordflow">if</span> (rerr*sqrt(1.0*n) < cgtol){ <span class="keywordflow">break</span>;} <a name="l00207"></a>00207 <span class="keywordflow">if</span> (rz*n < cgtol*cgtol){ <span class="keywordflow">break</span>;} <a name="l00208"></a>00208 <span class="keywordflow">if</span> (rerr/r0 < cgtol){ <span class="keywordflow">break</span>;} <a name="l00209"></a>00209 } <a name="l00210"></a>00210 <span class="keywordflow">if</span> (rerr<=0){ <span class="keywordflow">break</span>;} <a name="l00211"></a>00211 info=DSDPCGMatPre(B,R,Z);DSDPCHKERR(info); <a name="l00212"></a>00212 rzold=rz; <a name="l00213"></a>00213 info=DSDPVecDot(R,Z,&rz); DSDPCHKERR(info); <a name="l00214"></a>00214 beta=rz/rzold; <a name="l00215"></a>00215 info= DSDPVecAYPX(beta,Z,P); DSDPCHKERR(info); <a name="l00216"></a>00216 } <a name="l00217"></a>00217 DSDPLogInfo(0,2,<span class="stringliteral">"Conjugate Gradient, Initial ||r||: %4.2e, Final ||r||: %4.2e, Initial ||rz||: %4.4e, ||rz||: %4.4e, Iterates: %d\n"</span>,r0,rerr,rz0,rz,i+1); <a name="l00218"></a>00218 <a name="l00219"></a>00219 *iter=i+1; <a name="l00220"></a>00220 <a name="l00221"></a>00221 DSDPFunctionReturn(0); <a name="l00222"></a>00222 } <a name="l00223"></a>00223 <a name="l00237"></a>00237 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00238"></a>00238 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGSolve"</span> <a name="l00239"></a><a class="code" href="dsdpcg_8c.html#a6f1f6f35fc31a5bec5c06b152a4b4105">00239</a> <span class="preprocessor"></span><span class="keywordtype">int</span> <a class="code" href="dsdp_8h.html#a4e27e7414827a3737eb8a22456e0a9db" title="Apply CG to solve for the step directions.">DSDPCGSolve</a>(<a class="code" href="structDSDP__C.html" title="Internal structures for the DSDP solver.">DSDP</a> dsdp, <a class="code" href="structDSDPSchurMat__C.html" title="Schur complement matrix whose solution is the Newton direction.">DSDPSchurMat</a> MM, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> RHS, <a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X,<span class="keywordtype">double</span> cgtol, <a class="code" href="dsdpbasictypes_8h.html#ae667f2bff3ab5ab0bc109fa76dc4ed65" title="Boolean variables.">DSDPTruth</a> *success){ <a name="l00240"></a>00240 <span class="keywordtype">int</span> iter=0,n,info,maxit=10; <a name="l00241"></a>00241 <span class="keywordtype">double</span> dd,ymax; <a name="l00242"></a>00242 DSDPCG *sles=dsdp->sles; <a name="l00243"></a>00243 DSDPCGMat CGM; <a name="l00244"></a>00244 DSDPFunctionBegin; <a name="l00245"></a>00245 <a name="l00246"></a>00246 info=DSDPEventLogBegin(dsdp->cgtime); <a name="l00247"></a>00247 info=DSDPVecZero(X);DSDPCHKERR(info); <a name="l00248"></a>00248 info=DSDPVecGetSize(X,&n); DSDPCHKERR(info); <a name="l00249"></a>00249 *success=<a class="code" href="dsdpbasictypes_8h.html#ae667f2bff3ab5ab0bc109fa76dc4ed65ae6b2e53e51e94a5fa7204ceac78f824f">DSDP_TRUE</a>; <a name="l00250"></a>00250 <span class="keywordflow">if</span> (0){ <a name="l00251"></a>00251 maxit=0; <a name="l00252"></a>00252 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (dsdp->slestype==1){ <a name="l00253"></a>00253 <a name="l00254"></a>00254 CGM.type=DSDPNoMatrix; <a name="l00255"></a>00255 CGM.M=MM; <a name="l00256"></a>00256 CGM.dsdp=dsdp; <a name="l00257"></a>00257 cgtol*=1000; <a name="l00258"></a>00258 maxit=5; <a name="l00259"></a>00259 <a name="l00260"></a>00260 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (dsdp->slestype==2){ <a name="l00261"></a>00261 CGM.type=DSDPUnfactoredMatrix; <a name="l00262"></a>00262 CGM.M=MM; <a name="l00263"></a>00263 CGM.Diag=sles->Diag; <a name="l00264"></a>00264 CGM.dsdp=dsdp; <a name="l00265"></a>00265 cgtol*=100; <a name="l00266"></a>00266 maxit=(int)sqrt(1.0*n)+10; <a name="l00267"></a>00267 <span class="keywordflow">if</span> (maxit>20) maxit=20; <a name="l00268"></a>00268 info=DSDPVecSet(1.0,sles->Diag);DSDPCHKERR(info); <a name="l00269"></a>00269 <span class="comment">/*</span> <a name="l00270"></a>00270 <span class="comment"> info=DSDPSchurMatGetDiagonal(MM,sles->Diag);DSDPCHKERR(info);</span> <a name="l00271"></a>00271 <span class="comment"> info=DSDPVecReciprocalSqrt(sles->Diag); DSDPCHKERR(info);</span> <a name="l00272"></a>00272 <span class="comment"> */</span> <a name="l00273"></a>00273 <a name="l00274"></a>00274 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (dsdp->slestype==3){ <a name="l00275"></a>00275 CGM.type=DSDPFactoredMatrix; <a name="l00276"></a>00276 CGM.M=MM; <a name="l00277"></a>00277 CGM.dsdp=dsdp; <a name="l00278"></a>00278 maxit=0; <a name="l00279"></a>00279 info=<a class="code" href="group__DSDPSolution.html#gae2e8b2a5e1b71d699d753a045de95add" title="Copy the the infinity norm of the variables y.">DSDPGetMaxYElement</a>(dsdp,&ymax);DSDPCHKERR(info); <a name="l00280"></a>00280 <span class="keywordflow">if</span> (ymax > 1e5 && dsdp->rgap<1e-1) maxit=3; <a name="l00281"></a>00281 <span class="keywordflow">if</span> (0 && ymax > 1e5 && dsdp->rgap<1e-2){ <a name="l00282"></a>00282 maxit=6; <a name="l00283"></a>00283 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (dsdp->rgap<1e-5){ <a name="l00284"></a>00284 maxit=3; <a name="l00285"></a>00285 } <a name="l00286"></a>00286 <a name="l00287"></a>00287 info=<a class="code" href="dsdpschurmat_8c.html#a4d6f4dcdea100c12545aac3890eb0b9b" title="Solve the linear system.">DSDPSchurMatSolve</a>(MM,RHS,X);DSDPCHKERR(info); <a name="l00288"></a>00288 <a name="l00289"></a>00289 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (dsdp->slestype==4){ <a name="l00290"></a>00290 CGM.type=DSDPFactoredMatrix; <a name="l00291"></a>00291 CGM.M=MM; <a name="l00292"></a>00292 CGM.dsdp=dsdp; <a name="l00293"></a>00293 maxit=3; <a name="l00294"></a>00294 info=<a class="code" href="dsdpschurmat_8c.html#a4d6f4dcdea100c12545aac3890eb0b9b" title="Solve the linear system.">DSDPSchurMatSolve</a>(MM,RHS,X);DSDPCHKERR(info); <a name="l00295"></a>00295 } <a name="l00296"></a>00296 <span class="keywordflow">if</span> (n<maxit) maxit=n; <a name="l00297"></a>00297 <a name="l00298"></a>00298 info=DSDPConjugateGradient(CGM,X,RHS, <a name="l00299"></a>00299 sles->R,sles->BR,sles->P,sles->BP, <a name="l00300"></a>00300 sles->TTT,cgtol,maxit,&iter);DSDPCHKERR(info); <a name="l00301"></a>00301 <a name="l00302"></a>00302 <span class="keywordflow">if</span> (iter>=maxit) *success=<a class="code" href="dsdpbasictypes_8h.html#ae667f2bff3ab5ab0bc109fa76dc4ed65a74385569aa7a59059a8847e2d39b754c">DSDP_FALSE</a>; <a name="l00303"></a>00303 info=DSDPVecDot(RHS,X,&dd);DSDPCHKERR(info); <a name="l00304"></a>00304 <span class="keywordflow">if</span> (dd<0) *success=<a class="code" href="dsdpbasictypes_8h.html#ae667f2bff3ab5ab0bc109fa76dc4ed65a74385569aa7a59059a8847e2d39b754c">DSDP_FALSE</a>; <a name="l00305"></a>00305 info=DSDPEventLogEnd(dsdp->cgtime); <a name="l00306"></a>00306 DSDPFunctionReturn(0); <a name="l00307"></a>00307 } <a name="l00308"></a>00308 <a name="l00309"></a>00309 <a name="l00310"></a>00310 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00311"></a>00311 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGInitialize"</span> <a name="l00312"></a>00312 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPCGInitialize(DSDPCG **sles){ <a name="l00313"></a>00313 DSDPCG *ssles; <a name="l00314"></a>00314 <span class="keywordtype">int</span> info; <a name="l00315"></a>00315 DSDPFunctionBegin; <a name="l00316"></a>00316 DSDPCALLOC1(&ssles,DSDPCG,&info);DSDPCHKERR(info); <a name="l00317"></a>00317 ssles->setup2=0; <a name="l00318"></a>00318 *sles=ssles; <a name="l00319"></a>00319 DSDPFunctionReturn(0); <a name="l00320"></a>00320 } <a name="l00321"></a>00321 <a name="l00322"></a>00322 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00323"></a>00323 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGSetup"</span> <a name="l00324"></a>00324 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPCGSetup(DSDPCG *sles,<a class="code" href="dsdpvec_8h.html#a421243d35437ad44ded3c1e34198e8e9" title="This object hold m+2 variables: a scaling of C, the y variables, and r.">DSDPVec</a> X){ <a name="l00325"></a>00325 <span class="keywordtype">int</span> info; <a name="l00326"></a>00326 DSDPFunctionBegin; <a name="l00327"></a>00327 info=DSDPVecGetSize(X,&sles->m);DSDPCHKERR(info); <a name="l00328"></a>00328 <span class="keywordflow">if</span> (sles->setup2==0){ <a name="l00329"></a>00329 info = DSDPVecDuplicate(X,&sles->R); DSDPCHKERR(info); <a name="l00330"></a>00330 info = DSDPVecDuplicate(X,&sles->P); DSDPCHKERR(info); <a name="l00331"></a>00331 info = DSDPVecDuplicate(X,&sles->BP); DSDPCHKERR(info); <a name="l00332"></a>00332 info = DSDPVecDuplicate(X,&sles->BR); DSDPCHKERR(info); <a name="l00333"></a>00333 info = DSDPVecDuplicate(X,&sles->Diag); DSDPCHKERR(info); <a name="l00334"></a>00334 info = DSDPVecDuplicate(X,&sles->TTT); DSDPCHKERR(info); <a name="l00335"></a>00335 } <a name="l00336"></a>00336 sles->setup2=1; <a name="l00337"></a>00337 DSDPFunctionReturn(0); <a name="l00338"></a>00338 } <a name="l00339"></a>00339 <a name="l00340"></a>00340 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00341"></a>00341 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPCGDestroy"</span> <a name="l00342"></a>00342 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPCGDestroy(DSDPCG **ssles){ <a name="l00343"></a>00343 <span class="keywordtype">int</span> info; <a name="l00344"></a>00344 DSDPCG *sles=*ssles; <a name="l00345"></a>00345 DSDPFunctionBegin; <a name="l00346"></a>00346 <span class="keywordflow">if</span> (ssles==0 || sles==0){DSDPFunctionReturn(0);} <a name="l00347"></a>00347 <span class="keywordflow">if</span> (sles->setup2==1){ <a name="l00348"></a>00348 info = DSDPVecDestroy(&sles->R); DSDPCHKERR(info); <a name="l00349"></a>00349 info = DSDPVecDestroy(&sles->P); DSDPCHKERR(info); <a name="l00350"></a>00350 info = DSDPVecDestroy(&sles->BP); DSDPCHKERR(info); <a name="l00351"></a>00351 info = DSDPVecDestroy(&sles->BR); DSDPCHKERR(info); <a name="l00352"></a>00352 info = DSDPVecDestroy(&sles->Diag); DSDPCHKERR(info); <a name="l00353"></a>00353 info = DSDPVecDestroy(&sles->TTT); DSDPCHKERR(info); <a name="l00354"></a>00354 } <a name="l00355"></a>00355 DSDPFREE(ssles,&info);DSDPCHKERR(info); <a name="l00356"></a>00356 DSDPFunctionReturn(0); <a name="l00357"></a>00357 } </pre></div></div> </div> <hr class="footer"/><address class="footer"><small>Generated on Wed Jun 8 2011 for DSDP by  <a href="http://www.doxygen.org/index.html"> <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.7.4 </small></address> </body> </html>