|
123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217221822192220222122222223222422252226222722282229223022312232223322342235223622372238223922402241224222432244224522462247224822492250225122522253225422552256225722582259226022612262226322642265226622672268226922702271227222732274227522762277227822792280228122822283228422852286228722882289229022912292229322942295229622972298229923002301230223032304230523062307230823092310231123122313231423152316231723182319232023212322232323242325232623272328232923302331233223332334233523362337233823392340234123422343234423452346234723482349235023512352235323542355235623572358235923602361236223632364236523662367236823692370237123722373237423752376237723782379238023812382238323842385238623872388238923902391239223932394239523962397239823992400240124022403240424052406240724082409241024112412241324142415241624172418241924202421242224232424242524262427242824292430243124322433243424352436243724382439244024412442244324442445244624472448244924502451245224532454245524562457245824592460246124622463246424652466246724682469247024712472247324742475247624772478247924802481248224832484248524862487248824892490249124922493249424952496249724982499250025012502250325042505250625072508250925102511251225132514251525162517251825192520252125222523252425252526252725282529253025312532253325342535253625372538253925402541254225432544254525462547254825492550255125522553255425552556255725582559256025612562256325642565256625672568256925702571257225732574257525762577257825792580258125822583258425852586258725882589259025912592259325942595259625972598259926002601260226032604260526062607260826092610261126122613261426152616261726182619262026212622262326242625262626272628262926302631263226332634263526362637263826392640264126422643264426452646264726482649265026512652265326542655265626572658265926602661266226632664266526662667266826692670267126722673267426752676267726782679268026812682268326842685268626872688268926902691269226932694269526962697269826992700270127022703270427052706270727082709271027112712271327142715271627172718271927202721272227232724272527262727272827292730273127322733273427352736273727382739274027412742274327442745274627472748274927502751275227532754275527562757275827592760276127622763276427652766276727682769277027712772277327742775277627772778277927802781278227832784278527862787278827892790279127922793279427952796279727982799280028012802280328042805280628072808280928102811281228132814281528162817281828192820282128222823282428252826282728282829283028312832283328342835283628372838283928402841284228432844284528462847284828492850285128522853285428552856285728582859286028612862286328642865286628672868286928702871287228732874287528762877287828792880288128822883288428852886288728882889289028912892289328942895289628972898289929002901290229032904290529062907290829092910291129122913291429152916291729182919292029212922292329242925292629272928292929302931293229332934293529362937293829392940294129422943294429452946294729482949295029512952295329542955295629572958295929602961296229632964296529662967296829692970297129722973297429752976297729782979298029812982298329842985298629872988298929902991299229932994299529962997299829993000300130023003300430053006300730083009301030113012301330143015301630173018301930203021302230233024302530263027302830293030303130323033303430353036303730383039304030413042304330443045304630473048304930503051305230533054305530563057305830593060306130623063306430653066306730683069307030713072307330743075307630773078307930803081308230833084308530863087308830893090309130923093309430953096309730983099310031013102310331043105310631073108310931103111311231133114311531163117311831193120312131223123312431253126312731283129313031313132313331343135313631373138313931403141314231433144314531463147314831493150315131523153315431553156315731583159316031613162316331643165316631673168316931703171317231733174317531763177317831793180318131823183318431853186318731883189319031913192319331943195319631973198319932003201320232033204320532063207320832093210321132123213321432153216321732183219322032213222322332243225322632273228322932303231323232333234323532363237323832393240324132423243324432453246324732483249325032513252325332543255325632573258325932603261326232633264326532663267326832693270327132723273327432753276327732783279328032813282328332843285328632873288328932903291329232933294329532963297329832993300330133023303330433053306330733083309331033113312331333143315331633173318331933203321332233233324332533263327332833293330333133323333333433353336333733383339334033413342334333443345334633473348334933503351335233533354335533563357335833593360336133623363336433653366336733683369337033713372337333743375337633773378337933803381338233833384338533863387338833893390339133923393339433953396339733983399340034013402340334043405340634073408340934103411341234133414341534163417341834193420342134223423342434253426342734283429343034313432343334343435343634373438343934403441344234433444344534463447344834493450345134523453345434553456345734583459346034613462346334643465346634673468346934703471347234733474347534763477347834793480348134823483348434853486348734883489349034913492349334943495349634973498349935003501 |
- *> \brief <b> DGESVD computes the singular value decomposition (SVD) for GE matrices</b>
- *
- * =========== DOCUMENTATION ===========
- *
- * Online html documentation available at
- * http://www.netlib.org/lapack/explore-html/
- *
- *> \htmlonly
- *> Download DGESVD + dependencies
- *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgesvd.f">
- *> [TGZ]</a>
- *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgesvd.f">
- *> [ZIP]</a>
- *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgesvd.f">
- *> [TXT]</a>
- *> \endhtmlonly
- *
- * Definition:
- * ===========
- *
- * SUBROUTINE DGESVD( JOBU, JOBVT, M, N, A, LDA, S, U, LDU, VT, LDVT,
- * WORK, LWORK, INFO )
- *
- * .. Scalar Arguments ..
- * CHARACTER JOBU, JOBVT
- * INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
- * ..
- * .. Array Arguments ..
- * DOUBLE PRECISION A( LDA, * ), S( * ), U( LDU, * ),
- * $ VT( LDVT, * ), WORK( * )
- * ..
- *
- *
- *> \par Purpose:
- * =============
- *>
- *> \verbatim
- *>
- *> DGESVD computes the singular value decomposition (SVD) of a real
- *> M-by-N matrix A, optionally computing the left and/or right singular
- *> vectors. The SVD is written
- *>
- *> A = U * SIGMA * transpose(V)
- *>
- *> where SIGMA is an M-by-N matrix which is zero except for its
- *> min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and
- *> V is an N-by-N orthogonal matrix. The diagonal elements of SIGMA
- *> are the singular values of A; they are real and non-negative, and
- *> are returned in descending order. The first min(m,n) columns of
- *> U and V are the left and right singular vectors of A.
- *>
- *> Note that the routine returns V**T, not V.
- *> \endverbatim
- *
- * Arguments:
- * ==========
- *
- *> \param[in] JOBU
- *> \verbatim
- *> JOBU is CHARACTER*1
- *> Specifies options for computing all or part of the matrix U:
- *> = 'A': all M columns of U are returned in array U:
- *> = 'S': the first min(m,n) columns of U (the left singular
- *> vectors) are returned in the array U;
- *> = 'O': the first min(m,n) columns of U (the left singular
- *> vectors) are overwritten on the array A;
- *> = 'N': no columns of U (no left singular vectors) are
- *> computed.
- *> \endverbatim
- *>
- *> \param[in] JOBVT
- *> \verbatim
- *> JOBVT is CHARACTER*1
- *> Specifies options for computing all or part of the matrix
- *> V**T:
- *> = 'A': all N rows of V**T are returned in the array VT;
- *> = 'S': the first min(m,n) rows of V**T (the right singular
- *> vectors) are returned in the array VT;
- *> = 'O': the first min(m,n) rows of V**T (the right singular
- *> vectors) are overwritten on the array A;
- *> = 'N': no rows of V**T (no right singular vectors) are
- *> computed.
- *>
- *> JOBVT and JOBU cannot both be 'O'.
- *> \endverbatim
- *>
- *> \param[in] M
- *> \verbatim
- *> M is INTEGER
- *> The number of rows of the input matrix A. M >= 0.
- *> \endverbatim
- *>
- *> \param[in] N
- *> \verbatim
- *> N is INTEGER
- *> The number of columns of the input matrix A. N >= 0.
- *> \endverbatim
- *>
- *> \param[in,out] A
- *> \verbatim
- *> A is DOUBLE PRECISION array, dimension (LDA,N)
- *> On entry, the M-by-N matrix A.
- *> On exit,
- *> if JOBU = 'O', A is overwritten with the first min(m,n)
- *> columns of U (the left singular vectors,
- *> stored columnwise);
- *> if JOBVT = 'O', A is overwritten with the first min(m,n)
- *> rows of V**T (the right singular vectors,
- *> stored rowwise);
- *> if JOBU .ne. 'O' and JOBVT .ne. 'O', the contents of A
- *> are destroyed.
- *> \endverbatim
- *>
- *> \param[in] LDA
- *> \verbatim
- *> LDA is INTEGER
- *> The leading dimension of the array A. LDA >= max(1,M).
- *> \endverbatim
- *>
- *> \param[out] S
- *> \verbatim
- *> S is DOUBLE PRECISION array, dimension (min(M,N))
- *> The singular values of A, sorted so that S(i) >= S(i+1).
- *> \endverbatim
- *>
- *> \param[out] U
- *> \verbatim
- *> U is DOUBLE PRECISION array, dimension (LDU,UCOL)
- *> (LDU,M) if JOBU = 'A' or (LDU,min(M,N)) if JOBU = 'S'.
- *> If JOBU = 'A', U contains the M-by-M orthogonal matrix U;
- *> if JOBU = 'S', U contains the first min(m,n) columns of U
- *> (the left singular vectors, stored columnwise);
- *> if JOBU = 'N' or 'O', U is not referenced.
- *> \endverbatim
- *>
- *> \param[in] LDU
- *> \verbatim
- *> LDU is INTEGER
- *> The leading dimension of the array U. LDU >= 1; if
- *> JOBU = 'S' or 'A', LDU >= M.
- *> \endverbatim
- *>
- *> \param[out] VT
- *> \verbatim
- *> VT is DOUBLE PRECISION array, dimension (LDVT,N)
- *> If JOBVT = 'A', VT contains the N-by-N orthogonal matrix
- *> V**T;
- *> if JOBVT = 'S', VT contains the first min(m,n) rows of
- *> V**T (the right singular vectors, stored rowwise);
- *> if JOBVT = 'N' or 'O', VT is not referenced.
- *> \endverbatim
- *>
- *> \param[in] LDVT
- *> \verbatim
- *> LDVT is INTEGER
- *> The leading dimension of the array VT. LDVT >= 1; if
- *> JOBVT = 'A', LDVT >= N; if JOBVT = 'S', LDVT >= min(M,N).
- *> \endverbatim
- *>
- *> \param[out] WORK
- *> \verbatim
- *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
- *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK;
- *> if INFO > 0, WORK(2:MIN(M,N)) contains the unconverged
- *> superdiagonal elements of an upper bidiagonal matrix B
- *> whose diagonal is in S (not necessarily sorted). B
- *> satisfies A = U * B * VT, so it has the same singular values
- *> as A, and singular vectors related by U and VT.
- *> \endverbatim
- *>
- *> \param[in] LWORK
- *> \verbatim
- *> LWORK is INTEGER
- *> The dimension of the array WORK.
- *> LWORK >= MAX(1,5*MIN(M,N)) for the paths (see comments inside code):
- *> - PATH 1 (M much larger than N, JOBU='N')
- *> - PATH 1t (N much larger than M, JOBVT='N')
- *> LWORK >= MAX(1,3*MIN(M,N) + MAX(M,N),5*MIN(M,N)) for the other paths
- *> For good performance, LWORK should generally be larger.
- *>
- *> If LWORK = -1, then a workspace query is assumed; the routine
- *> only calculates the optimal size of the WORK array, returns
- *> this value as the first entry of the WORK array, and no error
- *> message related to LWORK is issued by XERBLA.
- *> \endverbatim
- *>
- *> \param[out] INFO
- *> \verbatim
- *> INFO is INTEGER
- *> = 0: successful exit.
- *> < 0: if INFO = -i, the i-th argument had an illegal value.
- *> > 0: if DBDSQR did not converge, INFO specifies how many
- *> superdiagonals of an intermediate bidiagonal form B
- *> did not converge to zero. See the description of WORK
- *> above for details.
- *> \endverbatim
- *
- * Authors:
- * ========
- *
- *> \author Univ. of Tennessee
- *> \author Univ. of California Berkeley
- *> \author Univ. of Colorado Denver
- *> \author NAG Ltd.
- *
- *> \ingroup doubleGEsing
- *
- * =====================================================================
- SUBROUTINE DGESVD( JOBU, JOBVT, M, N, A, LDA, S, U, LDU,
- $ VT, LDVT, WORK, LWORK, INFO )
- *
- * -- LAPACK driver routine --
- * -- LAPACK is a software package provided by Univ. of Tennessee, --
- * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
- *
- * .. Scalar Arguments ..
- CHARACTER JOBU, JOBVT
- INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
- * ..
- * .. Array Arguments ..
- DOUBLE PRECISION A( LDA, * ), S( * ), U( LDU, * ),
- $ VT( LDVT, * ), WORK( * )
- * ..
- *
- * =====================================================================
- *
- * .. Parameters ..
- DOUBLE PRECISION ZERO, ONE
- PARAMETER ( ZERO = 0.0D0, ONE = 1.0D0 )
- * ..
- * .. Local Scalars ..
- LOGICAL LQUERY, WNTUA, WNTUAS, WNTUN, WNTUO, WNTUS,
- $ WNTVA, WNTVAS, WNTVN, WNTVO, WNTVS
- INTEGER BDSPAC, BLK, CHUNK, I, IE, IERR, IR, ISCL,
- $ ITAU, ITAUP, ITAUQ, IU, IWORK, LDWRKR, LDWRKU,
- $ MAXWRK, MINMN, MINWRK, MNTHR, NCU, NCVT, NRU,
- $ NRVT, WRKBL
- INTEGER LWORK_DGEQRF, LWORK_DORGQR_N, LWORK_DORGQR_M,
- $ LWORK_DGEBRD, LWORK_DORGBR_P, LWORK_DORGBR_Q,
- $ LWORK_DGELQF, LWORK_DORGLQ_N, LWORK_DORGLQ_M
- DOUBLE PRECISION ANRM, BIGNUM, EPS, SMLNUM
- * ..
- * .. Local Arrays ..
- DOUBLE PRECISION DUM( 1 )
- * ..
- * .. External Subroutines ..
- EXTERNAL DBDSQR, DGEBRD, DGELQF, DGEMM, DGEQRF, DLACPY,
- $ DLASCL, DLASET, DORGBR, DORGLQ, DORGQR, DORMBR,
- $ XERBLA
- * ..
- * .. External Functions ..
- LOGICAL LSAME
- INTEGER ILAENV
- DOUBLE PRECISION DLAMCH, DLANGE
- EXTERNAL LSAME, ILAENV, DLAMCH, DLANGE
- * ..
- * .. Intrinsic Functions ..
- INTRINSIC MAX, MIN, SQRT
- * ..
- * .. Executable Statements ..
- *
- * Test the input arguments
- *
- INFO = 0
- MINMN = MIN( M, N )
- WNTUA = LSAME( JOBU, 'A' )
- WNTUS = LSAME( JOBU, 'S' )
- WNTUAS = WNTUA .OR. WNTUS
- WNTUO = LSAME( JOBU, 'O' )
- WNTUN = LSAME( JOBU, 'N' )
- WNTVA = LSAME( JOBVT, 'A' )
- WNTVS = LSAME( JOBVT, 'S' )
- WNTVAS = WNTVA .OR. WNTVS
- WNTVO = LSAME( JOBVT, 'O' )
- WNTVN = LSAME( JOBVT, 'N' )
- LQUERY = ( LWORK.EQ.-1 )
- *
- IF( .NOT.( WNTUA .OR. WNTUS .OR. WNTUO .OR. WNTUN ) ) THEN
- INFO = -1
- ELSE IF( .NOT.( WNTVA .OR. WNTVS .OR. WNTVO .OR. WNTVN ) .OR.
- $ ( WNTVO .AND. WNTUO ) ) THEN
- INFO = -2
- ELSE IF( M.LT.0 ) THEN
- INFO = -3
- ELSE IF( N.LT.0 ) THEN
- INFO = -4
- ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
- INFO = -6
- ELSE IF( LDU.LT.1 .OR. ( WNTUAS .AND. LDU.LT.M ) ) THEN
- INFO = -9
- ELSE IF( LDVT.LT.1 .OR. ( WNTVA .AND. LDVT.LT.N ) .OR.
- $ ( WNTVS .AND. LDVT.LT.MINMN ) ) THEN
- INFO = -11
- END IF
- *
- * Compute workspace
- * (Note: Comments in the code beginning "Workspace:" describe the
- * minimal amount of workspace needed at that point in the code,
- * as well as the preferred amount for good performance.
- * NB refers to the optimal block size for the immediately
- * following subroutine, as returned by ILAENV.)
- *
- IF( INFO.EQ.0 ) THEN
- MINWRK = 1
- MAXWRK = 1
- IF( M.GE.N .AND. MINMN.GT.0 ) THEN
- *
- * Compute space needed for DBDSQR
- *
- MNTHR = ILAENV( 6, 'DGESVD', JOBU // JOBVT, M, N, 0, 0 )
- BDSPAC = 5*N
- * Compute space needed for DGEQRF
- CALL DGEQRF( M, N, A, LDA, DUM(1), DUM(1), -1, IERR )
- LWORK_DGEQRF = INT( DUM(1) )
- * Compute space needed for DORGQR
- CALL DORGQR( M, N, N, A, LDA, DUM(1), DUM(1), -1, IERR )
- LWORK_DORGQR_N = INT( DUM(1) )
- CALL DORGQR( M, M, N, A, LDA, DUM(1), DUM(1), -1, IERR )
- LWORK_DORGQR_M = INT( DUM(1) )
- * Compute space needed for DGEBRD
- CALL DGEBRD( N, N, A, LDA, S, DUM(1), DUM(1),
- $ DUM(1), DUM(1), -1, IERR )
- LWORK_DGEBRD = INT( DUM(1) )
- * Compute space needed for DORGBR P
- CALL DORGBR( 'P', N, N, N, A, LDA, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_P = INT( DUM(1) )
- * Compute space needed for DORGBR Q
- CALL DORGBR( 'Q', N, N, N, A, LDA, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_Q = INT( DUM(1) )
- *
- IF( M.GE.MNTHR ) THEN
- IF( WNTUN ) THEN
- *
- * Path 1 (M much larger than N, JOBU='N')
- *
- MAXWRK = N + LWORK_DGEQRF
- MAXWRK = MAX( MAXWRK, 3*N + LWORK_DGEBRD )
- IF( WNTVO .OR. WNTVAS )
- $ MAXWRK = MAX( MAXWRK, 3*N + LWORK_DORGBR_P )
- MAXWRK = MAX( MAXWRK, BDSPAC )
- MINWRK = MAX( 4*N, BDSPAC )
- ELSE IF( WNTUO .AND. WNTVN ) THEN
- *
- * Path 2 (M much larger than N, JOBU='O', JOBVT='N')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_N )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = MAX( N*N + WRKBL, N*N + M*N + N )
- MINWRK = MAX( 3*N + M, BDSPAC )
- ELSE IF( WNTUO .AND. WNTVAS ) THEN
- *
- * Path 3 (M much larger than N, JOBU='O', JOBVT='S' or
- * 'A')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_N )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = MAX( N*N + WRKBL, N*N + M*N + N )
- MINWRK = MAX( 3*N + M, BDSPAC )
- ELSE IF( WNTUS .AND. WNTVN ) THEN
- *
- * Path 4 (M much larger than N, JOBU='S', JOBVT='N')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_N )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = N*N + WRKBL
- MINWRK = MAX( 3*N + M, BDSPAC )
- ELSE IF( WNTUS .AND. WNTVO ) THEN
- *
- * Path 5 (M much larger than N, JOBU='S', JOBVT='O')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_N )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = 2*N*N + WRKBL
- MINWRK = MAX( 3*N + M, BDSPAC )
- ELSE IF( WNTUS .AND. WNTVAS ) THEN
- *
- * Path 6 (M much larger than N, JOBU='S', JOBVT='S' or
- * 'A')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_N )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = N*N + WRKBL
- MINWRK = MAX( 3*N + M, BDSPAC )
- ELSE IF( WNTUA .AND. WNTVN ) THEN
- *
- * Path 7 (M much larger than N, JOBU='A', JOBVT='N')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_M )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = N*N + WRKBL
- MINWRK = MAX( 3*N + M, BDSPAC )
- ELSE IF( WNTUA .AND. WNTVO ) THEN
- *
- * Path 8 (M much larger than N, JOBU='A', JOBVT='O')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_M )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = 2*N*N + WRKBL
- MINWRK = MAX( 3*N + M, BDSPAC )
- ELSE IF( WNTUA .AND. WNTVAS ) THEN
- *
- * Path 9 (M much larger than N, JOBU='A', JOBVT='S' or
- * 'A')
- *
- WRKBL = N + LWORK_DGEQRF
- WRKBL = MAX( WRKBL, N + LWORK_DORGQR_M )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, 3*N + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = N*N + WRKBL
- MINWRK = MAX( 3*N + M, BDSPAC )
- END IF
- ELSE
- *
- * Path 10 (M at least N, but not much larger)
- *
- CALL DGEBRD( M, N, A, LDA, S, DUM(1), DUM(1),
- $ DUM(1), DUM(1), -1, IERR )
- LWORK_DGEBRD = INT( DUM(1) )
- MAXWRK = 3*N + LWORK_DGEBRD
- IF( WNTUS .OR. WNTUO ) THEN
- CALL DORGBR( 'Q', M, N, N, A, LDA, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_Q = INT( DUM(1) )
- MAXWRK = MAX( MAXWRK, 3*N + LWORK_DORGBR_Q )
- END IF
- IF( WNTUA ) THEN
- CALL DORGBR( 'Q', M, M, N, A, LDA, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_Q = INT( DUM(1) )
- MAXWRK = MAX( MAXWRK, 3*N + LWORK_DORGBR_Q )
- END IF
- IF( .NOT.WNTVN ) THEN
- MAXWRK = MAX( MAXWRK, 3*N + LWORK_DORGBR_P )
- END IF
- MAXWRK = MAX( MAXWRK, BDSPAC )
- MINWRK = MAX( 3*N + M, BDSPAC )
- END IF
- ELSE IF( MINMN.GT.0 ) THEN
- *
- * Compute space needed for DBDSQR
- *
- MNTHR = ILAENV( 6, 'DGESVD', JOBU // JOBVT, M, N, 0, 0 )
- BDSPAC = 5*M
- * Compute space needed for DGELQF
- CALL DGELQF( M, N, A, LDA, DUM(1), DUM(1), -1, IERR )
- LWORK_DGELQF = INT( DUM(1) )
- * Compute space needed for DORGLQ
- CALL DORGLQ( N, N, M, DUM(1), N, DUM(1), DUM(1), -1, IERR )
- LWORK_DORGLQ_N = INT( DUM(1) )
- CALL DORGLQ( M, N, M, A, LDA, DUM(1), DUM(1), -1, IERR )
- LWORK_DORGLQ_M = INT( DUM(1) )
- * Compute space needed for DGEBRD
- CALL DGEBRD( M, M, A, LDA, S, DUM(1), DUM(1),
- $ DUM(1), DUM(1), -1, IERR )
- LWORK_DGEBRD = INT( DUM(1) )
- * Compute space needed for DORGBR P
- CALL DORGBR( 'P', M, M, M, A, N, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_P = INT( DUM(1) )
- * Compute space needed for DORGBR Q
- CALL DORGBR( 'Q', M, M, M, A, N, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_Q = INT( DUM(1) )
- IF( N.GE.MNTHR ) THEN
- IF( WNTVN ) THEN
- *
- * Path 1t(N much larger than M, JOBVT='N')
- *
- MAXWRK = M + LWORK_DGELQF
- MAXWRK = MAX( MAXWRK, 3*M + LWORK_DGEBRD )
- IF( WNTUO .OR. WNTUAS )
- $ MAXWRK = MAX( MAXWRK, 3*M + LWORK_DORGBR_Q )
- MAXWRK = MAX( MAXWRK, BDSPAC )
- MINWRK = MAX( 4*M, BDSPAC )
- ELSE IF( WNTVO .AND. WNTUN ) THEN
- *
- * Path 2t(N much larger than M, JOBU='N', JOBVT='O')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_M )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = MAX( M*M + WRKBL, M*M + M*N + M )
- MINWRK = MAX( 3*M + N, BDSPAC )
- ELSE IF( WNTVO .AND. WNTUAS ) THEN
- *
- * Path 3t(N much larger than M, JOBU='S' or 'A',
- * JOBVT='O')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_M )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = MAX( M*M + WRKBL, M*M + M*N + M )
- MINWRK = MAX( 3*M + N, BDSPAC )
- ELSE IF( WNTVS .AND. WNTUN ) THEN
- *
- * Path 4t(N much larger than M, JOBU='N', JOBVT='S')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_M )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = M*M + WRKBL
- MINWRK = MAX( 3*M + N, BDSPAC )
- ELSE IF( WNTVS .AND. WNTUO ) THEN
- *
- * Path 5t(N much larger than M, JOBU='O', JOBVT='S')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_M )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = 2*M*M + WRKBL
- MINWRK = MAX( 3*M + N, BDSPAC )
- ELSE IF( WNTVS .AND. WNTUAS ) THEN
- *
- * Path 6t(N much larger than M, JOBU='S' or 'A',
- * JOBVT='S')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_M )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = M*M + WRKBL
- MINWRK = MAX( 3*M + N, BDSPAC )
- ELSE IF( WNTVA .AND. WNTUN ) THEN
- *
- * Path 7t(N much larger than M, JOBU='N', JOBVT='A')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_N )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = M*M + WRKBL
- MINWRK = MAX( 3*M + N, BDSPAC )
- ELSE IF( WNTVA .AND. WNTUO ) THEN
- *
- * Path 8t(N much larger than M, JOBU='O', JOBVT='A')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_N )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = 2*M*M + WRKBL
- MINWRK = MAX( 3*M + N, BDSPAC )
- ELSE IF( WNTVA .AND. WNTUAS ) THEN
- *
- * Path 9t(N much larger than M, JOBU='S' or 'A',
- * JOBVT='A')
- *
- WRKBL = M + LWORK_DGELQF
- WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_N )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_P )
- WRKBL = MAX( WRKBL, 3*M + LWORK_DORGBR_Q )
- WRKBL = MAX( WRKBL, BDSPAC )
- MAXWRK = M*M + WRKBL
- MINWRK = MAX( 3*M + N, BDSPAC )
- END IF
- ELSE
- *
- * Path 10t(N greater than M, but not much larger)
- *
- CALL DGEBRD( M, N, A, LDA, S, DUM(1), DUM(1),
- $ DUM(1), DUM(1), -1, IERR )
- LWORK_DGEBRD = INT( DUM(1) )
- MAXWRK = 3*M + LWORK_DGEBRD
- IF( WNTVS .OR. WNTVO ) THEN
- * Compute space needed for DORGBR P
- CALL DORGBR( 'P', M, N, M, A, N, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_P = INT( DUM(1) )
- MAXWRK = MAX( MAXWRK, 3*M + LWORK_DORGBR_P )
- END IF
- IF( WNTVA ) THEN
- CALL DORGBR( 'P', N, N, M, A, N, DUM(1),
- $ DUM(1), -1, IERR )
- LWORK_DORGBR_P = INT( DUM(1) )
- MAXWRK = MAX( MAXWRK, 3*M + LWORK_DORGBR_P )
- END IF
- IF( .NOT.WNTUN ) THEN
- MAXWRK = MAX( MAXWRK, 3*M + LWORK_DORGBR_Q )
- END IF
- MAXWRK = MAX( MAXWRK, BDSPAC )
- MINWRK = MAX( 3*M + N, BDSPAC )
- END IF
- END IF
- MAXWRK = MAX( MAXWRK, MINWRK )
- WORK( 1 ) = MAXWRK
- *
- IF( LWORK.LT.MINWRK .AND. .NOT.LQUERY ) THEN
- INFO = -13
- END IF
- END IF
- *
- IF( INFO.NE.0 ) THEN
- CALL XERBLA( 'DGESVD', -INFO )
- RETURN
- ELSE IF( LQUERY ) THEN
- RETURN
- END IF
- *
- * Quick return if possible
- *
- IF( M.EQ.0 .OR. N.EQ.0 ) THEN
- RETURN
- END IF
- *
- * Get machine constants
- *
- EPS = DLAMCH( 'P' )
- SMLNUM = SQRT( DLAMCH( 'S' ) ) / EPS
- BIGNUM = ONE / SMLNUM
- *
- * Scale A if max element outside range [SMLNUM,BIGNUM]
- *
- ANRM = DLANGE( 'M', M, N, A, LDA, DUM )
- ISCL = 0
- IF( ANRM.GT.ZERO .AND. ANRM.LT.SMLNUM ) THEN
- ISCL = 1
- CALL DLASCL( 'G', 0, 0, ANRM, SMLNUM, M, N, A, LDA, IERR )
- ELSE IF( ANRM.GT.BIGNUM ) THEN
- ISCL = 1
- CALL DLASCL( 'G', 0, 0, ANRM, BIGNUM, M, N, A, LDA, IERR )
- END IF
- *
- IF( M.GE.N ) THEN
- *
- * A has at least as many rows as columns. If A has sufficiently
- * more rows than columns, first reduce using the QR
- * decomposition (if sufficient workspace available)
- *
- IF( M.GE.MNTHR ) THEN
- *
- IF( WNTUN ) THEN
- *
- * Path 1 (M much larger than N, JOBU='N')
- * No left singular vectors to be computed
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Zero out below R
- *
- IF( N .GT. 1 ) THEN
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO, A( 2, 1 ),
- $ LDA )
- END IF
- IE = 1
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in A
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ), LWORK-IWORK+1,
- $ IERR )
- NCVT = 0
- IF( WNTVO .OR. WNTVAS ) THEN
- *
- * If right singular vectors desired, generate P'.
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, A, LDA, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- NCVT = N
- END IF
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing right
- * singular vectors of A in A if desired
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, NCVT, 0, 0, S, WORK( IE ), A, LDA,
- $ DUM, 1, DUM, 1, WORK( IWORK ), INFO )
- *
- * If right singular vectors desired in VT, copy them there
- *
- IF( WNTVAS )
- $ CALL DLACPY( 'F', N, N, A, LDA, VT, LDVT )
- *
- ELSE IF( WNTUO .AND. WNTVN ) THEN
- *
- * Path 2 (M much larger than N, JOBU='O', JOBVT='N')
- * N left singular vectors to be overwritten on A and
- * no right singular vectors to be computed
- *
- IF( LWORK.GE.N*N+MAX( 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.MAX( WRKBL, LDA*N + N ) + LDA*N ) THEN
- *
- * WORK(IU) is LDA by N, WORK(IR) is LDA by N
- *
- LDWRKU = LDA
- LDWRKR = LDA
- ELSE IF( LWORK.GE.MAX( WRKBL, LDA*N + N ) + N*N ) THEN
- *
- * WORK(IU) is LDA by N, WORK(IR) is N by N
- *
- LDWRKU = LDA
- LDWRKR = N
- ELSE
- *
- * WORK(IU) is LDWRKU by N, WORK(IR) is N by N
- *
- LDWRKU = ( LWORK-N*N-N ) / N
- LDWRKR = N
- END IF
- ITAU = IR + LDWRKR*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to WORK(IR) and zero out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, WORK( IR ), LDWRKR )
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO, WORK( IR+1 ),
- $ LDWRKR )
- *
- * Generate Q in A
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in WORK(IR)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, WORK( IR ), LDWRKR, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate left vectors bidiagonalizing R
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IR ), LDWRKR,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IR)
- * (Workspace: need N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, 0, N, 0, S, WORK( IE ), DUM, 1,
- $ WORK( IR ), LDWRKR, DUM, 1,
- $ WORK( IWORK ), INFO )
- IU = IE + N
- *
- * Multiply Q in A by left singular vectors of R in
- * WORK(IR), storing result in WORK(IU) and copying to A
- * (Workspace: need N*N + 2*N, prefer N*N + M*N + N)
- *
- DO 10 I = 1, M, LDWRKU
- CHUNK = MIN( M-I+1, LDWRKU )
- CALL DGEMM( 'N', 'N', CHUNK, N, N, ONE, A( I, 1 ),
- $ LDA, WORK( IR ), LDWRKR, ZERO,
- $ WORK( IU ), LDWRKU )
- CALL DLACPY( 'F', CHUNK, N, WORK( IU ), LDWRKU,
- $ A( I, 1 ), LDA )
- 10 CONTINUE
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- IE = 1
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize A
- * (Workspace: need 3*N + M, prefer 3*N + (M + N)*NB)
- *
- CALL DGEBRD( M, N, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate left vectors bidiagonalizing A
- * (Workspace: need 4*N, prefer 3*N + N*NB)
- *
- CALL DORGBR( 'Q', M, N, N, A, LDA, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in A
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, 0, M, 0, S, WORK( IE ), DUM, 1,
- $ A, LDA, DUM, 1, WORK( IWORK ), INFO )
- *
- END IF
- *
- ELSE IF( WNTUO .AND. WNTVAS ) THEN
- *
- * Path 3 (M much larger than N, JOBU='O', JOBVT='S' or 'A')
- * N left singular vectors to be overwritten on A and
- * N right singular vectors to be computed in VT
- *
- IF( LWORK.GE.N*N+MAX( 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.MAX( WRKBL, LDA*N + N ) + LDA*N ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is LDA by N
- *
- LDWRKU = LDA
- LDWRKR = LDA
- ELSE IF( LWORK.GE.MAX( WRKBL, LDA*N + N ) + N*N ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is N by N
- *
- LDWRKU = LDA
- LDWRKR = N
- ELSE
- *
- * WORK(IU) is LDWRKU by N and WORK(IR) is N by N
- *
- LDWRKU = ( LWORK-N*N-N ) / N
- LDWRKR = N
- END IF
- ITAU = IR + LDWRKR*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to VT, zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, VT, LDVT )
- IF( N.GT.1 )
- $ CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ VT( 2, 1 ), LDVT )
- *
- * Generate Q in A
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in VT, copying result to WORK(IR)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, VT, LDVT, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', N, N, VT, LDVT, WORK( IR ), LDWRKR )
- *
- * Generate left vectors bidiagonalizing R in WORK(IR)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IR ), LDWRKR,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right vectors bidiagonalizing R in VT
- * (Workspace: need N*N + 4*N-1, prefer N*N + 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IR) and computing right
- * singular vectors of R in VT
- * (Workspace: need N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, N, 0, S, WORK( IE ), VT, LDVT,
- $ WORK( IR ), LDWRKR, DUM, 1,
- $ WORK( IWORK ), INFO )
- IU = IE + N
- *
- * Multiply Q in A by left singular vectors of R in
- * WORK(IR), storing result in WORK(IU) and copying to A
- * (Workspace: need N*N + 2*N, prefer N*N + M*N + N)
- *
- DO 20 I = 1, M, LDWRKU
- CHUNK = MIN( M-I+1, LDWRKU )
- CALL DGEMM( 'N', 'N', CHUNK, N, N, ONE, A( I, 1 ),
- $ LDA, WORK( IR ), LDWRKR, ZERO,
- $ WORK( IU ), LDWRKU )
- CALL DLACPY( 'F', CHUNK, N, WORK( IU ), LDWRKU,
- $ A( I, 1 ), LDA )
- 20 CONTINUE
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to VT, zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, VT, LDVT )
- IF( N.GT.1 )
- $ CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ VT( 2, 1 ), LDVT )
- *
- * Generate Q in A
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in VT
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, VT, LDVT, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply Q in A by left vectors bidiagonalizing R
- * (Workspace: need 3*N + M, prefer 3*N + M*NB)
- *
- CALL DORMBR( 'Q', 'R', 'N', M, N, N, VT, LDVT,
- $ WORK( ITAUQ ), A, LDA, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right vectors bidiagonalizing R in VT
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in A and computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, M, 0, S, WORK( IE ), VT, LDVT,
- $ A, LDA, DUM, 1, WORK( IWORK ), INFO )
- *
- END IF
- *
- ELSE IF( WNTUS ) THEN
- *
- IF( WNTVN ) THEN
- *
- * Path 4 (M much larger than N, JOBU='S', JOBVT='N')
- * N left singular vectors to be computed in U and
- * no right singular vectors to be computed
- *
- IF( LWORK.GE.N*N+MAX( 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.WRKBL+LDA*N ) THEN
- *
- * WORK(IR) is LDA by N
- *
- LDWRKR = LDA
- ELSE
- *
- * WORK(IR) is N by N
- *
- LDWRKR = N
- END IF
- ITAU = IR + LDWRKR*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to WORK(IR), zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, WORK( IR ),
- $ LDWRKR )
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ WORK( IR+1 ), LDWRKR )
- *
- * Generate Q in A
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in WORK(IR)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, WORK( IR ), LDWRKR, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left vectors bidiagonalizing R in WORK(IR)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IR ), LDWRKR,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IR)
- * (Workspace: need N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, 0, N, 0, S, WORK( IE ), DUM,
- $ 1, WORK( IR ), LDWRKR, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply Q in A by left singular vectors of R in
- * WORK(IR), storing result in U
- * (Workspace: need N*N)
- *
- CALL DGEMM( 'N', 'N', M, N, N, ONE, A, LDA,
- $ WORK( IR ), LDWRKR, ZERO, U, LDU )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DORGQR( M, N, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Zero out below R in A
- *
- IF( N .GT. 1 ) THEN
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ A( 2, 1 ), LDA )
- END IF
- *
- * Bidiagonalize R in A
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply Q in U by left vectors bidiagonalizing R
- * (Workspace: need 3*N + M, prefer 3*N + M*NB)
- *
- CALL DORMBR( 'Q', 'R', 'N', M, N, N, A, LDA,
- $ WORK( ITAUQ ), U, LDU, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, 0, M, 0, S, WORK( IE ), DUM,
- $ 1, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTVO ) THEN
- *
- * Path 5 (M much larger than N, JOBU='S', JOBVT='O')
- * N left singular vectors to be computed in U and
- * N right singular vectors to be overwritten on A
- *
- IF( LWORK.GE.2*N*N+MAX( 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+2*LDA*N ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is LDA by N
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*N
- LDWRKR = LDA
- ELSE IF( LWORK.GE.WRKBL+( LDA + N )*N ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is N by N
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*N
- LDWRKR = N
- ELSE
- *
- * WORK(IU) is N by N and WORK(IR) is N by N
- *
- LDWRKU = N
- IR = IU + LDWRKU*N
- LDWRKR = N
- END IF
- ITAU = IR + LDWRKR*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R
- * (Workspace: need 2*N*N + 2*N, prefer 2*N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to WORK(IU), zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ WORK( IU+1 ), LDWRKU )
- *
- * Generate Q in A
- * (Workspace: need 2*N*N + 2*N, prefer 2*N*N + N + N*NB)
- *
- CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in WORK(IU), copying result to
- * WORK(IR)
- * (Workspace: need 2*N*N + 4*N,
- * prefer 2*N*N+3*N+2*N*NB)
- *
- CALL DGEBRD( N, N, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', N, N, WORK( IU ), LDWRKU,
- $ WORK( IR ), LDWRKR )
- *
- * Generate left bidiagonalizing vectors in WORK(IU)
- * (Workspace: need 2*N*N + 4*N, prefer 2*N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IU ), LDWRKU,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in WORK(IR)
- * (Workspace: need 2*N*N + 4*N-1,
- * prefer 2*N*N+3*N+(N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, WORK( IR ), LDWRKR,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IU) and computing
- * right singular vectors of R in WORK(IR)
- * (Workspace: need 2*N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, N, 0, S, WORK( IE ),
- $ WORK( IR ), LDWRKR, WORK( IU ),
- $ LDWRKU, DUM, 1, WORK( IWORK ), INFO )
- *
- * Multiply Q in A by left singular vectors of R in
- * WORK(IU), storing result in U
- * (Workspace: need N*N)
- *
- CALL DGEMM( 'N', 'N', M, N, N, ONE, A, LDA,
- $ WORK( IU ), LDWRKU, ZERO, U, LDU )
- *
- * Copy right singular vectors of R to A
- * (Workspace: need N*N)
- *
- CALL DLACPY( 'F', N, N, WORK( IR ), LDWRKR, A,
- $ LDA )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DORGQR( M, N, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Zero out below R in A
- *
- IF( N .GT. 1 ) THEN
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ A( 2, 1 ), LDA )
- END IF
- *
- * Bidiagonalize R in A
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply Q in U by left vectors bidiagonalizing R
- * (Workspace: need 3*N + M, prefer 3*N + M*NB)
- *
- CALL DORMBR( 'Q', 'R', 'N', M, N, N, A, LDA,
- $ WORK( ITAUQ ), U, LDU, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right vectors bidiagonalizing R in A
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, A, LDA, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U and computing right
- * singular vectors of A in A
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, M, 0, S, WORK( IE ), A,
- $ LDA, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTVAS ) THEN
- *
- * Path 6 (M much larger than N, JOBU='S', JOBVT='S'
- * or 'A')
- * N left singular vectors to be computed in U and
- * N right singular vectors to be computed in VT
- *
- IF( LWORK.GE.N*N+MAX( 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+LDA*N ) THEN
- *
- * WORK(IU) is LDA by N
- *
- LDWRKU = LDA
- ELSE
- *
- * WORK(IU) is N by N
- *
- LDWRKU = N
- END IF
- ITAU = IU + LDWRKU*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to WORK(IU), zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ WORK( IU+1 ), LDWRKU )
- *
- * Generate Q in A
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in WORK(IU), copying result to VT
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', N, N, WORK( IU ), LDWRKU, VT,
- $ LDVT )
- *
- * Generate left bidiagonalizing vectors in WORK(IU)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IU ), LDWRKU,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in VT
- * (Workspace: need N*N + 4*N-1,
- * prefer N*N+3*N+(N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IU) and computing
- * right singular vectors of R in VT
- * (Workspace: need N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, N, 0, S, WORK( IE ), VT,
- $ LDVT, WORK( IU ), LDWRKU, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply Q in A by left singular vectors of R in
- * WORK(IU), storing result in U
- * (Workspace: need N*N)
- *
- CALL DGEMM( 'N', 'N', M, N, N, ONE, A, LDA,
- $ WORK( IU ), LDWRKU, ZERO, U, LDU )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DORGQR( M, N, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to VT, zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, VT, LDVT )
- IF( N.GT.1 )
- $ CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ VT( 2, 1 ), LDVT )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in VT
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, VT, LDVT, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply Q in U by left bidiagonalizing vectors
- * in VT
- * (Workspace: need 3*N + M, prefer 3*N + M*NB)
- *
- CALL DORMBR( 'Q', 'R', 'N', M, N, N, VT, LDVT,
- $ WORK( ITAUQ ), U, LDU, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in VT
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U and computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, M, 0, S, WORK( IE ), VT,
- $ LDVT, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- END IF
- *
- ELSE IF( WNTUA ) THEN
- *
- IF( WNTVN ) THEN
- *
- * Path 7 (M much larger than N, JOBU='A', JOBVT='N')
- * M left singular vectors to be computed in U and
- * no right singular vectors to be computed
- *
- IF( LWORK.GE.N*N+MAX( N+M, 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.WRKBL+LDA*N ) THEN
- *
- * WORK(IR) is LDA by N
- *
- LDWRKR = LDA
- ELSE
- *
- * WORK(IR) is N by N
- *
- LDWRKR = N
- END IF
- ITAU = IR + LDWRKR*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Copy R to WORK(IR), zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, WORK( IR ),
- $ LDWRKR )
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ WORK( IR+1 ), LDWRKR )
- *
- * Generate Q in U
- * (Workspace: need N*N + N + M, prefer N*N + N + M*NB)
- *
- CALL DORGQR( M, M, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in WORK(IR)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, WORK( IR ), LDWRKR, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in WORK(IR)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IR ), LDWRKR,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IR)
- * (Workspace: need N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, 0, N, 0, S, WORK( IE ), DUM,
- $ 1, WORK( IR ), LDWRKR, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply Q in U by left singular vectors of R in
- * WORK(IR), storing result in A
- * (Workspace: need N*N)
- *
- CALL DGEMM( 'N', 'N', M, N, N, ONE, U, LDU,
- $ WORK( IR ), LDWRKR, ZERO, A, LDA )
- *
- * Copy left singular vectors of A from A to U
- *
- CALL DLACPY( 'F', M, N, A, LDA, U, LDU )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need N + M, prefer N + M*NB)
- *
- CALL DORGQR( M, M, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Zero out below R in A
- *
- IF( N .GT. 1 ) THEN
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ A( 2, 1 ), LDA )
- END IF
- *
- * Bidiagonalize R in A
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply Q in U by left bidiagonalizing vectors
- * in A
- * (Workspace: need 3*N + M, prefer 3*N + M*NB)
- *
- CALL DORMBR( 'Q', 'R', 'N', M, N, N, A, LDA,
- $ WORK( ITAUQ ), U, LDU, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, 0, M, 0, S, WORK( IE ), DUM,
- $ 1, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTVO ) THEN
- *
- * Path 8 (M much larger than N, JOBU='A', JOBVT='O')
- * M left singular vectors to be computed in U and
- * N right singular vectors to be overwritten on A
- *
- IF( LWORK.GE.2*N*N+MAX( N+M, 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+2*LDA*N ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is LDA by N
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*N
- LDWRKR = LDA
- ELSE IF( LWORK.GE.WRKBL+( LDA + N )*N ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is N by N
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*N
- LDWRKR = N
- ELSE
- *
- * WORK(IU) is N by N and WORK(IR) is N by N
- *
- LDWRKU = N
- IR = IU + LDWRKU*N
- LDWRKR = N
- END IF
- ITAU = IR + LDWRKR*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need 2*N*N + 2*N, prefer 2*N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need 2*N*N + N + M, prefer 2*N*N + N + M*NB)
- *
- CALL DORGQR( M, M, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to WORK(IU), zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ WORK( IU+1 ), LDWRKU )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in WORK(IU), copying result to
- * WORK(IR)
- * (Workspace: need 2*N*N + 4*N,
- * prefer 2*N*N+3*N+2*N*NB)
- *
- CALL DGEBRD( N, N, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', N, N, WORK( IU ), LDWRKU,
- $ WORK( IR ), LDWRKR )
- *
- * Generate left bidiagonalizing vectors in WORK(IU)
- * (Workspace: need 2*N*N + 4*N, prefer 2*N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IU ), LDWRKU,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in WORK(IR)
- * (Workspace: need 2*N*N + 4*N-1,
- * prefer 2*N*N+3*N+(N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, WORK( IR ), LDWRKR,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IU) and computing
- * right singular vectors of R in WORK(IR)
- * (Workspace: need 2*N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, N, 0, S, WORK( IE ),
- $ WORK( IR ), LDWRKR, WORK( IU ),
- $ LDWRKU, DUM, 1, WORK( IWORK ), INFO )
- *
- * Multiply Q in U by left singular vectors of R in
- * WORK(IU), storing result in A
- * (Workspace: need N*N)
- *
- CALL DGEMM( 'N', 'N', M, N, N, ONE, U, LDU,
- $ WORK( IU ), LDWRKU, ZERO, A, LDA )
- *
- * Copy left singular vectors of A from A to U
- *
- CALL DLACPY( 'F', M, N, A, LDA, U, LDU )
- *
- * Copy right singular vectors of R from WORK(IR) to A
- *
- CALL DLACPY( 'F', N, N, WORK( IR ), LDWRKR, A,
- $ LDA )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need N + M, prefer N + M*NB)
- *
- CALL DORGQR( M, M, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Zero out below R in A
- *
- IF( N .GT. 1 ) THEN
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ A( 2, 1 ), LDA )
- END IF
- *
- * Bidiagonalize R in A
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply Q in U by left bidiagonalizing vectors
- * in A
- * (Workspace: need 3*N + M, prefer 3*N + M*NB)
- *
- CALL DORMBR( 'Q', 'R', 'N', M, N, N, A, LDA,
- $ WORK( ITAUQ ), U, LDU, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in A
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, A, LDA, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U and computing right
- * singular vectors of A in A
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, M, 0, S, WORK( IE ), A,
- $ LDA, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTVAS ) THEN
- *
- * Path 9 (M much larger than N, JOBU='A', JOBVT='S'
- * or 'A')
- * M left singular vectors to be computed in U and
- * N right singular vectors to be computed in VT
- *
- IF( LWORK.GE.N*N+MAX( N+M, 4*N, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+LDA*N ) THEN
- *
- * WORK(IU) is LDA by N
- *
- LDWRKU = LDA
- ELSE
- *
- * WORK(IU) is N by N
- *
- LDWRKU = N
- END IF
- ITAU = IU + LDWRKU*N
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need N*N + 2*N, prefer N*N + N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need N*N + N + M, prefer N*N + N + M*NB)
- *
- CALL DORGQR( M, M, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R to WORK(IU), zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ WORK( IU+1 ), LDWRKU )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in WORK(IU), copying result to VT
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', N, N, WORK( IU ), LDWRKU, VT,
- $ LDVT )
- *
- * Generate left bidiagonalizing vectors in WORK(IU)
- * (Workspace: need N*N + 4*N, prefer N*N + 3*N + N*NB)
- *
- CALL DORGBR( 'Q', N, N, N, WORK( IU ), LDWRKU,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in VT
- * (Workspace: need N*N + 4*N-1,
- * prefer N*N+3*N+(N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of R in WORK(IU) and computing
- * right singular vectors of R in VT
- * (Workspace: need N*N + BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, N, 0, S, WORK( IE ), VT,
- $ LDVT, WORK( IU ), LDWRKU, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply Q in U by left singular vectors of R in
- * WORK(IU), storing result in A
- * (Workspace: need N*N)
- *
- CALL DGEMM( 'N', 'N', M, N, N, ONE, U, LDU,
- $ WORK( IU ), LDWRKU, ZERO, A, LDA )
- *
- * Copy left singular vectors of A from A to U
- *
- CALL DLACPY( 'F', M, N, A, LDA, U, LDU )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + N
- *
- * Compute A=Q*R, copying result to U
- * (Workspace: need 2*N, prefer N + N*NB)
- *
- CALL DGEQRF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- *
- * Generate Q in U
- * (Workspace: need N + M, prefer N + M*NB)
- *
- CALL DORGQR( M, M, N, U, LDU, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy R from A to VT, zeroing out below it
- *
- CALL DLACPY( 'U', N, N, A, LDA, VT, LDVT )
- IF( N.GT.1 )
- $ CALL DLASET( 'L', N-1, N-1, ZERO, ZERO,
- $ VT( 2, 1 ), LDVT )
- IE = ITAU
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize R in VT
- * (Workspace: need 4*N, prefer 3*N + 2*N*NB)
- *
- CALL DGEBRD( N, N, VT, LDVT, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply Q in U by left bidiagonalizing vectors
- * in VT
- * (Workspace: need 3*N + M, prefer 3*N + M*NB)
- *
- CALL DORMBR( 'Q', 'R', 'N', M, N, N, VT, LDVT,
- $ WORK( ITAUQ ), U, LDU, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in VT
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + N
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U and computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, N, M, 0, S, WORK( IE ), VT,
- $ LDVT, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- END IF
- *
- END IF
- *
- ELSE
- *
- * M .LT. MNTHR
- *
- * Path 10 (M at least N, but not much larger)
- * Reduce to bidiagonal form without QR decomposition
- *
- IE = 1
- ITAUQ = IE + N
- ITAUP = ITAUQ + N
- IWORK = ITAUP + N
- *
- * Bidiagonalize A
- * (Workspace: need 3*N + M, prefer 3*N + (M + N)*NB)
- *
- CALL DGEBRD( M, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ), LWORK-IWORK+1,
- $ IERR )
- IF( WNTUAS ) THEN
- *
- * If left singular vectors desired in U, copy result to U
- * and generate left bidiagonalizing vectors in U
- * (Workspace: need 3*N + NCU, prefer 3*N + NCU*NB)
- *
- CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
- IF( WNTUS )
- $ NCU = N
- IF( WNTUA )
- $ NCU = M
- CALL DORGBR( 'Q', M, NCU, N, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IF( WNTVAS ) THEN
- *
- * If right singular vectors desired in VT, copy result to
- * VT and generate right bidiagonalizing vectors in VT
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DLACPY( 'U', N, N, A, LDA, VT, LDVT )
- CALL DORGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IF( WNTUO ) THEN
- *
- * If left singular vectors desired in A, generate left
- * bidiagonalizing vectors in A
- * (Workspace: need 4*N, prefer 3*N + N*NB)
- *
- CALL DORGBR( 'Q', M, N, N, A, LDA, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IF( WNTVO ) THEN
- *
- * If right singular vectors desired in A, generate right
- * bidiagonalizing vectors in A
- * (Workspace: need 4*N-1, prefer 3*N + (N-1)*NB)
- *
- CALL DORGBR( 'P', N, N, N, A, LDA, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IWORK = IE + N
- IF( WNTUAS .OR. WNTUO )
- $ NRU = M
- IF( WNTUN )
- $ NRU = 0
- IF( WNTVAS .OR. WNTVO )
- $ NCVT = N
- IF( WNTVN )
- $ NCVT = 0
- IF( ( .NOT.WNTUO ) .AND. ( .NOT.WNTVO ) ) THEN
- *
- * Perform bidiagonal QR iteration, if desired, computing
- * left singular vectors in U and computing right singular
- * vectors in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, NCVT, NRU, 0, S, WORK( IE ), VT,
- $ LDVT, U, LDU, DUM, 1, WORK( IWORK ), INFO )
- ELSE IF( ( .NOT.WNTUO ) .AND. WNTVO ) THEN
- *
- * Perform bidiagonal QR iteration, if desired, computing
- * left singular vectors in U and computing right singular
- * vectors in A
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, NCVT, NRU, 0, S, WORK( IE ), A, LDA,
- $ U, LDU, DUM, 1, WORK( IWORK ), INFO )
- ELSE
- *
- * Perform bidiagonal QR iteration, if desired, computing
- * left singular vectors in A and computing right singular
- * vectors in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', N, NCVT, NRU, 0, S, WORK( IE ), VT,
- $ LDVT, A, LDA, DUM, 1, WORK( IWORK ), INFO )
- END IF
- *
- END IF
- *
- ELSE
- *
- * A has more columns than rows. If A has sufficiently more
- * columns than rows, first reduce using the LQ decomposition (if
- * sufficient workspace available)
- *
- IF( N.GE.MNTHR ) THEN
- *
- IF( WNTVN ) THEN
- *
- * Path 1t(N much larger than M, JOBVT='N')
- * No right singular vectors to be computed
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Zero out above L
- *
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ), LDA )
- IE = 1
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in A
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ), LWORK-IWORK+1,
- $ IERR )
- IF( WNTUO .OR. WNTUAS ) THEN
- *
- * If left singular vectors desired, generate Q
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, A, LDA, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IWORK = IE + M
- NRU = 0
- IF( WNTUO .OR. WNTUAS )
- $ NRU = M
- *
- * Perform bidiagonal QR iteration, computing left singular
- * vectors of A in A if desired
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, 0, NRU, 0, S, WORK( IE ), DUM, 1, A,
- $ LDA, DUM, 1, WORK( IWORK ), INFO )
- *
- * If left singular vectors desired in U, copy them there
- *
- IF( WNTUAS )
- $ CALL DLACPY( 'F', M, M, A, LDA, U, LDU )
- *
- ELSE IF( WNTVO .AND. WNTUN ) THEN
- *
- * Path 2t(N much larger than M, JOBU='N', JOBVT='O')
- * M right singular vectors to be overwritten on A and
- * no left singular vectors to be computed
- *
- IF( LWORK.GE.M*M+MAX( 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.MAX( WRKBL, LDA*N + M ) + LDA*M ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is LDA by M
- *
- LDWRKU = LDA
- CHUNK = N
- LDWRKR = LDA
- ELSE IF( LWORK.GE.MAX( WRKBL, LDA*N + M ) + M*M ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is M by M
- *
- LDWRKU = LDA
- CHUNK = N
- LDWRKR = M
- ELSE
- *
- * WORK(IU) is M by CHUNK and WORK(IR) is M by M
- *
- LDWRKU = M
- CHUNK = ( LWORK-M*M-M ) / M
- LDWRKR = M
- END IF
- ITAU = IR + LDWRKR*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to WORK(IR) and zero out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, WORK( IR ), LDWRKR )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO,
- $ WORK( IR+LDWRKR ), LDWRKR )
- *
- * Generate Q in A
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in WORK(IR)
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, WORK( IR ), LDWRKR, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate right vectors bidiagonalizing L
- * (Workspace: need M*M + 4*M-1, prefer M*M + 3*M + (M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IR ), LDWRKR,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing right
- * singular vectors of L in WORK(IR)
- * (Workspace: need M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, 0, 0, S, WORK( IE ),
- $ WORK( IR ), LDWRKR, DUM, 1, DUM, 1,
- $ WORK( IWORK ), INFO )
- IU = IE + M
- *
- * Multiply right singular vectors of L in WORK(IR) by Q
- * in A, storing result in WORK(IU) and copying to A
- * (Workspace: need M*M + 2*M, prefer M*M + M*N + M)
- *
- DO 30 I = 1, N, CHUNK
- BLK = MIN( N-I+1, CHUNK )
- CALL DGEMM( 'N', 'N', M, BLK, M, ONE, WORK( IR ),
- $ LDWRKR, A( 1, I ), LDA, ZERO,
- $ WORK( IU ), LDWRKU )
- CALL DLACPY( 'F', M, BLK, WORK( IU ), LDWRKU,
- $ A( 1, I ), LDA )
- 30 CONTINUE
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- IE = 1
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize A
- * (Workspace: need 3*M + N, prefer 3*M + (M + N)*NB)
- *
- CALL DGEBRD( M, N, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate right vectors bidiagonalizing A
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'P', M, N, M, A, LDA, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing right
- * singular vectors of A in A
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'L', M, N, 0, 0, S, WORK( IE ), A, LDA,
- $ DUM, 1, DUM, 1, WORK( IWORK ), INFO )
- *
- END IF
- *
- ELSE IF( WNTVO .AND. WNTUAS ) THEN
- *
- * Path 3t(N much larger than M, JOBU='S' or 'A', JOBVT='O')
- * M right singular vectors to be overwritten on A and
- * M left singular vectors to be computed in U
- *
- IF( LWORK.GE.M*M+MAX( 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.MAX( WRKBL, LDA*N + M ) + LDA*M ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is LDA by M
- *
- LDWRKU = LDA
- CHUNK = N
- LDWRKR = LDA
- ELSE IF( LWORK.GE.MAX( WRKBL, LDA*N + M ) + M*M ) THEN
- *
- * WORK(IU) is LDA by N and WORK(IR) is M by M
- *
- LDWRKU = LDA
- CHUNK = N
- LDWRKR = M
- ELSE
- *
- * WORK(IU) is M by CHUNK and WORK(IR) is M by M
- *
- LDWRKU = M
- CHUNK = ( LWORK-M*M-M ) / M
- LDWRKR = M
- END IF
- ITAU = IR + LDWRKR*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to U, zeroing about above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, U, LDU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, U( 1, 2 ),
- $ LDU )
- *
- * Generate Q in A
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in U, copying result to WORK(IR)
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, U, LDU, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, M, U, LDU, WORK( IR ), LDWRKR )
- *
- * Generate right vectors bidiagonalizing L in WORK(IR)
- * (Workspace: need M*M + 4*M-1, prefer M*M + 3*M + (M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IR ), LDWRKR,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left vectors bidiagonalizing L in U
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of L in U, and computing right
- * singular vectors of L in WORK(IR)
- * (Workspace: need M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, M, 0, S, WORK( IE ),
- $ WORK( IR ), LDWRKR, U, LDU, DUM, 1,
- $ WORK( IWORK ), INFO )
- IU = IE + M
- *
- * Multiply right singular vectors of L in WORK(IR) by Q
- * in A, storing result in WORK(IU) and copying to A
- * (Workspace: need M*M + 2*M, prefer M*M + M*N + M))
- *
- DO 40 I = 1, N, CHUNK
- BLK = MIN( N-I+1, CHUNK )
- CALL DGEMM( 'N', 'N', M, BLK, M, ONE, WORK( IR ),
- $ LDWRKR, A( 1, I ), LDA, ZERO,
- $ WORK( IU ), LDWRKU )
- CALL DLACPY( 'F', M, BLK, WORK( IU ), LDWRKU,
- $ A( 1, I ), LDA )
- 40 CONTINUE
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to U, zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, U, LDU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, U( 1, 2 ),
- $ LDU )
- *
- * Generate Q in A
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in U
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, U, LDU, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply right vectors bidiagonalizing L by Q in A
- * (Workspace: need 3*M + N, prefer 3*M + N*NB)
- *
- CALL DORMBR( 'P', 'L', 'T', M, N, M, U, LDU,
- $ WORK( ITAUP ), A, LDA, WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left vectors bidiagonalizing L in U
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U and computing right
- * singular vectors of A in A
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, N, M, 0, S, WORK( IE ), A, LDA,
- $ U, LDU, DUM, 1, WORK( IWORK ), INFO )
- *
- END IF
- *
- ELSE IF( WNTVS ) THEN
- *
- IF( WNTUN ) THEN
- *
- * Path 4t(N much larger than M, JOBU='N', JOBVT='S')
- * M right singular vectors to be computed in VT and
- * no left singular vectors to be computed
- *
- IF( LWORK.GE.M*M+MAX( 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.WRKBL+LDA*M ) THEN
- *
- * WORK(IR) is LDA by M
- *
- LDWRKR = LDA
- ELSE
- *
- * WORK(IR) is M by M
- *
- LDWRKR = M
- END IF
- ITAU = IR + LDWRKR*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to WORK(IR), zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, WORK( IR ),
- $ LDWRKR )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO,
- $ WORK( IR+LDWRKR ), LDWRKR )
- *
- * Generate Q in A
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in WORK(IR)
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, WORK( IR ), LDWRKR, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right vectors bidiagonalizing L in
- * WORK(IR)
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + (M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IR ), LDWRKR,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing right
- * singular vectors of L in WORK(IR)
- * (Workspace: need M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, 0, 0, S, WORK( IE ),
- $ WORK( IR ), LDWRKR, DUM, 1, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply right singular vectors of L in WORK(IR) by
- * Q in A, storing result in VT
- * (Workspace: need M*M)
- *
- CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IR ),
- $ LDWRKR, A, LDA, ZERO, VT, LDVT )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy result to VT
- *
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DORGLQ( M, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Zero out above L in A
- *
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ),
- $ LDA )
- *
- * Bidiagonalize L in A
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply right vectors bidiagonalizing L by Q in VT
- * (Workspace: need 3*M + N, prefer 3*M + N*NB)
- *
- CALL DORMBR( 'P', 'L', 'T', M, N, M, A, LDA,
- $ WORK( ITAUP ), VT, LDVT,
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, N, 0, 0, S, WORK( IE ), VT,
- $ LDVT, DUM, 1, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTUO ) THEN
- *
- * Path 5t(N much larger than M, JOBU='O', JOBVT='S')
- * M right singular vectors to be computed in VT and
- * M left singular vectors to be overwritten on A
- *
- IF( LWORK.GE.2*M*M+MAX( 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+2*LDA*M ) THEN
- *
- * WORK(IU) is LDA by M and WORK(IR) is LDA by M
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*M
- LDWRKR = LDA
- ELSE IF( LWORK.GE.WRKBL+( LDA + M )*M ) THEN
- *
- * WORK(IU) is LDA by M and WORK(IR) is M by M
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*M
- LDWRKR = M
- ELSE
- *
- * WORK(IU) is M by M and WORK(IR) is M by M
- *
- LDWRKU = M
- IR = IU + LDWRKU*M
- LDWRKR = M
- END IF
- ITAU = IR + LDWRKR*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need 2*M*M + 2*M, prefer 2*M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to WORK(IU), zeroing out below it
- *
- CALL DLACPY( 'L', M, M, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO,
- $ WORK( IU+LDWRKU ), LDWRKU )
- *
- * Generate Q in A
- * (Workspace: need 2*M*M + 2*M, prefer 2*M*M + M + M*NB)
- *
- CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in WORK(IU), copying result to
- * WORK(IR)
- * (Workspace: need 2*M*M + 4*M,
- * prefer 2*M*M+3*M+2*M*NB)
- *
- CALL DGEBRD( M, M, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, M, WORK( IU ), LDWRKU,
- $ WORK( IR ), LDWRKR )
- *
- * Generate right bidiagonalizing vectors in WORK(IU)
- * (Workspace: need 2*M*M + 4*M-1,
- * prefer 2*M*M+3*M+(M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IU ), LDWRKU,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in WORK(IR)
- * (Workspace: need 2*M*M + 4*M, prefer 2*M*M + 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, WORK( IR ), LDWRKR,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of L in WORK(IR) and computing
- * right singular vectors of L in WORK(IU)
- * (Workspace: need 2*M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, M, 0, S, WORK( IE ),
- $ WORK( IU ), LDWRKU, WORK( IR ),
- $ LDWRKR, DUM, 1, WORK( IWORK ), INFO )
- *
- * Multiply right singular vectors of L in WORK(IU) by
- * Q in A, storing result in VT
- * (Workspace: need M*M)
- *
- CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IU ),
- $ LDWRKU, A, LDA, ZERO, VT, LDVT )
- *
- * Copy left singular vectors of L to A
- * (Workspace: need M*M)
- *
- CALL DLACPY( 'F', M, M, WORK( IR ), LDWRKR, A,
- $ LDA )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DORGLQ( M, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Zero out above L in A
- *
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ),
- $ LDA )
- *
- * Bidiagonalize L in A
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply right vectors bidiagonalizing L by Q in VT
- * (Workspace: need 3*M + N, prefer 3*M + N*NB)
- *
- CALL DORMBR( 'P', 'L', 'T', M, N, M, A, LDA,
- $ WORK( ITAUP ), VT, LDVT,
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors of L in A
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, A, LDA, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, compute left
- * singular vectors of A in A and compute right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, N, M, 0, S, WORK( IE ), VT,
- $ LDVT, A, LDA, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTUAS ) THEN
- *
- * Path 6t(N much larger than M, JOBU='S' or 'A',
- * JOBVT='S')
- * M right singular vectors to be computed in VT and
- * M left singular vectors to be computed in U
- *
- IF( LWORK.GE.M*M+MAX( 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+LDA*M ) THEN
- *
- * WORK(IU) is LDA by N
- *
- LDWRKU = LDA
- ELSE
- *
- * WORK(IU) is LDA by M
- *
- LDWRKU = M
- END IF
- ITAU = IU + LDWRKU*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to WORK(IU), zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO,
- $ WORK( IU+LDWRKU ), LDWRKU )
- *
- * Generate Q in A
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in WORK(IU), copying result to U
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, M, WORK( IU ), LDWRKU, U,
- $ LDU )
- *
- * Generate right bidiagonalizing vectors in WORK(IU)
- * (Workspace: need M*M + 4*M-1,
- * prefer M*M+3*M+(M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IU ), LDWRKU,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in U
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of L in U and computing right
- * singular vectors of L in WORK(IU)
- * (Workspace: need M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, M, 0, S, WORK( IE ),
- $ WORK( IU ), LDWRKU, U, LDU, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply right singular vectors of L in WORK(IU) by
- * Q in A, storing result in VT
- * (Workspace: need M*M)
- *
- CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IU ),
- $ LDWRKU, A, LDA, ZERO, VT, LDVT )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DORGLQ( M, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to U, zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, U, LDU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, U( 1, 2 ),
- $ LDU )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in U
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, U, LDU, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply right bidiagonalizing vectors in U by Q
- * in VT
- * (Workspace: need 3*M + N, prefer 3*M + N*NB)
- *
- CALL DORMBR( 'P', 'L', 'T', M, N, M, U, LDU,
- $ WORK( ITAUP ), VT, LDVT,
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in U
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U and computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, N, M, 0, S, WORK( IE ), VT,
- $ LDVT, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- END IF
- *
- ELSE IF( WNTVA ) THEN
- *
- IF( WNTUN ) THEN
- *
- * Path 7t(N much larger than M, JOBU='N', JOBVT='A')
- * N right singular vectors to be computed in VT and
- * no left singular vectors to be computed
- *
- IF( LWORK.GE.M*M+MAX( N + M, 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IR = 1
- IF( LWORK.GE.WRKBL+LDA*M ) THEN
- *
- * WORK(IR) is LDA by M
- *
- LDWRKR = LDA
- ELSE
- *
- * WORK(IR) is M by M
- *
- LDWRKR = M
- END IF
- ITAU = IR + LDWRKR*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Copy L to WORK(IR), zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, WORK( IR ),
- $ LDWRKR )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO,
- $ WORK( IR+LDWRKR ), LDWRKR )
- *
- * Generate Q in VT
- * (Workspace: need M*M + M + N, prefer M*M + M + N*NB)
- *
- CALL DORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in WORK(IR)
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, WORK( IR ), LDWRKR, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate right bidiagonalizing vectors in WORK(IR)
- * (Workspace: need M*M + 4*M-1,
- * prefer M*M+3*M+(M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IR ), LDWRKR,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing right
- * singular vectors of L in WORK(IR)
- * (Workspace: need M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, 0, 0, S, WORK( IE ),
- $ WORK( IR ), LDWRKR, DUM, 1, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply right singular vectors of L in WORK(IR) by
- * Q in VT, storing result in A
- * (Workspace: need M*M)
- *
- CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IR ),
- $ LDWRKR, VT, LDVT, ZERO, A, LDA )
- *
- * Copy right singular vectors of A from A to VT
- *
- CALL DLACPY( 'F', M, N, A, LDA, VT, LDVT )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need M + N, prefer M + N*NB)
- *
- CALL DORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Zero out above L in A
- *
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ),
- $ LDA )
- *
- * Bidiagonalize L in A
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply right bidiagonalizing vectors in A by Q
- * in VT
- * (Workspace: need 3*M + N, prefer 3*M + N*NB)
- *
- CALL DORMBR( 'P', 'L', 'T', M, N, M, A, LDA,
- $ WORK( ITAUP ), VT, LDVT,
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, N, 0, 0, S, WORK( IE ), VT,
- $ LDVT, DUM, 1, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTUO ) THEN
- *
- * Path 8t(N much larger than M, JOBU='O', JOBVT='A')
- * N right singular vectors to be computed in VT and
- * M left singular vectors to be overwritten on A
- *
- IF( LWORK.GE.2*M*M+MAX( N + M, 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+2*LDA*M ) THEN
- *
- * WORK(IU) is LDA by M and WORK(IR) is LDA by M
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*M
- LDWRKR = LDA
- ELSE IF( LWORK.GE.WRKBL+( LDA + M )*M ) THEN
- *
- * WORK(IU) is LDA by M and WORK(IR) is M by M
- *
- LDWRKU = LDA
- IR = IU + LDWRKU*M
- LDWRKR = M
- ELSE
- *
- * WORK(IU) is M by M and WORK(IR) is M by M
- *
- LDWRKU = M
- IR = IU + LDWRKU*M
- LDWRKR = M
- END IF
- ITAU = IR + LDWRKR*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need 2*M*M + 2*M, prefer 2*M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need 2*M*M + M + N, prefer 2*M*M + M + N*NB)
- *
- CALL DORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to WORK(IU), zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO,
- $ WORK( IU+LDWRKU ), LDWRKU )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in WORK(IU), copying result to
- * WORK(IR)
- * (Workspace: need 2*M*M + 4*M,
- * prefer 2*M*M+3*M+2*M*NB)
- *
- CALL DGEBRD( M, M, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, M, WORK( IU ), LDWRKU,
- $ WORK( IR ), LDWRKR )
- *
- * Generate right bidiagonalizing vectors in WORK(IU)
- * (Workspace: need 2*M*M + 4*M-1,
- * prefer 2*M*M+3*M+(M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IU ), LDWRKU,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in WORK(IR)
- * (Workspace: need 2*M*M + 4*M, prefer 2*M*M + 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, WORK( IR ), LDWRKR,
- $ WORK( ITAUQ ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of L in WORK(IR) and computing
- * right singular vectors of L in WORK(IU)
- * (Workspace: need 2*M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, M, 0, S, WORK( IE ),
- $ WORK( IU ), LDWRKU, WORK( IR ),
- $ LDWRKR, DUM, 1, WORK( IWORK ), INFO )
- *
- * Multiply right singular vectors of L in WORK(IU) by
- * Q in VT, storing result in A
- * (Workspace: need M*M)
- *
- CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IU ),
- $ LDWRKU, VT, LDVT, ZERO, A, LDA )
- *
- * Copy right singular vectors of A from A to VT
- *
- CALL DLACPY( 'F', M, N, A, LDA, VT, LDVT )
- *
- * Copy left singular vectors of A from WORK(IR) to A
- *
- CALL DLACPY( 'F', M, M, WORK( IR ), LDWRKR, A,
- $ LDA )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need M + N, prefer M + N*NB)
- *
- CALL DORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Zero out above L in A
- *
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ),
- $ LDA )
- *
- * Bidiagonalize L in A
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, A, LDA, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply right bidiagonalizing vectors in A by Q
- * in VT
- * (Workspace: need 3*M + N, prefer 3*M + N*NB)
- *
- CALL DORMBR( 'P', 'L', 'T', M, N, M, A, LDA,
- $ WORK( ITAUP ), VT, LDVT,
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in A
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, A, LDA, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in A and computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, N, M, 0, S, WORK( IE ), VT,
- $ LDVT, A, LDA, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- ELSE IF( WNTUAS ) THEN
- *
- * Path 9t(N much larger than M, JOBU='S' or 'A',
- * JOBVT='A')
- * N right singular vectors to be computed in VT and
- * M left singular vectors to be computed in U
- *
- IF( LWORK.GE.M*M+MAX( N + M, 4*M, BDSPAC ) ) THEN
- *
- * Sufficient workspace for a fast algorithm
- *
- IU = 1
- IF( LWORK.GE.WRKBL+LDA*M ) THEN
- *
- * WORK(IU) is LDA by M
- *
- LDWRKU = LDA
- ELSE
- *
- * WORK(IU) is M by M
- *
- LDWRKU = M
- END IF
- ITAU = IU + LDWRKU*M
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need M*M + 2*M, prefer M*M + M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need M*M + M + N, prefer M*M + M + N*NB)
- *
- CALL DORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to WORK(IU), zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, WORK( IU ),
- $ LDWRKU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO,
- $ WORK( IU+LDWRKU ), LDWRKU )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in WORK(IU), copying result to U
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, WORK( IU ), LDWRKU, S,
- $ WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- CALL DLACPY( 'L', M, M, WORK( IU ), LDWRKU, U,
- $ LDU )
- *
- * Generate right bidiagonalizing vectors in WORK(IU)
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + (M-1)*NB)
- *
- CALL DORGBR( 'P', M, M, M, WORK( IU ), LDWRKU,
- $ WORK( ITAUP ), WORK( IWORK ),
- $ LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in U
- * (Workspace: need M*M + 4*M, prefer M*M + 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of L in U and computing right
- * singular vectors of L in WORK(IU)
- * (Workspace: need M*M + BDSPAC)
- *
- CALL DBDSQR( 'U', M, M, M, 0, S, WORK( IE ),
- $ WORK( IU ), LDWRKU, U, LDU, DUM, 1,
- $ WORK( IWORK ), INFO )
- *
- * Multiply right singular vectors of L in WORK(IU) by
- * Q in VT, storing result in A
- * (Workspace: need M*M)
- *
- CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IU ),
- $ LDWRKU, VT, LDVT, ZERO, A, LDA )
- *
- * Copy right singular vectors of A from A to VT
- *
- CALL DLACPY( 'F', M, N, A, LDA, VT, LDVT )
- *
- ELSE
- *
- * Insufficient workspace for a fast algorithm
- *
- ITAU = 1
- IWORK = ITAU + M
- *
- * Compute A=L*Q, copying result to VT
- * (Workspace: need 2*M, prefer M + M*NB)
- *
- CALL DGELQF( M, N, A, LDA, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- *
- * Generate Q in VT
- * (Workspace: need M + N, prefer M + N*NB)
- *
- CALL DORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Copy L to U, zeroing out above it
- *
- CALL DLACPY( 'L', M, M, A, LDA, U, LDU )
- CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, U( 1, 2 ),
- $ LDU )
- IE = ITAU
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize L in U
- * (Workspace: need 4*M, prefer 3*M + 2*M*NB)
- *
- CALL DGEBRD( M, M, U, LDU, S, WORK( IE ),
- $ WORK( ITAUQ ), WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Multiply right bidiagonalizing vectors in U by Q
- * in VT
- * (Workspace: need 3*M + N, prefer 3*M + N*NB)
- *
- CALL DORMBR( 'P', 'L', 'T', M, N, M, U, LDU,
- $ WORK( ITAUP ), VT, LDVT,
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- *
- * Generate left bidiagonalizing vectors in U
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'Q', M, M, M, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- IWORK = IE + M
- *
- * Perform bidiagonal QR iteration, computing left
- * singular vectors of A in U and computing right
- * singular vectors of A in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'U', M, N, M, 0, S, WORK( IE ), VT,
- $ LDVT, U, LDU, DUM, 1, WORK( IWORK ),
- $ INFO )
- *
- END IF
- *
- END IF
- *
- END IF
- *
- ELSE
- *
- * N .LT. MNTHR
- *
- * Path 10t(N greater than M, but not much larger)
- * Reduce to bidiagonal form without LQ decomposition
- *
- IE = 1
- ITAUQ = IE + M
- ITAUP = ITAUQ + M
- IWORK = ITAUP + M
- *
- * Bidiagonalize A
- * (Workspace: need 3*M + N, prefer 3*M + (M + N)*NB)
- *
- CALL DGEBRD( M, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
- $ WORK( ITAUP ), WORK( IWORK ), LWORK-IWORK+1,
- $ IERR )
- IF( WNTUAS ) THEN
- *
- * If left singular vectors desired in U, copy result to U
- * and generate left bidiagonalizing vectors in U
- * (Workspace: need 4*M-1, prefer 3*M + (M-1)*NB)
- *
- CALL DLACPY( 'L', M, M, A, LDA, U, LDU )
- CALL DORGBR( 'Q', M, M, N, U, LDU, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IF( WNTVAS ) THEN
- *
- * If right singular vectors desired in VT, copy result to
- * VT and generate right bidiagonalizing vectors in VT
- * (Workspace: need 3*M + NRVT, prefer 3*M + NRVT*NB)
- *
- CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
- IF( WNTVA )
- $ NRVT = N
- IF( WNTVS )
- $ NRVT = M
- CALL DORGBR( 'P', NRVT, N, M, VT, LDVT, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IF( WNTUO ) THEN
- *
- * If left singular vectors desired in A, generate left
- * bidiagonalizing vectors in A
- * (Workspace: need 4*M-1, prefer 3*M + (M-1)*NB)
- *
- CALL DORGBR( 'Q', M, M, N, A, LDA, WORK( ITAUQ ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IF( WNTVO ) THEN
- *
- * If right singular vectors desired in A, generate right
- * bidiagonalizing vectors in A
- * (Workspace: need 4*M, prefer 3*M + M*NB)
- *
- CALL DORGBR( 'P', M, N, M, A, LDA, WORK( ITAUP ),
- $ WORK( IWORK ), LWORK-IWORK+1, IERR )
- END IF
- IWORK = IE + M
- IF( WNTUAS .OR. WNTUO )
- $ NRU = M
- IF( WNTUN )
- $ NRU = 0
- IF( WNTVAS .OR. WNTVO )
- $ NCVT = N
- IF( WNTVN )
- $ NCVT = 0
- IF( ( .NOT.WNTUO ) .AND. ( .NOT.WNTVO ) ) THEN
- *
- * Perform bidiagonal QR iteration, if desired, computing
- * left singular vectors in U and computing right singular
- * vectors in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'L', M, NCVT, NRU, 0, S, WORK( IE ), VT,
- $ LDVT, U, LDU, DUM, 1, WORK( IWORK ), INFO )
- ELSE IF( ( .NOT.WNTUO ) .AND. WNTVO ) THEN
- *
- * Perform bidiagonal QR iteration, if desired, computing
- * left singular vectors in U and computing right singular
- * vectors in A
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'L', M, NCVT, NRU, 0, S, WORK( IE ), A, LDA,
- $ U, LDU, DUM, 1, WORK( IWORK ), INFO )
- ELSE
- *
- * Perform bidiagonal QR iteration, if desired, computing
- * left singular vectors in A and computing right singular
- * vectors in VT
- * (Workspace: need BDSPAC)
- *
- CALL DBDSQR( 'L', M, NCVT, NRU, 0, S, WORK( IE ), VT,
- $ LDVT, A, LDA, DUM, 1, WORK( IWORK ), INFO )
- END IF
- *
- END IF
- *
- END IF
- *
- * If DBDSQR failed to converge, copy unconverged superdiagonals
- * to WORK( 2:MINMN )
- *
- IF( INFO.NE.0 ) THEN
- IF( IE.GT.2 ) THEN
- DO 50 I = 1, MINMN - 1
- WORK( I+1 ) = WORK( I+IE-1 )
- 50 CONTINUE
- END IF
- IF( IE.LT.2 ) THEN
- DO 60 I = MINMN - 1, 1, -1
- WORK( I+1 ) = WORK( I+IE-1 )
- 60 CONTINUE
- END IF
- END IF
- *
- * Undo scaling if necessary
- *
- IF( ISCL.EQ.1 ) THEN
- IF( ANRM.GT.BIGNUM )
- $ CALL DLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN, 1, S, MINMN,
- $ IERR )
- IF( INFO.NE.0 .AND. ANRM.GT.BIGNUM )
- $ CALL DLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN-1, 1, WORK( 2 ),
- $ MINMN, IERR )
- IF( ANRM.LT.SMLNUM )
- $ CALL DLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN, 1, S, MINMN,
- $ IERR )
- IF( INFO.NE.0 .AND. ANRM.LT.SMLNUM )
- $ CALL DLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN-1, 1, WORK( 2 ),
- $ MINMN, IERR )
- END IF
- *
- * Return optimal workspace in WORK(1)
- *
- WORK( 1 ) = MAXWRK
- *
- RETURN
- *
- * End of DGESVD
- *
- END
|