<!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/sdp/dsdpstep.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_2c21778722bc8593635d5ebd154ac3f7.html">sdp</a> </li> </ul> </div> </div> <div class="header"> <div class="headertitle"> <div class="title">dsdpstep.c</div> </div> </div> <div class="contents"> <a href="dsdpstep_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="dsdpdualmat_8h.html" title="The interface between the SDPCone and the matrix S.">dsdpdualmat.h</a>"</span> <a name="l00002"></a>00002 <span class="preprocessor">#include "<a class="code" href="dsdpdsmat_8h.html" title="The interface between the SDPCone and the Delta S matrix.">dsdpdsmat.h</a>"</span> <a name="l00003"></a>00003 <span class="preprocessor">#include "<a class="code" href="dsdpxmat_8h.html" title="The interface between the SDPCone and the dense matrix array.">dsdpxmat.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="dsdplanczos_8h.html" title="Lanczos procedure determines the maximum step length.">dsdplanczos.h</a>"</span> <a name="l00006"></a>00006 <span class="preprocessor">#include "<a class="code" href="dsdplapack_8h.html" title="DSDP uses BLAS and LAPACK for many of its operations.">dsdplapack.h</a>"</span> <a name="l00007"></a>00007 <a name="l00013"></a>00013 <span class="keyword">typedef</span> <span class="keyword">struct </span>_P_Mat3* Mat3; <a name="l00014"></a>00014 <a name="l00015"></a>00015 <span class="keyword">static</span> <span class="keywordtype">int</span> MatMult3(Mat3 A, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> x, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> y); <a name="l00016"></a>00016 <span class="keyword">static</span> <span class="keywordtype">int</span> ComputeStepROBUST(Mat3 A, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> *Q, <span class="keywordtype">int</span> m, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> R,<span class="keywordtype">double</span>*, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> QAQTv, <span class="keywordtype">double</span> *dwork, <span class="keywordtype">double</span> *maxstep, <span class="keywordtype">double</span> *mineig); <a name="l00017"></a>00017 <span class="keyword">static</span> <span class="keywordtype">int</span> ComputeStepFAST(Mat3 A, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> *Q, <span class="keywordtype">int</span> m, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W, <span class="keywordtype">double</span> *dwork, <span class="keywordtype">int</span>*iwork,<span class="keywordtype">double</span> *maxstep, <span class="keywordtype">double</span> *mineig); <a name="l00018"></a>00018 <a name="l00019"></a>00019 <span class="keyword">extern</span> <span class="keywordtype">int</span> DSDPGetEigsSTEP(<span class="keywordtype">double</span>[],<span class="keywordtype">int</span>,<span class="keywordtype">double</span>[],<span class="keywordtype">int</span>,<span class="keywordtype">long</span> <span class="keywordtype">int</span>[],<span class="keywordtype">int</span>, <a name="l00020"></a>00020 <span class="keywordtype">double</span>[],<span class="keywordtype">int</span>,<span class="keywordtype">double</span>[],<span class="keywordtype">int</span>,<span class="keywordtype">int</span>[],<span class="keywordtype">int</span>); <a name="l00021"></a>00021 <a name="l00022"></a>00022 <span class="keywordtype">int</span> DSDPGetTriDiagonalEigs(<span class="keywordtype">int</span> n,<span class="keywordtype">double</span> *D,<span class="keywordtype">double</span> *E, <span class="keywordtype">double</span>*WORK2N,<span class="keywordtype">int</span>*); <a name="l00023"></a>00023 <a name="l00024"></a>00024 <span class="keyword">struct </span>_P_Mat3{ <a name="l00025"></a>00025 <span class="keywordtype">int</span> type; <a name="l00026"></a>00026 <a class="code" href="structDSDPDualMat__C.html" title="Represents an S matrix for one block in the semidefinite cone.">DSDPDualMat</a> ss; <a name="l00027"></a>00027 <a class="code" href="structDSDPDSMat__C.html" title="Symmetric Delta S matrix for one block in the semidefinite cone.">DSDPDSMat</a> ds; <a name="l00028"></a>00028 <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> V; <a name="l00029"></a>00029 <a class="code" href="structDSDPVMat__C.html" title="Dense symmetric matrix for one block in the semidefinite cone.">DSDPVMat</a> x; <a name="l00030"></a>00030 }; <a name="l00031"></a>00031 <a name="l00032"></a>00032 <a name="l00033"></a>00033 <span class="keywordtype">int</span> DSDPGetTriDiagonalEigs(<span class="keywordtype">int</span> N,<span class="keywordtype">double</span> D[],<span class="keywordtype">double</span> E[], <span class="keywordtype">double</span> WORK[],<span class="keywordtype">int</span> IIWORK[]){ <a name="l00034"></a>00034 ffinteger LDZ=DSDPMax(1,N),INFO,NN=N; <a name="l00035"></a>00035 ffinteger M,IL=N-1,IU=N,*ISUPPZ=0; <a name="l00036"></a>00036 ffinteger *IWORK=(ffinteger*)IIWORK,LIWORK,LWORK; <a name="l00037"></a>00037 <span class="keywordtype">double</span> WW[2],VL=-1e10,VU=1e10,*Z=0,ABSTOL=0; <a name="l00038"></a>00038 <span class="keywordtype">char</span> JOBZ=<span class="charliteral">'N'</span>, RANGE=<span class="charliteral">'I'</span>; <a name="l00039"></a>00039 <span class="keywordflow">if</span> (N<50){ <a name="l00040"></a>00040 dstev(&JOBZ,&NN,D,E,Z,&LDZ,WORK,&INFO); <a name="l00041"></a>00041 } <span class="keywordflow">else</span> { <a name="l00042"></a>00042 <a name="l00043"></a>00043 <span class="keywordflow">if</span> (N<=1) IL=1; <a name="l00044"></a>00044 <span class="keywordflow">if</span> (N<=1) IU=1; <a name="l00045"></a>00045 <a name="l00046"></a>00046 LWORK=20*N+1; <a name="l00047"></a>00047 LIWORK=10*N+1; <a name="l00048"></a>00048 <a name="l00049"></a>00049 dstevr(&JOBZ,&RANGE,&NN,D,E,&VL,&VU,&IL,&IU,&ABSTOL,&M,WW,Z,&LDZ,ISUPPZ,WORK,&LWORK,IWORK,&LIWORK,&INFO); <a name="l00050"></a>00050 <a name="l00051"></a>00051 <span class="keywordflow">if</span> (N==1){ <a name="l00052"></a>00052 D[0]=WW[0]; <a name="l00053"></a>00053 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (N>=2){ <a name="l00054"></a>00054 D[N-2]=WW[0]; <a name="l00055"></a>00055 D[N-1]=WW[1]; <a name="l00056"></a>00056 } <a name="l00057"></a>00057 <a name="l00058"></a>00058 } <a name="l00059"></a>00059 <span class="keywordflow">return</span> INFO; <a name="l00060"></a>00060 } <a name="l00061"></a>00061 <a name="l00062"></a>00062 <span class="comment">/* static int id1=0, id2=0; */</span> <a name="l00063"></a>00063 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00064"></a>00064 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "MatMult3"</span> <a name="l00065"></a>00065 <span class="preprocessor"></span><span class="keyword">static</span> <span class="keywordtype">int</span> MatMult3(Mat3 A, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> X, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> Y){ <a name="l00066"></a>00066 <a name="l00067"></a>00067 <span class="keywordtype">int</span> info=0; <a name="l00068"></a>00068 <span class="keywordtype">double</span> minus_one=-1.0; <a name="l00069"></a>00069 <a name="l00070"></a>00070 DSDPFunctionBegin; <a name="l00071"></a>00071 <span class="comment">/* DSDPEventLogBegin(id2); */</span> <a name="l00072"></a>00072 <span class="keywordflow">if</span> (A->type==2){ <a name="l00073"></a>00073 info=<a class="code" href="dsdpxmat_8c.html#aaeb6e1bdb3b88f647f81e2319d614543" title="Multiply X by a vector.">DSDPVMatMult</a>(A->x,X,Y);DSDPCHKERR(info); <a name="l00074"></a>00074 } <span class="keywordflow">else</span> { <a name="l00075"></a>00075 info=<a class="code" href="dsdpdualmat_8c.html#a853fa8c51c83aa5f2e0c6923f76e0ee4" title="Backward triangular solve.">DSDPDualMatCholeskySolveBackward</a>(A->ss,X,Y); DSDPCHKERR(info); <a name="l00076"></a>00076 info=<a class="code" href="dsdpdsmat_8c.html#ad5c54b432eb405cd469a48b885cf47f5" title="Set values into the matrix.">DSDPDSMatMult</a>(A->ds,Y,A->V); DSDPCHKERR(info); <a name="l00077"></a>00077 info=<a class="code" href="dsdpdualmat_8c.html#aa64bc6161dd47646962c7226982f6d62" title="Forward triangular solve.">DSDPDualMatCholeskySolveForward</a>(A->ss,A->V,Y); DSDPCHKERR(info); <a name="l00078"></a>00078 info=<a class="code" href="sdpconevec_8c.html#aa96c836f6ec003af73fd3fc9d24db24e" title="Compute the Euclidean norm.">SDPConeVecScale</a>(minus_one,Y); DSDPCHKERR(info); <a name="l00079"></a>00079 } <a name="l00080"></a>00080 <span class="comment">/* DSDPEventLogEnd(id2);*/</span> <a name="l00081"></a>00081 DSDPFunctionReturn(0); <a name="l00082"></a>00082 } <a name="l00083"></a>00083 <a name="l00084"></a>00084 <a name="l00085"></a>00085 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00086"></a>00086 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPLanczosInitialize"</span> <a name="l00087"></a>00087 <span class="preprocessor"></span> <a name="l00092"></a><a class="code" href="dsdpstep_8c.html#a1186e99d92695e2380ed73eea701c95e">00092</a> <span class="keywordtype">int</span> <a class="code" href="dsdplanczos_8h.html#a5f6b8919f166ba6149cf12e71b8d3bca" title="Initialize Lanczos structure.">DSDPLanczosInitialize</a>( <a class="code" href="structDSDPLanczosStepLength.html" title="Apply Lanczos prodedure to find distance to boundary.">DSDPLanczosStepLength</a> *LZ ){ <a name="l00093"></a>00093 DSDPFunctionBegin; <a name="l00094"></a>00094 <span class="comment">/* Build Lanczos structures */</span> <a name="l00095"></a>00095 LZ->n=0; <a name="l00096"></a>00096 LZ->lanczosm=0; <a name="l00097"></a>00097 LZ->maxlanczosm=20; <a name="l00098"></a>00098 LZ->type=0; <a name="l00099"></a>00099 LZ->dwork4n=0; <a name="l00100"></a>00100 LZ->Q=0; <a name="l00101"></a>00101 LZ->darray=0; <a name="l00102"></a>00102 <span class="comment">/*</span> <a name="l00103"></a>00103 <span class="comment"> if (id1==0 && id2==0){</span> <a name="l00104"></a>00104 <span class="comment"> DSDPEventLogRegister("STEP EIGS",&id1); printf("ID1: %d\n",id1);</span> <a name="l00105"></a>00105 <span class="comment"> DSDPEventLogRegister("STEP MULT",&id2); printf("ID2: %d\n",id2);</span> <a name="l00106"></a>00106 <span class="comment"> }</span> <a name="l00107"></a>00107 <span class="comment"> */</span> <a name="l00108"></a>00108 DSDPFunctionReturn(0); <a name="l00109"></a>00109 } <a name="l00110"></a>00110 <a name="l00111"></a>00111 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00112"></a>00112 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPSetMaximumLanczosIterations"</span> <a name="l00113"></a>00113 <span class="preprocessor"></span> <a name="l00119"></a><a class="code" href="dsdpstep_8c.html#a0a64d376a4373536b179f65b92e7c79a">00119</a> <span class="keywordtype">int</span> <a class="code" href="dsdplanczos_8h.html#a94aaca7c5feb571597a089bc7302a021" title="Set parameter.">DSDPSetMaximumLanczosIterations</a>( <a class="code" href="structDSDPLanczosStepLength.html" title="Apply Lanczos prodedure to find distance to boundary.">DSDPLanczosStepLength</a> *LZ, <span class="keywordtype">int</span> maxlanczos ){ <a name="l00120"></a>00120 DSDPFunctionBegin; <a name="l00121"></a>00121 <span class="keywordflow">if</span> (maxlanczos>0) LZ->lanczosm=maxlanczos; <a name="l00122"></a>00122 DSDPFunctionReturn(0); <a name="l00123"></a>00123 } <a name="l00124"></a>00124 <a name="l00125"></a>00125 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00126"></a>00126 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPFastLanczosSetup"</span> <a name="l00127"></a>00127 <span class="preprocessor"></span> <a name="l00133"></a><a class="code" href="dsdpstep_8c.html#a35e73b100b4385a69d3fb00496dc9dd2">00133</a> <span class="keywordtype">int</span> <a class="code" href="dsdplanczos_8h.html#a0db98a9eb4ae08a45498bf56de769636" title="Use Lanczos procedure. Assume off tridiagonal entries are zero.">DSDPFastLanczosSetup</a>( <a class="code" href="structDSDPLanczosStepLength.html" title="Apply Lanczos prodedure to find distance to boundary.">DSDPLanczosStepLength</a> *LZ, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> V ){ <a name="l00134"></a>00134 <span class="keywordtype">int</span> i,n,info; <a name="l00135"></a>00135 DSDPFunctionBegin; <a name="l00136"></a>00136 <span class="comment">/* Build Lanczos structures */</span> <a name="l00137"></a>00137 info=SDPConeVecGetSize(V,&n);DSDPCHKERR(info); <a name="l00138"></a>00138 LZ->lanczosm=DSDPMin(LZ->maxlanczosm,n); <a name="l00139"></a>00139 LZ->type=1; <a name="l00140"></a>00140 LZ->n=n; <a name="l00141"></a>00141 <span class="keywordflow">if</span> (LZ->lanczosm<50){ <a name="l00142"></a>00142 DSDPCALLOC2(&LZ->dwork4n,<span class="keywordtype">double</span>,4*(LZ->lanczosm)+2,&info); DSDPCHKERR(info); <a name="l00143"></a>00143 DSDPCALLOC2(&LZ->iwork10n,<span class="keywordtype">int</span>,1,&info); DSDPCHKERR(info); <a name="l00144"></a>00144 } <span class="keywordflow">else</span> { <a name="l00145"></a>00145 DSDPCALLOC2(&LZ->dwork4n,<span class="keywordtype">double</span>,23*(LZ->lanczosm)+2,&info); DSDPCHKERR(info); <a name="l00146"></a>00146 DSDPCALLOC2(&LZ->iwork10n,<span class="keywordtype">int</span>,10*(LZ->lanczosm),&info); DSDPCHKERR(info); <a name="l00147"></a>00147 } <a name="l00148"></a>00148 DSDPCALLOC2(&LZ->Q,<a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a>,2,&info);DSDPCHKERR(info); <a name="l00149"></a>00149 <span class="keywordflow">for</span> (i=0;i<2;i++){ <a name="l00150"></a>00150 info = <a class="code" href="sdpconevec_8c.html#a0b5c0b92c582546c95a5d676a3ed06c4" title="Allocate another vector with the same structure as the first.">SDPConeVecDuplicate</a>(V,&LZ->Q[i]);DSDPCHKERR(info); <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__ "DSDPRobustLanczosSetup"</span> <a name="l00157"></a>00157 <span class="preprocessor"></span> <a name="l00163"></a><a class="code" href="dsdpstep_8c.html#a31f54eb95d3afec9b543db50c32673d3">00163</a> <span class="keywordtype">int</span> <a class="code" href="dsdplanczos_8h.html#af016e5b9270f774d081faff7296f0259" title="Use slowerer but more robust method.">DSDPRobustLanczosSetup</a>( <a class="code" href="structDSDPLanczosStepLength.html" title="Apply Lanczos prodedure to find distance to boundary.">DSDPLanczosStepLength</a> *LZ, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> V ){ <a name="l00164"></a>00164 <span class="keywordtype">int</span> i,n,leig,info; <a name="l00165"></a>00165 DSDPFunctionBegin; <a name="l00166"></a>00166 <span class="comment">/* Build Lanczos structures */</span> <a name="l00167"></a>00167 info=SDPConeVecGetSize(V,&n);DSDPCHKERR(info); <a name="l00168"></a>00168 leig=DSDPMin(LZ->maxlanczosm,n); <a name="l00169"></a>00169 LZ->n=n; <a name="l00170"></a>00170 LZ->lanczosm=leig; <a name="l00171"></a>00171 LZ->type=2; <a name="l00172"></a>00172 <a name="l00173"></a>00173 DSDPCALLOC2(&LZ->dwork4n,<span class="keywordtype">double</span>,3*(LZ->lanczosm)+1,&info); DSDPCHKERR(info); <a name="l00174"></a>00174 DSDPCALLOC2(&LZ->darray,<span class="keywordtype">double</span>,(leig*leig),&info); DSDPCHKERR(info); <a name="l00175"></a>00175 DSDPCALLOC2(&LZ->Q,<a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a>,(leig+1),&info);DSDPCHKERR(info); <a name="l00176"></a>00176 <a name="l00177"></a>00177 <span class="keywordflow">for</span> (i=0;i<=leig;i++){ <a name="l00178"></a>00178 info = <a class="code" href="sdpconevec_8c.html#a0b5c0b92c582546c95a5d676a3ed06c4" title="Allocate another vector with the same structure as the first.">SDPConeVecDuplicate</a>(V,&LZ->Q[i]);DSDPCHKERR(info); <a name="l00179"></a>00179 } <a name="l00180"></a>00180 info = SDPConeVecCreate(leig,&LZ->Tv);DSDPCHKERR(info); <a name="l00181"></a>00181 DSDPFunctionReturn(0); <a name="l00182"></a>00182 } <a name="l00183"></a>00183 <a name="l00184"></a>00184 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00185"></a>00185 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPLanczosDestroy"</span> <a name="l00186"></a>00186 <span class="preprocessor"></span> <a name="l00191"></a><a class="code" href="dsdpstep_8c.html#a1a5f06be4c3ba6959de6ae9a76c32a69">00191</a> <span class="keywordtype">int</span> <a class="code" href="dsdplanczos_8h.html#a13450b87544302fa3e6ceda031189a7a" title="Free data structure.">DSDPLanczosDestroy</a>( <a class="code" href="structDSDPLanczosStepLength.html" title="Apply Lanczos prodedure to find distance to boundary.">DSDPLanczosStepLength</a> *LZ){ <a name="l00192"></a>00192 <span class="keywordtype">int</span> i,info; <a name="l00193"></a>00193 DSDPFunctionBegin; <a name="l00194"></a>00194 <span class="keywordflow">if</span> (LZ->type==2){ <a name="l00195"></a>00195 <span class="keywordflow">for</span> (i=0;i<=LZ->lanczosm;i++){ <a name="l00196"></a>00196 info = SDPConeVecDestroy(&LZ->Q[i]);DSDPCHKERR(info); <a name="l00197"></a>00197 } <a name="l00198"></a>00198 info=SDPConeVecDestroy(&LZ->Tv);DSDPCHKERR(info); <a name="l00199"></a>00199 DSDPFREE(&LZ->darray,&info);DSDPCHKERR(info); <a name="l00200"></a>00200 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (LZ->type==1){ <a name="l00201"></a>00201 info = SDPConeVecDestroy(&LZ->Q[1]);DSDPCHKERR(info); <a name="l00202"></a>00202 info = SDPConeVecDestroy(&LZ->Q[0]);DSDPCHKERR(info); <a name="l00203"></a>00203 DSDPFREE(&LZ->iwork10n,&info);DSDPCHKERR(info); <a name="l00204"></a>00204 } <a name="l00205"></a>00205 DSDPFREE(&LZ->Q,&info);DSDPCHKERR(info); <a name="l00206"></a>00206 DSDPFREE(&LZ->dwork4n,&info);DSDPCHKERR(info); <a name="l00207"></a>00207 info=<a class="code" href="dsdplanczos_8h.html#a5f6b8919f166ba6149cf12e71b8d3bca" title="Initialize Lanczos structure.">DSDPLanczosInitialize</a>(LZ);DSDPCHKERR(info); <a name="l00208"></a>00208 DSDPFunctionReturn(0); <a name="l00209"></a>00209 } <a name="l00210"></a>00210 <a name="l00211"></a>00211 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00212"></a>00212 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPLanczosMinXEig"</span> <a name="l00213"></a>00213 <span class="preprocessor"></span><span class="keywordtype">int</span> DSDPLanczosMinXEig( <a class="code" href="structDSDPLanczosStepLength.html" title="Apply Lanczos prodedure to find distance to boundary.">DSDPLanczosStepLength</a> *LZ, <a class="code" href="structDSDPVMat__C.html" title="Dense symmetric matrix for one block in the semidefinite cone.">DSDPVMat</a> X, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W1, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W2, <span class="keywordtype">double</span> *mineig ){ <a name="l00214"></a>00214 <span class="keywordtype">int</span> info,m; <a name="l00215"></a>00215 <span class="keywordtype">double</span> smaxstep; <a name="l00216"></a>00216 <span class="keyword">struct </span>_P_Mat3 PP; <a name="l00217"></a>00217 Mat3 A=&PP; <a name="l00218"></a>00218 <a name="l00219"></a>00219 DSDPFunctionBegin; <a name="l00220"></a>00220 A->type=2; <a name="l00221"></a>00221 A->x=X; <a name="l00222"></a>00222 A->V=W2; <a name="l00223"></a>00223 m=LZ->lanczosm; <a name="l00224"></a>00224 <a name="l00225"></a>00225 <span class="keywordflow">if</span> (LZ->type==1){ <a name="l00226"></a>00226 info = ComputeStepFAST(A,LZ->Q,m,W1,LZ->dwork4n,LZ->iwork10n,&smaxstep,mineig);DSDPCHKERR(info); <a name="l00227"></a>00227 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (LZ->type==2){ <a name="l00228"></a>00228 info = ComputeStepROBUST(A,LZ->Q,m,LZ->Q[m],W1,LZ->darray<span class="comment">/*LZ->TT*/</span>,LZ->Tv,LZ->dwork4n,&smaxstep,mineig);DSDPCHKERR(info); <a name="l00229"></a>00229 } <span class="keywordflow">else</span> { <a name="l00230"></a>00230 DSDPSETERR1(1,<span class="stringliteral">"Lanczos Step Length Has not been SetUp. Type: %d\n"</span>,LZ->type); <a name="l00231"></a>00231 } <a name="l00232"></a>00232 DSDPFunctionReturn(0); <a name="l00233"></a>00233 } <a name="l00234"></a>00234 <a name="l00235"></a>00235 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00236"></a>00236 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "DSDPLanczosStepSize"</span> <a name="l00237"></a>00237 <span class="preprocessor"></span> <a name="l00247"></a><a class="code" href="dsdpstep_8c.html#a3a47fc7937cf0aba8111ac60767a9ba8">00247</a> <span class="keywordtype">int</span> <a class="code" href="dsdplanczos_8h.html#a9f710224114c584d949d0fe7b4e4529a" title="Compute distance to boundary.">DSDPLanczosStepSize</a>( <a class="code" href="structDSDPLanczosStepLength.html" title="Apply Lanczos prodedure to find distance to boundary.">DSDPLanczosStepLength</a> *LZ, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W1, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W2, <a class="code" href="structDSDPDualMat__C.html" title="Represents an S matrix for one block in the semidefinite cone.">DSDPDualMat</a> S, <a class="code" href="structDSDPDSMat__C.html" title="Symmetric Delta S matrix for one block in the semidefinite cone.">DSDPDSMat</a> DS, <span class="keywordtype">double</span> *maxstep ){ <a name="l00248"></a>00248 <span class="keywordtype">int</span> info,m; <a name="l00249"></a>00249 <span class="keywordtype">double</span> smaxstep,mineig; <a name="l00250"></a>00250 <span class="keyword">struct </span>_P_Mat3 PP; <a name="l00251"></a>00251 Mat3 A=&PP; <a name="l00252"></a>00252 <a name="l00253"></a>00253 DSDPFunctionBegin; <a name="l00254"></a>00254 A->ss=S; <a name="l00255"></a>00255 A->ds=DS; A->V=W2; <a name="l00256"></a>00256 A->type=1; <a name="l00257"></a>00257 m=LZ->lanczosm; <a name="l00258"></a>00258 <a name="l00259"></a>00259 <span class="keywordflow">if</span> (LZ->type==1){ <a name="l00260"></a>00260 info = ComputeStepFAST(A,LZ->Q,m,W1,LZ->dwork4n,LZ->iwork10n,&smaxstep,&mineig);DSDPCHKERR(info); <a name="l00261"></a>00261 *maxstep=smaxstep; <a name="l00262"></a>00262 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (LZ->type==2){ <a name="l00263"></a>00263 info = ComputeStepROBUST(A,LZ->Q,m,LZ->Q[m],W1,LZ->darray<span class="comment">/*LZ->TT*/</span>,LZ->Tv,LZ->dwork4n,&smaxstep,&mineig);DSDPCHKERR(info); <a name="l00264"></a>00264 *maxstep=smaxstep; <a name="l00265"></a>00265 } <span class="keywordflow">else</span> { <a name="l00266"></a>00266 DSDPSETERR1(1,<span class="stringliteral">"Lanczos Step Length Has not been SetUp. Type: %d\n"</span>,LZ->type); <a name="l00267"></a>00267 } <a name="l00268"></a>00268 DSDPFunctionReturn(0); <a name="l00269"></a>00269 } <a name="l00270"></a>00270 <a name="l00271"></a>00271 <a name="l00272"></a>00272 <a name="l00273"></a>00273 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00274"></a>00274 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "ComputeStepROBUST"</span> <a name="l00275"></a>00275 <span class="preprocessor"></span><span class="keyword">static</span> <span class="keywordtype">int</span> ComputeStepROBUST(Mat3 A, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> *Q, <span class="keywordtype">int</span> m, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> R, <span class="keywordtype">double</span>*darray, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> QAQTv, <span class="keywordtype">double</span> *dwork, <span class="keywordtype">double</span> *maxstep , <span class="keywordtype">double</span> *mineig){ <a name="l00276"></a>00276 <a name="l00277"></a>00277 <span class="keywordtype">int</span> i,j,n,info; <a name="l00278"></a>00278 <span class="keywordtype">double</span> tt,wnorm, phi; <a name="l00279"></a>00279 <span class="keywordtype">double</span> lambda1=0,lambda2=0,delta=0; <a name="l00280"></a>00280 <span class="keywordtype">double</span> res1,res2,beta; <a name="l00281"></a>00281 <span class="keywordtype">double</span> one=1.0; <a name="l00282"></a>00282 <span class="keywordtype">double</span> *eigvec; <a name="l00283"></a>00283 <span class="keywordtype">int</span> N, LDA, LWORK; <a name="l00284"></a>00284 <a name="l00285"></a>00285 DSDPFunctionBegin; <a name="l00286"></a>00286 <a name="l00287"></a>00287 memset((<span class="keywordtype">void</span>*)darray,0,m*m*<span class="keyword">sizeof</span>(<span class="keywordtype">double</span>)); <a name="l00288"></a>00288 <span class="keywordflow">if</span> (A->type==1){ <a name="l00289"></a>00289 <span class="keywordflow">for</span> (i=0; i<m; i++){ darray[i*m+i]=-1.0;} <a name="l00290"></a>00290 } <span class="keywordflow">else</span> { <a name="l00291"></a>00291 <span class="keywordflow">for</span> (i=0; i<m; i++){ darray[i*m+i]=1.0;} <a name="l00292"></a>00292 } <a name="l00293"></a>00293 <a name="l00294"></a>00294 info = <a class="code" href="sdpconevec_8c.html#aecbde49bc8a9f21a21f8910ee3605753" title="Set each element of vector to this number.">SDPConeVecSet</a>(one,Q[0]);DSDPCHKERR(info); <a name="l00295"></a>00295 info = <a class="code" href="sdpconevec_8c.html#ad4e291595a68e1f476f047a1983eb9d5" title="Scale the vector to norm of 1.">SDPConeVecNormalize</a>(Q[0]);DSDPCHKERR(info); <a name="l00296"></a>00296 <a name="l00297"></a>00297 <span class="keywordflow">for</span> (i=0; i<m; i++){ <a name="l00298"></a>00298 info = MatMult3(A,Q[i],W);DSDPCHKERR(info); <a name="l00299"></a>00299 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(W,&phi);DSDPCHKERR(info); <a name="l00300"></a>00300 <span class="keywordflow">if</span> (phi!=phi){ *maxstep = 0.0; <span class="keywordflow">return</span> 0;} <a name="l00301"></a>00301 <span class="keywordflow">if</span> (i>0){ <a name="l00302"></a>00302 tt=-darray[i*m+i-1]; <a name="l00303"></a>00303 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[i-1],W);DSDPCHKERR(info); <a name="l00304"></a>00304 } <a name="l00305"></a>00305 info = <a class="code" href="sdpconevec_8c.html#aaf2acb6a6221592957980b795972776f" title="Inner product of two vectors.">SDPConeVecDot</a>(W,Q[i],&tt);DSDPCHKERR(info); <a name="l00306"></a>00306 darray[i*m+i]=tt; <a name="l00307"></a>00307 tt*=-1.0; <a name="l00308"></a>00308 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[i],W);DSDPCHKERR(info); <a name="l00309"></a>00309 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(W,&wnorm);DSDPCHKERR(info); <a name="l00310"></a>00310 <span class="keywordflow">if</span> (wnorm <= 0.8 * phi){ <a name="l00311"></a>00311 <span class="keywordflow">for</span> (j=0;j<=i;j++){ <a name="l00312"></a>00312 info = <a class="code" href="sdpconevec_8c.html#aaf2acb6a6221592957980b795972776f" title="Inner product of two vectors.">SDPConeVecDot</a>(W,Q[j],&tt);DSDPCHKERR(info); <a name="l00313"></a>00313 <span class="keywordflow">if</span> (tt==tt){tt*=-1.0;} <span class="keywordflow">else</span> {tt=0;} <a name="l00314"></a>00314 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[j],W);DSDPCHKERR(info); <a name="l00315"></a>00315 darray[j*m+i]-=tt; <a name="l00316"></a>00316 <span class="keywordflow">if</span> (i!=j){ darray[i*m+j]-=tt; } <a name="l00317"></a>00317 } <a name="l00318"></a>00318 } <a name="l00319"></a>00319 <a name="l00320"></a>00320 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(W,&wnorm);DSDPCHKERR(info); <a name="l00321"></a>00321 <span class="keywordflow">if</span> (i<m-1){ <a name="l00322"></a>00322 darray[i*m+i+1]=wnorm; <a name="l00323"></a>00323 darray[i*m+m+i]=wnorm; <a name="l00324"></a>00324 } <a name="l00325"></a>00325 <span class="keywordflow">if</span> (fabs(wnorm)<=1.0e-14) <span class="keywordflow">break</span>; <a name="l00326"></a>00326 info=<a class="code" href="sdpconevec_8c.html#a1e15e8a3438013cc83b22e3d09df5adc" title="Copy v1 to v2.">SDPConeVecCopy</a>(W,Q[i+1]);DSDPCHKERR(info); <a name="l00327"></a>00327 info=<a class="code" href="sdpconevec_8c.html#ad4e291595a68e1f476f047a1983eb9d5" title="Scale the vector to norm of 1.">SDPConeVecNormalize</a>(Q[i+1]); DSDPCHKERR(info); <a name="l00328"></a>00328 <a name="l00329"></a>00329 } <a name="l00330"></a>00330 <span class="comment">/*</span> <a name="l00331"></a>00331 <span class="comment"> DSDPLogInfo("Step Length: Lanczos Iterates: %2.0f, ",i);</span> <a name="l00332"></a>00332 <span class="comment"> DSDPLogInfo("VNorm: %3.1e, ",wnorm);</span> <a name="l00333"></a>00333 <span class="comment"> */</span> <a name="l00334"></a>00334 <span class="comment">/*</span> <a name="l00335"></a>00335 <span class="comment"> printf(" --- TRI DIAGONAL MATRIX ---- \n");</span> <a name="l00336"></a>00336 <span class="comment"> */</span> <a name="l00337"></a>00337 <a name="l00338"></a>00338 <a name="l00339"></a>00339 LWORK=DSDPMax(3*m,1); LDA=DSDPMax(1,m); N=m; <a name="l00340"></a>00340 info=SDPConeVecGetArray(QAQTv,&eigvec);DSDPCHKERR(info); <a name="l00341"></a>00341 info=DSDPGetEigsSTEP(darray,m,0,0,0,0,eigvec,m,dwork,LWORK,0,0);DSDPCHKERR(info); <a name="l00342"></a>00342 info=SDPConeVecRestoreArray(QAQTv,&eigvec);DSDPCHKERR(info); <a name="l00343"></a>00343 <a name="l00344"></a>00344 <span class="keywordflow">if</span> (N==0){ <a name="l00345"></a>00345 lambda1=-0.0; <a name="l00346"></a>00346 delta=1.0e-20; <a name="l00347"></a>00347 *mineig=0; <a name="l00348"></a>00348 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (N==1){ <a name="l00349"></a>00349 info=SDPConeVecGetArray(QAQTv,&eigvec);DSDPCHKERR(info); <a name="l00350"></a>00350 lambda1=-eigvec[0]; <a name="l00351"></a>00351 info=SDPConeVecRestoreArray(QAQTv,&eigvec);DSDPCHKERR(info); <a name="l00352"></a>00352 delta=1.0e-20; <a name="l00353"></a>00353 *mineig=lambda1; <a name="l00354"></a>00354 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (N>1){ <a name="l00355"></a>00355 info=SDPConeVecGetArray(QAQTv,&eigvec);DSDPCHKERR(info); <a name="l00356"></a>00356 *mineig=eigvec[0]; <a name="l00357"></a>00357 lambda1=-eigvec[N-1]; <a name="l00358"></a>00358 lambda2=-eigvec[N-2]; <a name="l00359"></a>00359 info=SDPConeVecRestoreArray(QAQTv,&eigvec);DSDPCHKERR(info); <a name="l00360"></a>00360 <a name="l00361"></a>00361 info = <a class="code" href="sdpconevec_8c.html#aabe353b42ee6b9c2f4db24b3996d2baf" title="Zero the elements of the vector.">SDPConeVecZero</a>(W);DSDPCHKERR(info); <a name="l00362"></a>00362 <span class="keywordflow">for</span> (i=0;i<m;i++){ <a name="l00363"></a>00363 tt=darray[(N-1)*m+i]; <a name="l00364"></a>00364 info=<a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[i],W);DSDPCHKERR(info); <a name="l00365"></a>00365 } <a name="l00366"></a>00366 info = MatMult3(A,W,R);DSDPCHKERR(info); <a name="l00367"></a>00367 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(lambda1,W,R);DSDPCHKERR(info); <a name="l00368"></a>00368 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(R,&res1);DSDPCHKERR(info); <a name="l00369"></a>00369 <a name="l00370"></a>00370 info = <a class="code" href="sdpconevec_8c.html#aabe353b42ee6b9c2f4db24b3996d2baf" title="Zero the elements of the vector.">SDPConeVecZero</a>(W);DSDPCHKERR(info); <a name="l00371"></a>00371 <span class="keywordflow">for</span> (i=0;i<m;i++){ <a name="l00372"></a>00372 tt=darray[(N-2)*m+i]; <a name="l00373"></a>00373 info=<a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[i],W);DSDPCHKERR(info); <a name="l00374"></a>00374 } <a name="l00375"></a>00375 info = MatMult3(A,W,R);DSDPCHKERR(info); <a name="l00376"></a>00376 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(lambda2,W,R);DSDPCHKERR(info); <a name="l00377"></a>00377 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(R,&res2);DSDPCHKERR(info); <a name="l00378"></a>00378 <a name="l00379"></a>00379 tt = -lambda1 + lambda2 - res2; <a name="l00380"></a>00380 <span class="keywordflow">if</span> (tt>0) beta=tt; <a name="l00381"></a>00381 <span class="keywordflow">else</span> beta=1.0e-20; <a name="l00382"></a>00382 delta = DSDPMin(res1,sqrt(res1)/beta); <a name="l00383"></a>00383 <a name="l00384"></a>00384 } <a name="l00385"></a>00385 <a name="l00386"></a>00386 <a name="l00387"></a>00387 <span class="keywordflow">if</span> (delta-lambda1>0) <a name="l00388"></a>00388 *maxstep = 1.0/(delta-lambda1); <a name="l00389"></a>00389 <span class="keywordflow">else</span> <a name="l00390"></a>00390 *maxstep = 1.0e+30; <a name="l00391"></a>00391 <a name="l00392"></a>00392 info=SDPConeVecGetSize(W,&n);DSDPCHKERR(info); <a name="l00393"></a>00393 DSDPLogInfo(0,19,<span class="stringliteral">"Robust Lanczos StepLength: Iterates %d, Max: %d, BlockSize: %d, Lambda1: %4.2e, Res1: %4.2e, Lambda2: %4.2e, Res2: %4.2e, Delta: %4.2e, MaxStep: %4.2e\n"</span>,i,m,n,lambda1,res1*res1,lambda2,res2*res2,delta,*maxstep); <a name="l00394"></a>00394 <a name="l00395"></a>00395 DSDPFunctionReturn(0); <a name="l00396"></a>00396 } <a name="l00397"></a>00397 <a name="l00398"></a>00398 <span class="preprocessor">#undef __FUNCT__ </span> <a name="l00399"></a>00399 <span class="preprocessor"></span><span class="preprocessor">#define __FUNCT__ "ComputeStepFAST"</span> <a name="l00400"></a>00400 <span class="preprocessor"></span><span class="keyword">static</span> <span class="keywordtype">int</span> ComputeStepFAST(Mat3 A, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> *Q, <span class="keywordtype">int</span> m, <a class="code" href="structSDPConeVec__C.html" title="Vector whose length corresponds to dimension of a block in a cone.">SDPConeVec</a> W, <span class="keywordtype">double</span> *dwork, <span class="keywordtype">int</span> *iwork,<span class="keywordtype">double</span> *maxstep ,<span class="keywordtype">double</span> *mineig){ <a name="l00401"></a>00401 <a name="l00402"></a>00402 <span class="keywordtype">int</span> i,j,n,info; <a name="l00403"></a>00403 <span class="keywordtype">double</span> tt,wnorm, phi; <a name="l00404"></a>00404 <span class="keywordtype">double</span> lambda1=0,lambda2=0,delta=0; <a name="l00405"></a>00405 <span class="keywordtype">double</span> res1,res2,beta; <a name="l00406"></a>00406 <span class="keywordtype">double</span> one=1.0; <a name="l00407"></a>00407 <span class="keywordtype">int</span> N=m; <a name="l00408"></a>00408 <span class="keywordtype">double</span> *diag,*subdiag,*ddwork; <a name="l00409"></a>00409 <a name="l00410"></a>00410 DSDPFunctionBegin; <a name="l00411"></a>00411 diag=dwork; <a name="l00412"></a>00412 subdiag=dwork+m; <a name="l00413"></a>00413 ddwork=dwork+2*m; <a name="l00414"></a>00414 <a name="l00415"></a>00415 <span class="keywordflow">if</span> (A->type==1){ <a name="l00416"></a>00416 <span class="keywordflow">for</span> (i=0; i<m; i++){ diag[i]=-1; subdiag[i]=0;} <a name="l00417"></a>00417 } <span class="keywordflow">else</span> { <a name="l00418"></a>00418 <span class="keywordflow">for</span> (i=0; i<m; i++){ diag[i]=1.0; subdiag[i]=0;} <a name="l00419"></a>00419 } <a name="l00420"></a>00420 info = <a class="code" href="sdpconevec_8c.html#aecbde49bc8a9f21a21f8910ee3605753" title="Set each element of vector to this number.">SDPConeVecSet</a>(one,Q[0]);DSDPCHKERR(info); <a name="l00421"></a>00421 info = <a class="code" href="sdpconevec_8c.html#ad4e291595a68e1f476f047a1983eb9d5" title="Scale the vector to norm of 1.">SDPConeVecNormalize</a>(Q[0]);DSDPCHKERR(info); <a name="l00422"></a>00422 <a name="l00423"></a>00423 <span class="keywordflow">for</span> (i=0; i<m; i++){ <a name="l00424"></a>00424 info = MatMult3(A,Q[0],W);DSDPCHKERR(info); <a name="l00425"></a>00425 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(W,&phi);DSDPCHKERR(info); <a name="l00426"></a>00426 <span class="keywordflow">if</span> (phi!=phi){ *maxstep = 0.0; <span class="keywordflow">return</span> 0;} <a name="l00427"></a>00427 <span class="keywordflow">if</span> (i>0){ <a name="l00428"></a>00428 tt=-subdiag[i-1]; <a name="l00429"></a>00429 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[1],W);DSDPCHKERR(info); <a name="l00430"></a>00430 } <a name="l00431"></a>00431 info = <a class="code" href="sdpconevec_8c.html#aaf2acb6a6221592957980b795972776f" title="Inner product of two vectors.">SDPConeVecDot</a>(W,Q[0],&tt);DSDPCHKERR(info); <a name="l00432"></a>00432 diag[i]=tt; <a name="l00433"></a>00433 tt*=-1.0; <a name="l00434"></a>00434 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[0],W);DSDPCHKERR(info); <a name="l00435"></a>00435 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(W,&wnorm);DSDPCHKERR(info); <a name="l00436"></a>00436 <span class="keywordflow">if</span> (wnorm <= 1.0 * phi){ <a name="l00437"></a>00437 <span class="keywordflow">for</span> (j=0;j<=i;j++){ <a name="l00438"></a>00438 <span class="keywordflow">if</span> (j==i-1){ <a name="l00439"></a>00439 info = <a class="code" href="sdpconevec_8c.html#aaf2acb6a6221592957980b795972776f" title="Inner product of two vectors.">SDPConeVecDot</a>(W,Q[1],&tt);DSDPCHKERR(info); <a name="l00440"></a>00440 <span class="keywordflow">if</span> (tt==tt){tt*=-1.0;} <span class="keywordflow">else</span> {tt=0;} <a name="l00441"></a>00441 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[1],W);DSDPCHKERR(info); <a name="l00442"></a>00442 subdiag[i-1]-=tt; <a name="l00443"></a>00443 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (j==i){ <a name="l00444"></a>00444 info = <a class="code" href="sdpconevec_8c.html#aaf2acb6a6221592957980b795972776f" title="Inner product of two vectors.">SDPConeVecDot</a>(W,Q[0],&tt);DSDPCHKERR(info); <a name="l00445"></a>00445 <span class="keywordflow">if</span> (tt==tt){tt*=-1.0;} <span class="keywordflow">else</span> {tt=0;} <a name="l00446"></a>00446 info = <a class="code" href="sdpconevec_8c.html#a09b3931fbac3a916e530016400213d71" title="Add a multiple of X to Y.">SDPConeVecAXPY</a>(tt,Q[0],W);DSDPCHKERR(info); <a name="l00447"></a>00447 diag[i]-=tt; <a name="l00448"></a>00448 } <a name="l00449"></a>00449 <a name="l00450"></a>00450 } <a name="l00451"></a>00451 } <a name="l00452"></a>00452 <a name="l00453"></a>00453 info = <a class="code" href="sdpconevec_8c.html#a0f0dc0f70e8bd7cbdde3dbcf4c925f06" title="Compute the Euclidean norm.">SDPConeVecNorm2</a>(W,&wnorm);DSDPCHKERR(info); <a name="l00454"></a>00454 <span class="comment">/* printf("PHI: %4.4e, VNORM: %4.2e Diag: %4.2e\n",phi,wnorm,diag[i]); */</span> <a name="l00455"></a>00455 <span class="keywordflow">if</span> (i<m-1){ <a name="l00456"></a>00456 subdiag[i]=wnorm; <a name="l00457"></a>00457 } <a name="l00458"></a>00458 <span class="keywordflow">if</span> (fabs(wnorm)<=1.0e-10){i++;<span class="keywordflow">break</span>;} <a name="l00459"></a>00459 info=<a class="code" href="sdpconevec_8c.html#a1e15e8a3438013cc83b22e3d09df5adc" title="Copy v1 to v2.">SDPConeVecCopy</a>(Q[0],Q[1]);DSDPCHKERR(info); <a name="l00460"></a>00460 info=<a class="code" href="sdpconevec_8c.html#a1e15e8a3438013cc83b22e3d09df5adc" title="Copy v1 to v2.">SDPConeVecCopy</a>(W,Q[0]);DSDPCHKERR(info); <a name="l00461"></a>00461 info=<a class="code" href="sdpconevec_8c.html#ad4e291595a68e1f476f047a1983eb9d5" title="Scale the vector to norm of 1.">SDPConeVecNormalize</a>(Q[0]); DSDPCHKERR(info); <a name="l00462"></a>00462 <a name="l00463"></a>00463 } <a name="l00464"></a>00464 <a name="l00465"></a>00465 <span class="comment">/* DSDPEventLogBegin(id1); */</span> <a name="l00466"></a>00466 info=DSDPGetTriDiagonalEigs(m,diag,subdiag,ddwork,iwork); DSDPCHKERR(info); <a name="l00467"></a>00467 <span class="comment">/* DSDPEventLogEnd(id1); */</span> <a name="l00468"></a>00468 <span class="keywordflow">if</span> (N==0){ <a name="l00469"></a>00469 lambda1=-0.0; <a name="l00470"></a>00470 delta=1.0e-20; <a name="l00471"></a>00471 *mineig=0; <a name="l00472"></a>00472 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (N==1){ <a name="l00473"></a>00473 lambda1=-diag[0]; <a name="l00474"></a>00474 delta=1.0e-20; <a name="l00475"></a>00475 *mineig=diag[0]; <a name="l00476"></a>00476 } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (N>1){ <a name="l00477"></a>00477 lambda1=-diag[N-1]; <a name="l00478"></a>00478 lambda2=-diag[N-2]; <a name="l00479"></a>00479 <a name="l00480"></a>00480 res1=1.0e-8; <a name="l00481"></a>00481 res2=1.0e-8; <a name="l00482"></a>00482 <a name="l00483"></a>00483 tt = -lambda1 + lambda2 - res2; <a name="l00484"></a>00484 <span class="keywordflow">if</span> (tt>0) beta=tt; <a name="l00485"></a>00485 <span class="keywordflow">else</span> beta=1.0e-20; <a name="l00486"></a>00486 delta = DSDPMin(res1,sqrt(res1)/beta); <a name="l00487"></a>00487 <a name="l00488"></a>00488 *mineig=diag[0]; <a name="l00489"></a>00489 } <a name="l00490"></a>00490 <a name="l00491"></a>00491 <a name="l00492"></a>00492 <span class="keywordflow">if</span> (delta-lambda1>0) <a name="l00493"></a>00493 *maxstep = 1.0/(delta-lambda1); <a name="l00494"></a>00494 <span class="keywordflow">else</span> <a name="l00495"></a>00495 *maxstep = 1.0e+30; <a name="l00496"></a>00496 <a name="l00497"></a>00497 info=SDPConeVecGetSize(W,&n);DSDPCHKERR(info); <a name="l00498"></a>00498 DSDPLogInfo(0,19,<span class="stringliteral">"Step Length: Fast Lanczos Iterates: %2d, Max: %d, Block Size: %d, VNorm: %3.1e, Lambda1: %4.4e, Lambda2: %4.4e, Delta: %4.2e, Maxstep: %4.2e\n"</span>, <a name="l00499"></a>00499 i,m,n,wnorm,lambda1,lambda2,delta,*maxstep); <a name="l00500"></a>00500 <a name="l00501"></a>00501 DSDPFunctionReturn(0); <a name="l00502"></a>00502 } </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>