1 /******************************************************************************
2  * Copyright 1998-2019 Lawrence Livermore National Security, LLC and other
3  * HYPRE Project Developers. See the top-level COPYRIGHT file for details.
4  *
5  * SPDX-License-Identifier: (Apache-2.0 OR MIT)
6  ******************************************************************************/
7 
8 /******************************************************************************
9  *
10  * HYPRE_LOBPCG interface
11  *
12  *****************************************************************************/
13 
14 #include "_hypre_utilities.h"
15 
16 #include "HYPRE_config.h"
17 
18 #include "HYPRE_lobpcg.h"
19 #include "lobpcg.h"
20 
21 #include "interpreter.h"
22 #include "HYPRE_MatvecFunctions.h"
23 
24 #include "_hypre_lapack.h"
25 
26 typedef struct
27 {
28 	HYPRE_Int    (*Precond)(void*,void*,void*,void*);
29    HYPRE_Int    (*PrecondSetup)(void*,void*,void*,void*);
30 
31 } hypre_LOBPCGPrecond;
32 
33 typedef struct
34 {
35    lobpcg_Tolerance              tolerance;
36    HYPRE_Int                           maxIterations;
37    HYPRE_Int                           verbosityLevel;
38    HYPRE_Int                           precondUsageMode;
39    HYPRE_Int                           iterationNumber;
40    utilities_FortranMatrix*      eigenvaluesHistory;
41    utilities_FortranMatrix*      residualNorms;
42    utilities_FortranMatrix*      residualNormsHistory;
43 
44 } lobpcg_Data;
45 
46 #define lobpcg_tolerance(data)            ((data).tolerance)
47 #define lobpcg_absoluteTolerance(data)    ((data).tolerance.absolute)
48 #define lobpcg_relativeTolerance(data)    ((data).tolerance.relative)
49 #define lobpcg_maxIterations(data)        ((data).maxIterations)
50 #define lobpcg_verbosityLevel(data)       ((data).verbosityLevel)
51 #define lobpcg_precondUsageMode(data)     ((data).precondUsageMode)
52 #define lobpcg_iterationNumber(data)      ((data).iterationNumber)
53 #define lobpcg_eigenvaluesHistory(data)   ((data).eigenvaluesHistory)
54 #define lobpcg_residualNorms(data)        ((data).residualNorms)
55 #define lobpcg_residualNormsHistory(data) ((data).residualNormsHistory)
56 
57 typedef struct
58 {
59 
60    lobpcg_Data                   lobpcgData;
61 
62    mv_InterfaceInterpreter*      interpreter;
63 
64    void*                         A;
65    void*                         matvecData;
66    void*                         precondData;
67 
68    void*                         B;
69    void*                         matvecDataB;
70    void*                         T;
71    void*                         matvecDataT;
72 
73    hypre_LOBPCGPrecond           precondFunctions;
74 
75    HYPRE_MatvecFunctions*        matvecFunctions;
76 
77 } hypre_LOBPCGData;
78 
dsygv_interface(HYPRE_Int * itype,char * jobz,char * uplo,HYPRE_Int * n,HYPRE_Real * a,HYPRE_Int * lda,HYPRE_Real * b,HYPRE_Int * ldb,HYPRE_Real * w,HYPRE_Real * work,HYPRE_Int * lwork,HYPRE_Int * info)79 static HYPRE_Int dsygv_interface (HYPRE_Int *itype, char *jobz, char *uplo, HYPRE_Int *
80                             n, HYPRE_Real *a, HYPRE_Int *lda, HYPRE_Real *b, HYPRE_Int *ldb,
81                             HYPRE_Real *w, HYPRE_Real *work, HYPRE_Int *lwork, HYPRE_Int *info)
82 {
83    hypre_dsygv(itype, jobz, uplo, n, a, lda, b, ldb, w, work, lwork, info);
84    return 0;
85 }
86 
dpotrf_interface(const char * uplo,HYPRE_Int * n,HYPRE_Real * a,HYPRE_Int * lda,HYPRE_Int * info)87 static HYPRE_Int dpotrf_interface (const char *uplo, HYPRE_Int *n, HYPRE_Real *a, HYPRE_Int *
88                              lda, HYPRE_Int *info)
89 {
90    hypre_dpotrf(uplo, n, a, lda, info);
91    return 0;
92 }
93 
94 
95 HYPRE_Int
lobpcg_initialize(lobpcg_Data * data)96 lobpcg_initialize( lobpcg_Data* data )
97 {
98    (data->tolerance).absolute    = 1.0e-06;
99    (data->tolerance).relative    = 1.0e-06;
100    (data->maxIterations)         = 500;
101    (data->precondUsageMode)      = 0;
102    (data->verbosityLevel)        = 0;
103    (data->eigenvaluesHistory)    = utilities_FortranMatrixCreate();
104    (data->residualNorms)         = utilities_FortranMatrixCreate();
105    (data->residualNormsHistory)  = utilities_FortranMatrixCreate();
106 
107    return 0;
108 }
109 
110 HYPRE_Int
lobpcg_clean(lobpcg_Data * data)111 lobpcg_clean( lobpcg_Data* data )
112 {
113    utilities_FortranMatrixDestroy( data->eigenvaluesHistory );
114    utilities_FortranMatrixDestroy( data->residualNorms );
115    utilities_FortranMatrixDestroy( data->residualNormsHistory );
116 
117    return 0;
118 }
119 
120 HYPRE_Int
hypre_LOBPCGDestroy(void * pcg_vdata)121 hypre_LOBPCGDestroy( void *pcg_vdata )
122 {
123 	hypre_LOBPCGData      *pcg_data      = (hypre_LOBPCGData*)pcg_vdata;
124 
125    if (pcg_data) {
126       HYPRE_MatvecFunctions * mv = pcg_data->matvecFunctions;
127       if ( pcg_data->matvecData != NULL ) {
128          (*(mv->MatvecDestroy))(pcg_data->matvecData);
129          pcg_data->matvecData = NULL;
130       }
131       if ( pcg_data->matvecDataB != NULL ) {
132          (*(mv->MatvecDestroy))(pcg_data->matvecDataB);
133          pcg_data->matvecDataB = NULL;
134       }
135       if ( pcg_data->matvecDataT != NULL ) {
136          (*(mv->MatvecDestroy))(pcg_data->matvecDataT);
137          pcg_data->matvecDataT = NULL;
138       }
139 
140       lobpcg_clean( &(pcg_data->lobpcgData) );
141 
142       hypre_TFree( pcg_vdata , HYPRE_MEMORY_HOST);
143    }
144 
145    return hypre_error_flag;
146 }
147 
148 HYPRE_Int
hypre_LOBPCGSetup(void * pcg_vdata,void * A,void * b,void * x)149 hypre_LOBPCGSetup( void *pcg_vdata, void *A, void *b, void *x )
150 {
151    hypre_LOBPCGData *pcg_data = (hypre_LOBPCGData*)pcg_vdata;
152    HYPRE_MatvecFunctions * mv = pcg_data->matvecFunctions;
153    HYPRE_Int  (*precond_setup)(void*,void*,void*,void*) = (pcg_data->precondFunctions).PrecondSetup;
154    void *precond_data = (pcg_data->precondData);
155 
156    (pcg_data->A) = A;
157 
158    if ( pcg_data->matvecData != NULL )
159       (*(mv->MatvecDestroy))(pcg_data->matvecData);
160    (pcg_data->matvecData) = (*(mv->MatvecCreate))(A, x);
161 
162    if ( precond_setup != NULL ) {
163       if ( pcg_data->T == NULL )
164          precond_setup(precond_data, A, b, x);
165       else
166          precond_setup(precond_data, pcg_data->T, b, x);
167    }
168 
169    return hypre_error_flag;
170 }
171 
172 HYPRE_Int
hypre_LOBPCGSetupB(void * pcg_vdata,void * B,void * x)173 hypre_LOBPCGSetupB( void *pcg_vdata, void *B, void *x )
174 {
175    hypre_LOBPCGData *pcg_data = (hypre_LOBPCGData*)pcg_vdata;
176    HYPRE_MatvecFunctions * mv = pcg_data->matvecFunctions;
177 
178    (pcg_data->B) = B;
179 
180    if ( pcg_data->matvecDataB != NULL )
181       (*(mv->MatvecDestroy))(pcg_data -> matvecDataB);
182    (pcg_data->matvecDataB) = (*(mv->MatvecCreate))(B, x);
183    if ( B != NULL )
184       (pcg_data->matvecDataB) = (*(mv->MatvecCreate))(B, x);
185    else
186       (pcg_data->matvecDataB) = NULL;
187 
188    return hypre_error_flag;
189 }
190 
191 HYPRE_Int
hypre_LOBPCGSetupT(void * pcg_vdata,void * T,void * x)192 hypre_LOBPCGSetupT( void *pcg_vdata, void *T, void *x )
193 {
194    hypre_LOBPCGData *pcg_data = (hypre_LOBPCGData*)pcg_vdata;
195    HYPRE_MatvecFunctions * mv = pcg_data->matvecFunctions;
196 
197    (pcg_data -> T) = T;
198 
199    if ( pcg_data->matvecDataT != NULL )
200       (*(mv->MatvecDestroy))(pcg_data->matvecDataT);
201    if ( T != NULL )
202       (pcg_data->matvecDataT) = (*(mv->MatvecCreate))(T, x);
203    else
204       (pcg_data->matvecDataT) = NULL;
205 
206    return hypre_error_flag;
207 }
208 
209 HYPRE_Int
hypre_LOBPCGSetTol(void * pcg_vdata,HYPRE_Real tol)210 hypre_LOBPCGSetTol( void* pcg_vdata, HYPRE_Real tol )
211 {
212    hypre_LOBPCGData *pcg_data	= (hypre_LOBPCGData*)pcg_vdata;
213 
214    lobpcg_absoluteTolerance(pcg_data->lobpcgData) = tol;
215 
216    return hypre_error_flag;
217 }
218 
219 HYPRE_Int
hypre_LOBPCGSetRTol(void * pcg_vdata,HYPRE_Real tol)220 hypre_LOBPCGSetRTol( void* pcg_vdata, HYPRE_Real tol )
221 {
222    hypre_LOBPCGData *pcg_data	= (hypre_LOBPCGData*) pcg_vdata;
223 
224    lobpcg_relativeTolerance(pcg_data->lobpcgData) = tol;
225 
226    return hypre_error_flag;
227 }
228 
229 HYPRE_Int
hypre_LOBPCGSetMaxIter(void * pcg_vdata,HYPRE_Int max_iter)230 hypre_LOBPCGSetMaxIter( void* pcg_vdata, HYPRE_Int max_iter  )
231 {
232    hypre_LOBPCGData *pcg_data	= (hypre_LOBPCGData*)pcg_vdata;
233 
234    lobpcg_maxIterations(pcg_data->lobpcgData) = max_iter;
235 
236    return hypre_error_flag;
237 }
238 
239 HYPRE_Int
hypre_LOBPCGSetPrecondUsageMode(void * pcg_vdata,HYPRE_Int mode)240 hypre_LOBPCGSetPrecondUsageMode( void* pcg_vdata, HYPRE_Int mode  )
241 {
242    hypre_LOBPCGData *pcg_data	= (hypre_LOBPCGData*)pcg_vdata;
243 
244    lobpcg_precondUsageMode(pcg_data->lobpcgData) = mode;
245 
246    return hypre_error_flag;
247 }
248 
249 HYPRE_Int
hypre_LOBPCGGetPrecond(void * pcg_vdata,HYPRE_Solver * precond_data_ptr)250 hypre_LOBPCGGetPrecond( void         *pcg_vdata,
251 			HYPRE_Solver *precond_data_ptr )
252 {
253    hypre_LOBPCGData*	pcg_data	= (hypre_LOBPCGData*)pcg_vdata;
254 
255    *precond_data_ptr = (HYPRE_Solver)(pcg_data -> precondData);
256 
257    return hypre_error_flag;
258 }
259 
260 HYPRE_Int
hypre_LOBPCGSetPrecond(void * pcg_vdata,HYPRE_Int (* precond)(void *,void *,void *,void *),HYPRE_Int (* precond_setup)(void *,void *,void *,void *),void * precond_data)261 hypre_LOBPCGSetPrecond( void  *pcg_vdata,
262 						HYPRE_Int  (*precond)(void*,void*,void*,void*),
263 						HYPRE_Int  (*precond_setup)(void*,void*,void*,void*),
264 			void  *precond_data )
265 {
266    hypre_LOBPCGData* pcg_data = (hypre_LOBPCGData*)pcg_vdata;
267 
268    (pcg_data->precondFunctions).Precond      = precond;
269    (pcg_data->precondFunctions).PrecondSetup = precond_setup;
270    (pcg_data->precondData)                   = precond_data;
271 
272    return hypre_error_flag;
273 }
274 
275 HYPRE_Int
hypre_LOBPCGSetPrintLevel(void * pcg_vdata,HYPRE_Int level)276 hypre_LOBPCGSetPrintLevel( void *pcg_vdata, HYPRE_Int level )
277 {
278    hypre_LOBPCGData *pcg_data = (hypre_LOBPCGData*)pcg_vdata;
279 
280    lobpcg_verbosityLevel(pcg_data->lobpcgData) = level;
281 
282    return hypre_error_flag;
283 }
284 
285 void
hypre_LOBPCGPreconditioner(void * vdata,void * x,void * y)286 hypre_LOBPCGPreconditioner( void *vdata, void* x, void* y )
287 {
288    hypre_LOBPCGData *data = (hypre_LOBPCGData*)vdata;
289    mv_InterfaceInterpreter* ii = data->interpreter;
290    HYPRE_Int (*precond)(void*,void*,void*,void*) = (data->precondFunctions).Precond;
291 
292    if ( precond == NULL ) {
293       (*(ii->CopyVector))(x,y);
294       return;
295    }
296 
297    if ( lobpcg_precondUsageMode(data->lobpcgData) == 0 )
298       (*(ii->ClearVector))(y);
299    else
300       (*(ii->CopyVector))(x,y);
301 
302    if ( data->T == NULL )
303       precond(data->precondData, data->A, x, y);
304    else
305       precond(data->precondData, data->T, x, y);
306 }
307 
308 void
hypre_LOBPCGOperatorA(void * pcg_vdata,void * x,void * y)309 hypre_LOBPCGOperatorA( void *pcg_vdata, void* x, void* y )
310 {
311    hypre_LOBPCGData*           pcg_data    = (hypre_LOBPCGData*)pcg_vdata;
312    HYPRE_MatvecFunctions * mv = pcg_data->matvecFunctions;
313    void*	              	      matvec_data = (pcg_data -> matvecData);
314 
315    (*(mv->Matvec))(matvec_data, 1.0, pcg_data->A, x, 0.0, y);
316 }
317 
318 void
hypre_LOBPCGOperatorB(void * pcg_vdata,void * x,void * y)319 hypre_LOBPCGOperatorB( void *pcg_vdata, void* x, void* y )
320 {
321    hypre_LOBPCGData*           pcg_data    = (hypre_LOBPCGData*)pcg_vdata;
322    mv_InterfaceInterpreter* ii          = pcg_data->interpreter;
323    HYPRE_MatvecFunctions * mv = pcg_data->matvecFunctions;
324    void*                       matvec_data = (pcg_data -> matvecDataB);
325 
326    if ( pcg_data->B == NULL ) {
327       (*(ii->CopyVector))(x, y);
328 
329       /* a test */
330       /*
331         (*(ii->ScaleVector))(2.0, y);
332       */
333 
334       return;
335    }
336 
337    (*(mv->Matvec))(matvec_data, 1.0, pcg_data->B, x, 0.0, y);
338 }
339 
340 void
hypre_LOBPCGMultiPreconditioner(void * data,void * x,void * y)341 hypre_LOBPCGMultiPreconditioner( void *data, void * x, void*  y )
342 {
343    hypre_LOBPCGData *pcg_data = (hypre_LOBPCGData*)data;
344    mv_InterfaceInterpreter* ii = pcg_data->interpreter;
345 
346    ii->Eval( hypre_LOBPCGPreconditioner, data, x, y );
347 }
348 
349 void
hypre_LOBPCGMultiOperatorA(void * data,void * x,void * y)350 hypre_LOBPCGMultiOperatorA( void *data, void * x, void*  y )
351 {
352    hypre_LOBPCGData *pcg_data = (hypre_LOBPCGData*)data;
353    mv_InterfaceInterpreter* ii = pcg_data->interpreter;
354 
355    ii->Eval( hypre_LOBPCGOperatorA, data, x, y );
356 }
357 
358 void
hypre_LOBPCGMultiOperatorB(void * data,void * x,void * y)359 hypre_LOBPCGMultiOperatorB( void *data, void * x, void*  y )
360 {
361    hypre_LOBPCGData *pcg_data = (hypre_LOBPCGData*)data;
362    mv_InterfaceInterpreter* ii = pcg_data->interpreter;
363 
364    ii->Eval( hypre_LOBPCGOperatorB, data, x, y );
365 }
366 
367 HYPRE_Int
hypre_LOBPCGSolve(void * vdata,mv_MultiVectorPtr con,mv_MultiVectorPtr vec,HYPRE_Real * val)368 hypre_LOBPCGSolve( void *vdata,
369 		   mv_MultiVectorPtr con,
370 		   mv_MultiVectorPtr vec,
371 		   HYPRE_Real* val )
372 {
373    hypre_LOBPCGData* data = (hypre_LOBPCGData*)vdata;
374    HYPRE_Int (*precond)(void*,void*,void*,void*) = (data->precondFunctions).Precond;
375    void* opB = data->B;
376 
377    void (*prec)( void*, void*, void* );
378    void (*operatorA)( void*, void*, void* );
379    void (*operatorB)( void*, void*, void* );
380 
381    HYPRE_Int maxit = lobpcg_maxIterations(data->lobpcgData);
382    HYPRE_Int verb  = lobpcg_verbosityLevel(data->lobpcgData);
383 
384    HYPRE_Int n	= mv_MultiVectorWidth( vec );
385    lobpcg_BLASLAPACKFunctions blap_fn;
386 
387    utilities_FortranMatrix* lambdaHistory;
388    utilities_FortranMatrix* residuals;
389    utilities_FortranMatrix* residualsHistory;
390 
391    lambdaHistory	= lobpcg_eigenvaluesHistory(data->lobpcgData);
392    residuals = lobpcg_residualNorms(data->lobpcgData);
393    residualsHistory = lobpcg_residualNormsHistory(data->lobpcgData);
394 
395    utilities_FortranMatrixAllocateData( n, maxit + 1,	lambdaHistory );
396    utilities_FortranMatrixAllocateData( n, 1,		residuals );
397    utilities_FortranMatrixAllocateData( n, maxit + 1,	residualsHistory );
398 
399    if ( precond != NULL )
400       prec = hypre_LOBPCGMultiPreconditioner;
401    else
402       prec = NULL;
403 
404    operatorA = hypre_LOBPCGMultiOperatorA;
405 
406    if ( opB != NULL )
407       operatorB = hypre_LOBPCGMultiOperatorB;
408    else
409       operatorB = NULL;
410 
411    blap_fn.dsygv = dsygv_interface;
412    blap_fn.dpotrf = dpotrf_interface;
413 
414    lobpcg_solve( vec,
415                  vdata, operatorA,
416                  vdata, operatorB,
417                  vdata, prec,
418                  con,
419                  blap_fn,
420                  lobpcg_tolerance(data->lobpcgData), maxit, verb,
421                  &(lobpcg_iterationNumber(data->lobpcgData)),
422                  val,
423                  utilities_FortranMatrixValues(lambdaHistory),
424                  utilities_FortranMatrixGlobalHeight(lambdaHistory),
425                  utilities_FortranMatrixValues(residuals),
426                  utilities_FortranMatrixValues(residualsHistory),
427                  utilities_FortranMatrixGlobalHeight(residualsHistory)
428       );
429 
430    return hypre_error_flag;
431 }
432 
433 utilities_FortranMatrix*
hypre_LOBPCGResidualNorms(void * vdata)434 hypre_LOBPCGResidualNorms( void *vdata )
435 {
436    hypre_LOBPCGData *data = (hypre_LOBPCGData*)vdata;
437    return (lobpcg_residualNorms(data->lobpcgData));
438 }
439 
440 utilities_FortranMatrix*
hypre_LOBPCGResidualNormsHistory(void * vdata)441 hypre_LOBPCGResidualNormsHistory( void *vdata )
442 {
443    hypre_LOBPCGData *data = (hypre_LOBPCGData*)vdata;
444    return (lobpcg_residualNormsHistory(data->lobpcgData));
445 }
446 
447 utilities_FortranMatrix*
hypre_LOBPCGEigenvaluesHistory(void * vdata)448 hypre_LOBPCGEigenvaluesHistory( void *vdata )
449 {
450    hypre_LOBPCGData *data = (hypre_LOBPCGData*)vdata;
451    return (lobpcg_eigenvaluesHistory(data->lobpcgData));
452 }
453 
454 HYPRE_Int
hypre_LOBPCGIterations(void * vdata)455 hypre_LOBPCGIterations( void* vdata )
456 {
457    hypre_LOBPCGData *data = (hypre_LOBPCGData*)vdata;
458    return (lobpcg_iterationNumber(data->lobpcgData));
459 }
460 
461 
462 HYPRE_Int
HYPRE_LOBPCGCreate(mv_InterfaceInterpreter * ii,HYPRE_MatvecFunctions * mv,HYPRE_Solver * solver)463 HYPRE_LOBPCGCreate( mv_InterfaceInterpreter* ii, HYPRE_MatvecFunctions* mv,
464                     HYPRE_Solver* solver )
465 {
466    hypre_LOBPCGData *pcg_data;
467 
468    pcg_data = hypre_CTAlloc(hypre_LOBPCGData, 1, HYPRE_MEMORY_HOST);
469 
470    (pcg_data->precondFunctions).Precond = NULL;
471    (pcg_data->precondFunctions).PrecondSetup = NULL;
472 
473    /* set defaults */
474 
475    (pcg_data->interpreter)               = ii;
476    pcg_data->matvecFunctions             = mv;
477 
478    (pcg_data->matvecData)	       	= NULL;
479    (pcg_data->B)	       			= NULL;
480    (pcg_data->matvecDataB)	       	= NULL;
481    (pcg_data->T)	       			= NULL;
482    (pcg_data->matvecDataT)	       	= NULL;
483    (pcg_data->precondData)	       	= NULL;
484 
485    lobpcg_initialize( &(pcg_data->lobpcgData) );
486 
487    *solver = (HYPRE_Solver)pcg_data;
488 
489    return hypre_error_flag;
490 }
491 
492 HYPRE_Int
HYPRE_LOBPCGDestroy(HYPRE_Solver solver)493 HYPRE_LOBPCGDestroy( HYPRE_Solver solver )
494 {
495    return( hypre_LOBPCGDestroy( (void *) solver ) );
496 }
497 
498 HYPRE_Int
HYPRE_LOBPCGSetup(HYPRE_Solver solver,HYPRE_Matrix A,HYPRE_Vector b,HYPRE_Vector x)499 HYPRE_LOBPCGSetup( HYPRE_Solver solver,
500                    HYPRE_Matrix A,
501                    HYPRE_Vector b,
502                    HYPRE_Vector x      )
503 {
504    return( hypre_LOBPCGSetup( solver, A, b, x ) );
505 }
506 
507 HYPRE_Int
HYPRE_LOBPCGSetupB(HYPRE_Solver solver,HYPRE_Matrix B,HYPRE_Vector x)508 HYPRE_LOBPCGSetupB( HYPRE_Solver solver,
509                     HYPRE_Matrix B,
510                     HYPRE_Vector x      )
511 {
512    return( hypre_LOBPCGSetupB( solver, B, x ) );
513 }
514 
515 HYPRE_Int
HYPRE_LOBPCGSetupT(HYPRE_Solver solver,HYPRE_Matrix T,HYPRE_Vector x)516 HYPRE_LOBPCGSetupT( HYPRE_Solver solver,
517                     HYPRE_Matrix T,
518                     HYPRE_Vector x      )
519 {
520    return( hypre_LOBPCGSetupT( solver, T, x ) );
521 }
522 
523 HYPRE_Int
HYPRE_LOBPCGSolve(HYPRE_Solver solver,mv_MultiVectorPtr con,mv_MultiVectorPtr vec,HYPRE_Real * val)524 HYPRE_LOBPCGSolve( HYPRE_Solver solver, mv_MultiVectorPtr con,
525 		   mv_MultiVectorPtr vec, HYPRE_Real* val )
526 {
527    return( hypre_LOBPCGSolve( (void *) solver, con, vec, val ) );
528 }
529 
530 HYPRE_Int
HYPRE_LOBPCGSetTol(HYPRE_Solver solver,HYPRE_Real tol)531 HYPRE_LOBPCGSetTol( HYPRE_Solver solver, HYPRE_Real tol )
532 {
533    return( hypre_LOBPCGSetTol( (void *) solver, tol ) );
534 }
535 
536 HYPRE_Int
HYPRE_LOBPCGSetRTol(HYPRE_Solver solver,HYPRE_Real tol)537 HYPRE_LOBPCGSetRTol( HYPRE_Solver solver, HYPRE_Real tol )
538 {
539    return( hypre_LOBPCGSetRTol( (void *) solver, tol ) );
540 }
541 
542 HYPRE_Int
HYPRE_LOBPCGSetMaxIter(HYPRE_Solver solver,HYPRE_Int max_iter)543 HYPRE_LOBPCGSetMaxIter( HYPRE_Solver solver, HYPRE_Int max_iter )
544 {
545    return( hypre_LOBPCGSetMaxIter( (void *) solver, max_iter ) );
546 }
547 
548 HYPRE_Int
HYPRE_LOBPCGSetPrecondUsageMode(HYPRE_Solver solver,HYPRE_Int mode)549 HYPRE_LOBPCGSetPrecondUsageMode( HYPRE_Solver solver, HYPRE_Int mode )
550 {
551    return( hypre_LOBPCGSetPrecondUsageMode( (void *) solver, mode ) );
552 }
553 
554 HYPRE_Int
HYPRE_LOBPCGSetPrecond(HYPRE_Solver solver,HYPRE_PtrToSolverFcn precond,HYPRE_PtrToSolverFcn precond_setup,HYPRE_Solver precond_solver)555 HYPRE_LOBPCGSetPrecond( HYPRE_Solver         solver,
556                         HYPRE_PtrToSolverFcn precond,
557                         HYPRE_PtrToSolverFcn precond_setup,
558                         HYPRE_Solver         precond_solver )
559 {
560    return( hypre_LOBPCGSetPrecond( (void *) solver,
561                                    (HYPRE_Int (*)(void*, void*, void*, void*))precond,
562 								   (HYPRE_Int (*)(void*, void*, void*, void*))precond_setup,
563                                    (void *) precond_solver ) );
564 }
565 
566 HYPRE_Int
HYPRE_LOBPCGGetPrecond(HYPRE_Solver solver,HYPRE_Solver * precond_data_ptr)567 HYPRE_LOBPCGGetPrecond( HYPRE_Solver  solver,
568                         HYPRE_Solver *precond_data_ptr )
569 {
570    return( hypre_LOBPCGGetPrecond( (void *)     solver,
571                                    (HYPRE_Solver *) precond_data_ptr ) );
572 }
573 
574 HYPRE_Int
HYPRE_LOBPCGSetPrintLevel(HYPRE_Solver solver,HYPRE_Int level)575 HYPRE_LOBPCGSetPrintLevel( HYPRE_Solver solver, HYPRE_Int level )
576 {
577    return( hypre_LOBPCGSetPrintLevel( (void*)solver, level ) );
578 }
579 
580 utilities_FortranMatrix*
HYPRE_LOBPCGResidualNorms(HYPRE_Solver solver)581 HYPRE_LOBPCGResidualNorms( HYPRE_Solver solver )
582 {
583    return ( hypre_LOBPCGResidualNorms( (void*)solver ) );
584 }
585 
586 utilities_FortranMatrix*
HYPRE_LOBPCGResidualNormsHistory(HYPRE_Solver solver)587 HYPRE_LOBPCGResidualNormsHistory( HYPRE_Solver solver )
588 {
589    return ( hypre_LOBPCGResidualNormsHistory( (void*)solver ) );
590 }
591 
592 utilities_FortranMatrix*
HYPRE_LOBPCGEigenvaluesHistory(HYPRE_Solver solver)593 HYPRE_LOBPCGEigenvaluesHistory( HYPRE_Solver solver )
594 {
595    return ( hypre_LOBPCGEigenvaluesHistory( (void*)solver ) );
596 }
597 
598 HYPRE_Int
HYPRE_LOBPCGIterations(HYPRE_Solver solver)599 HYPRE_LOBPCGIterations( HYPRE_Solver solver )
600 {
601    return ( hypre_LOBPCGIterations( (void*)solver ) );
602 }
603 
604 void
lobpcg_MultiVectorByMultiVector(mv_MultiVectorPtr x,mv_MultiVectorPtr y,utilities_FortranMatrix * xy)605 lobpcg_MultiVectorByMultiVector( mv_MultiVectorPtr x,
606                                  mv_MultiVectorPtr y,
607                                  utilities_FortranMatrix* xy )
608 {
609    mv_MultiVectorByMultiVector( x, y,
610                                 utilities_FortranMatrixGlobalHeight( xy ),
611                                 utilities_FortranMatrixHeight( xy ),
612                                 utilities_FortranMatrixWidth( xy ),
613                                 utilities_FortranMatrixValues( xy ) );
614 }
615 
616