void lis_esolver_destroy_f(LIS_ESOLVER_F *esolver, LIS_INT *ierr)
{
	LIS_DEBUG_FUNC_IN;

	*ierr = lis_esolver_destroy((LIS_ESOLVER)LIS_V2P(esolver));
	if( *ierr )	return;

	LIS_DEBUG_FUNC_OUT;
	return;
}
예제 #2
0
int main(int argc, char* argv[])
{
  int               err;
  int               nprocs,mtype,my_rank;
  int               nsol;
  LIS_MATRIX        A,A0;
  LIS_VECTOR        x;
  LIS_REAL          evalue0;
  LIS_SCALAR        shift, *tmpa;
  LIS_ESOLVER       esolver;
  LIS_SOLVER        solver;
  LIS_REAL          residual;
  int               iters;
  double            times;
  double	      itimes,ptimes,p_c_times,p_i_times;
  LIS_PRECON        precon;
  int		      *ptr,*index;
  LIS_SCALAR	      *value;
  int               nsolver;
  char	      esolvername[128];
  int               mode;

  LIS_DEBUG_FUNC_IN;
    
  lis_initialize(&argc, &argv);

#ifdef USE_MPI
  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
  MPI_Comm_rank(MPI_COMM_WORLD,&my_rank);
#else
  nprocs  = 1;
  my_rank = 0;
#endif
    
  if( argc < 4 )
    {
      if( my_rank==0 ) printf("Usage: etest1 matrix_filename solution_filename residual_filename [options]\n");
      lis_finalize();
      exit(0);
    }

  if( my_rank==0 )
    {
      printf("\n");
      printf("number of processes = %d\n",nprocs);
    }

#ifdef _OPENMP
  if( my_rank==0 )
    {
      printf("max number of threads = %d\n",omp_get_num_procs());
      printf("number of threads = %d\n",omp_get_max_threads());
    }
#endif
		
  /* create matrix and vectors */
  lis_matrix_create(LIS_COMM_WORLD,&A);
  lis_input_matrix(A,argv[1]);
  lis_vector_duplicate(A,&x);
  lis_esolver_create(&esolver);
  lis_esolver_set_option("-eprint mem",esolver);
  lis_esolver_set_optionC(esolver);
  lis_esolve(A, x, &evalue0, esolver);
  lis_esolver_get_esolver(esolver,&nsol);
  lis_get_esolvername(nsol,esolvername);
  lis_esolver_get_residualnorm(esolver, &residual);
  lis_esolver_get_iters(esolver, &iters);
  lis_esolver_get_timeex(esolver,&times,&itimes,&ptimes,&p_c_times,&p_i_times);
  if( my_rank==0 ) {
    printf("%s: mode number              = %d\n", esolvername, esolver->options[LIS_EOPTIONS_MODE]);
    printf("%s: eigenvalue               = %e\n", esolvername, evalue0);
    printf("%s: number of iterations     = %d\n",esolvername, iters);
    printf("%s: elapsed time             = %e sec.\n", esolvername, times);
    printf("%s:   preconditioner         = %e sec.\n", esolvername, ptimes);
    printf("%s:     matrix creation      = %e sec.\n", esolvername, p_c_times);
    printf("%s:   linear solver          = %e sec.\n", esolvername, itimes);
    printf("%s: relative residual 2-norm = %e\n\n",esolvername, residual);
  }

  /* write solution */
  lis_output_vector(x,LIS_FMT_MM,argv[2]);

  /* write residual */
  lis_esolver_output_rhistory(esolver, argv[3]);

  lis_esolver_destroy(esolver);
  lis_matrix_destroy(A);
  lis_vector_destroy(x);

  lis_finalize();

  LIS_DEBUG_FUNC_OUT;

  return 0;
}
예제 #3
0
파일: etest4.c 프로젝트: florianl/lis
LIS_INT main(LIS_INT argc, char* argv[])
{
    LIS_INT i,n,gn,is,ie;
    LIS_INT nprocs,my_rank;
    int int_nprocs,int_my_rank;
    LIS_INT nesol;
    LIS_MATRIX A;
    LIS_VECTOR x;
    LIS_REAL evalue0;
    LIS_ESOLVER esolver;
    LIS_REAL residual;
    LIS_INT iter;
    double time;
    double itime,ptime,p_c_time,p_i_time;
    char esolvername[128];
    
    LIS_DEBUG_FUNC_IN;

    lis_initialize(&argc, &argv);

#ifdef USE_MPI
    MPI_Comm_size(MPI_COMM_WORLD,&int_nprocs);
    MPI_Comm_rank(MPI_COMM_WORLD,&int_my_rank);
    nprocs = int_nprocs;
    my_rank = int_my_rank;
#else
    nprocs  = 1;
    my_rank = 0;
#endif
    
    if( argc < 2 )
      {
	if( my_rank==0 ) 
	  {
	      printf("Usage: %s n [eoptions]\n", argv[0]);
	  }
	CHKERR(1);
      }

  if( my_rank==0 )
    {
      printf("\n");
      printf("number of processes = %d\n",nprocs);
    }

#ifdef _OPENMP
  if( my_rank==0 )
    {
#ifdef _LONG__LONG
      printf("max number of threads = %lld\n",omp_get_num_procs());
      printf("number of threads = %lld\n",omp_get_max_threads());
#else
      printf("max number of threads = %d\n",omp_get_num_procs());
      printf("number of threads = %d\n",omp_get_max_threads());
#endif
    }
#endif
		
    /* generate coefficient matrix for one dimensional Poisson equation */
    n = atoi(argv[1]);
    lis_matrix_create(LIS_COMM_WORLD,&A);
    lis_matrix_set_size(A,0,n);
    lis_matrix_get_size(A,&n,&gn);
    lis_matrix_get_range(A,&is,&ie);
    for(i=is;i<ie;i++)
    {
      if( i>0   )  lis_matrix_set_value(LIS_INS_VALUE,i,i-1,-1.0,A);
      if( i<gn-1 ) lis_matrix_set_value(LIS_INS_VALUE,i,i+1,-1.0,A);
      lis_matrix_set_value(LIS_INS_VALUE,i,i,2.0,A);
    }
    lis_matrix_set_type(A,LIS_MATRIX_CSR);
    lis_matrix_assemble(A);
    lis_vector_duplicate(A,&x);

    lis_esolver_create(&esolver);
    lis_esolver_set_option("-eprint mem",esolver);
    lis_esolver_set_optionC(esolver);
    lis_esolve(A, x, &evalue0, esolver);
    lis_esolver_get_esolver(esolver,&nesol);
    lis_esolver_get_esolvername(nesol,esolvername);
    lis_esolver_get_residualnorm(esolver, &residual);
    lis_esolver_get_iter(esolver, &iter);
    lis_esolver_get_timeex(esolver,&time,&itime,&ptime,&p_c_time,&p_i_time);
    if( my_rank==0 ) {
      printf("%s: mode number          = %d\n", esolvername, 0);
#ifdef _LONG__DOUBLE
      printf("%s: eigenvalue           = %Le\n", esolvername, evalue0);
#else
      printf("%s: eigenvalue           = %e\n", esolvername, evalue0);
#endif
#ifdef _LONG__LONG
      printf("%s: number of iterations = %lld\n",esolvername, iter);
#else
      printf("%s: number of iterations = %d\n",esolvername, iter);
#endif
      printf("%s: elapsed time         = %e sec.\n", esolvername, time);
      printf("%s:   preconditioner     = %e sec.\n", esolvername, ptime);
      printf("%s:     matrix creation  = %e sec.\n", esolvername, p_c_time);
      printf("%s:   linear solver      = %e sec.\n", esolvername, itime);
#ifdef _LONG__DOUBLE
      printf("%s: relative residual    = %Le\n\n",esolvername, residual);
#else
      printf("%s: relative residual    = %e\n\n",esolvername, residual);
#endif
  }

    /*
    lis_vector_nrm2(x, &xnrm2);
    lis_vector_scale((1/xnrm2*sqrt(n)), x);
    lis_vector_print(x);
    */

    /*
    lis_vector_create(LIS_COMM_WORLD,&y);
    lis_matrix_create(LIS_COMM_WORLD,&B);
    lis_esolver_get_evalues(esolver,y);
    lis_esolver_get_evectors(esolver,B);
    lis_output_vector(y,LIS_FMT_MM,"evalues.out");
    lis_output_matrix(B,LIS_FMT_MM,"evectors.out");
    lis_vector_destroy(y);
    lis_matrix_destroy(B);
    */

    lis_esolver_destroy(esolver);
    lis_matrix_destroy(A);
    lis_vector_destroy(x);

    lis_finalize();

    LIS_DEBUG_FUNC_OUT;

    return 0;
}
예제 #4
0
파일: getest1.c 프로젝트: anishida/lis
LIS_INT main(int argc, char* argv[])
{
  LIS_Comm comm;
  LIS_INT err;
  int nprocs,my_rank;
  LIS_INT nesol;
  LIS_MATRIX A,B;
  LIS_VECTOR x;
  LIS_SCALAR evalue0;
  LIS_ESOLVER esolver;
  LIS_REAL residual;
  LIS_INT iter;
  double time;
  double itime,ptime,p_c_time,p_i_time;
  char esolvername[128];

  LIS_DEBUG_FUNC_IN;
    
  lis_initialize(&argc, &argv);

  comm = LIS_COMM_WORLD;

#ifdef USE_MPI
  MPI_Comm_size(comm,&nprocs);
  MPI_Comm_rank(comm,&my_rank);
#else
  nprocs  = 1;
  my_rank = 0;
#endif
    
  if( argc < 5 )
    {
      lis_printf(comm,"Usage: %s matrix_a_filename matrix_b_filename evector_filename rhistory_filename [options]\n", argv[0]);
      CHKERR(1);
    }

  lis_printf(comm,"\n");
  lis_printf(comm,"number of processes = %d\n",nprocs);

#ifdef _OPENMP
  lis_printf(comm,"max number of threads = %d\n",omp_get_num_procs());
  lis_printf(comm,"number of threads = %d\n",omp_get_max_threads());
#endif
		
  /* create matrix and vectors */
  lis_matrix_create(comm,&A);  
  lis_matrix_create(comm,&B);
  lis_printf(comm,"\nmatrix A:\n");  
  lis_input_matrix(A,argv[1]);
  lis_printf(comm,"matrix B:\n");    
  lis_input_matrix(B,argv[2]);  
  lis_vector_duplicate(A,&x);
  lis_esolver_create(&esolver);
  lis_esolver_set_option("-e gii -eprint mem",esolver);
  err = lis_esolver_set_optionC(esolver);
  CHKERR(err);      

  err = lis_gesolve(A,B,x,&evalue0,esolver);
  CHKERR(err);
  lis_esolver_get_esolver(esolver,&nesol);
  lis_esolver_get_esolvername(nesol,esolvername);
  lis_esolver_get_residualnorm(esolver,&residual);
  lis_esolver_get_iter(esolver,&iter);
  lis_esolver_get_timeex(esolver,&time,&itime,&ptime,&p_c_time,&p_i_time);

  lis_printf(comm,"%s: mode number          = %d\n", esolvername, 0);
#ifdef _COMPLEX      
  lis_printf(comm,"%s: eigenvalue           = (%e, %e)\n", esolvername, (double)creal(evalue0), (double)cimag(evalue0));
#else
  lis_printf(comm,"%s: eigenvalue           = %e\n", esolvername, (double)evalue0);
#endif      
  lis_printf(comm,"%s: number of iterations = %D\n",esolvername, iter);
  lis_printf(comm,"%s: elapsed time         = %e sec.\n", esolvername, time);
  lis_printf(comm,"%s:   preconditioner     = %e sec.\n", esolvername, ptime);
  lis_printf(comm,"%s:     matrix creation  = %e sec.\n", esolvername, p_c_time);
  lis_printf(comm,"%s:   linear solver      = %e sec.\n", esolvername, itime);
  lis_printf(comm,"%s: relative residual    = %e\n\n",esolvername, (double)residual);

  /* write eigenvector */
  lis_output_vector(x,LIS_FMT_MM,argv[3]);

  /* write residual history */
  lis_esolver_output_rhistory(esolver,argv[4]);

  lis_esolver_destroy(esolver);
  lis_matrix_destroy(A);
  lis_matrix_destroy(B);  
  lis_vector_destroy(x);

  lis_finalize();

  LIS_DEBUG_FUNC_OUT;

  return 0;
}
예제 #5
0
파일: etest5.c 프로젝트: florianl/lis
LIS_INT main(LIS_INT argc, char* argv[])
{
  LIS_INT nprocs,my_rank;
  int int_nprocs,int_my_rank;
  LIS_INT nsol;
  LIS_MATRIX A,B;
  LIS_VECTOR x,y,z,w;
  LIS_SCALAR evalue0;
  LIS_ESOLVER esolver;
  LIS_REAL residual;
  LIS_INT iter;
  double time;
  double itime,ptime,p_c_time,p_i_time;
  char esolvername[128];

  LIS_DEBUG_FUNC_IN;
    
  lis_initialize(&argc, &argv);

#ifdef USE_MPI
  MPI_Comm_size(MPI_COMM_WORLD,&int_nprocs);
  MPI_Comm_rank(MPI_COMM_WORLD,&int_my_rank);
  nprocs = int_nprocs;
  my_rank = int_my_rank;
#else
  nprocs  = 1;
  my_rank = 0;
#endif
    
  if( argc < 6 )
    {
      if( my_rank==0 ) 
	{
	  printf("Usage: %s matrix_filename evalues_filename evectors_filename residuals_filename iters_filename [options]\n", argv[0]);
	}
      CHKERR(1);
    }

  if( my_rank==0 )
    {
      printf("\n");
      printf("number of processes = %d\n",nprocs);
    }

#ifdef _OPENMP
  if( my_rank==0 )
    {
      printf("max number of threads = %d\n",omp_get_num_procs());
      printf("number of threads = %d\n",omp_get_max_threads());
    }
#endif
		
  /* create matrix and vectors */
  lis_matrix_create(LIS_COMM_WORLD,&A);
  lis_input_matrix(A,argv[1]);
  lis_vector_duplicate(A,&x);
  lis_esolver_create(&esolver);
  lis_esolver_set_option("-e si -ss 1 -eprint mem",esolver);
  lis_esolver_set_optionC(esolver);
  lis_esolve(A, x, &evalue0, esolver);
  lis_esolver_get_esolver(esolver,&nsol);
  lis_esolver_get_esolvername(nsol,esolvername);
  lis_esolver_get_residualnorm(esolver, &residual);
  lis_esolver_get_iter(esolver, &iter);
  lis_esolver_get_timeex(esolver,&time,&itime,&ptime,&p_c_time,&p_i_time);
  if( my_rank==0 ) {
    printf("%s: mode number          = %d\n", esolvername, 0);
#ifdef _LONG__DOUBLE
    printf("%s: eigenvalue           = %Le\n", esolvername, evalue0);
#else
    printf("%s: eigenvalue           = %e\n", esolvername, evalue0);
#endif
#ifdef _LONG__LONG
    printf("%s: number of iterations = %lld\n",esolvername, iter);
#else
    printf("%s: number of iterations = %d\n",esolvername, iter);
#endif
    printf("%s: elapsed time         = %e sec.\n", esolvername, time);
    printf("%s:   preconditioner     = %e sec.\n", esolvername, ptime);
    printf("%s:     matrix creation  = %e sec.\n", esolvername, p_c_time);
    printf("%s:   linear solver      = %e sec.\n", esolvername, itime);
#ifdef _LONG__DOUBLE
    printf("%s: relative residual    = %Le\n\n",esolvername, residual);
#else
    printf("%s: relative residual    = %e\n\n",esolvername, residual);
#endif
  }

  lis_vector_create(LIS_COMM_WORLD,&y);
  lis_vector_create(LIS_COMM_WORLD,&z);
  lis_vector_create(LIS_COMM_WORLD,&w);
  lis_matrix_create(LIS_COMM_WORLD,&B);
  lis_esolver_get_evalues(esolver,y);
  lis_esolver_get_residualnorms(esolver,z);
  lis_esolver_get_iters(esolver,w);
  lis_esolver_get_evectors(esolver,B);

  /* write eigenvalues */
  lis_output_vector(y,LIS_FMT_MM,argv[2]);

  /* write eigenvectors */
  lis_output_matrix(B,LIS_FMT_MM,argv[3]);

  /* write residual norms */
  lis_output_vector(z,LIS_FMT_MM,argv[4]);

  /* write numbers of iterations */
  lis_output_vector(w,LIS_FMT_MM,argv[5]);

  lis_esolver_destroy(esolver);
  lis_matrix_destroy(A);
  lis_vector_destroy(x);
  lis_matrix_destroy(B);
  lis_vector_destroy(y);
  lis_vector_destroy(z);
  lis_vector_destroy(w);

  lis_finalize();

  LIS_DEBUG_FUNC_OUT;

  return 0;
}
예제 #6
0
LIS_INT main(LIS_INT argc, char* argv[])
{
    LIS_INT           err,i,n,nnz,is,ie,maxiter,nn,ii,jj;
    LIS_INT           j,k,m;
    LIS_INT           nprocs,mtype,my_rank;
    int               int_nprocs,int_my_rank;
    LIS_INT           nesol;
    LIS_MATRIX        A,A0;
    LIS_VECTOR        x;
    LIS_REAL          evalue0;
    LIS_SCALAR        *tmpa;
    LIS_ESOLVER       esolver;
    LIS_REAL          residual;
    LIS_INT           iters;
    double            times;
    double        itimes,ptimes,p_c_times,p_i_times;
    LIS_INT        *ptr,*index;
    LIS_SCALAR        *value;
    LIS_INT           nesolver;
    char        esolvername[128], solution_filename[128];
    
    LIS_DEBUG_FUNC_IN;
    
    lis_initialize(&argc, &argv);

#ifdef USE_MPI
    MPI_Comm_size(MPI_COMM_WORLD,&int_nprocs);
    MPI_Comm_rank(MPI_COMM_WORLD,&int_my_rank);
    nprocs = int_nprocs;
    my_rank = int_my_rank;
#else
    nprocs  = 1;
    my_rank = 0;
#endif
    
    if( argc < 6 )
      {
  if( my_rank==0 ) 
    {
        printf("Usage: %s m n matrix_type solution_filename residual_filename [options]\n", argv[0]);
    }
  CHKERR(1);
      }

    if( my_rank==0 )
      {
  printf("\n");
#ifdef _LONGLONG
  printf("number of processes = %lld\n",nprocs);
#else
  printf("number of processes = %d\n",nprocs);
#endif
      }

#ifdef _OPENMP
    if( my_rank==0 )
      {
#ifdef _LONGLONG
  printf("max number of threads = %lld\n",omp_get_num_procs());
  printf("number of threads = %lld\n",omp_get_max_threads());
#else
  printf("max number of threads = %d\n",omp_get_num_procs());
  printf("number of threads = %d\n",omp_get_max_threads());
#endif
      }
#endif
    
    /* generate matrix */
    m  = atoi(argv[1]);
    n  = atoi(argv[2]);
    mtype  = atoi(argv[3]);
    if( m<=0 || n<=0 )
      {
#ifdef _LONGLONG
  if( my_rank==0 ) printf("m=%lld <=0 or n=%lld <=0\n",m,n);
#else
  if( my_rank==0 ) printf("m=%d <=0 or n=%d <=0\n",m,n);
#endif
  CHKERR(1);
      }
    nn = m*n;
    err = lis_matrix_create(LIS_COMM_WORLD,&A);
    err = lis_matrix_set_size(A,0,nn);
    CHKERR(err);

    ptr   = (LIS_INT *)malloc((A->n+1)*sizeof(LIS_INT));
    if( ptr==NULL ) CHKERR(1);
    index = (LIS_INT *)malloc(5*A->n*sizeof(LIS_INT));
    if( index==NULL ) CHKERR(1);
    value = (LIS_SCALAR *)malloc(5*A->n*sizeof(LIS_SCALAR));
    if( value==NULL ) CHKERR(1);

    lis_matrix_get_range(A,&is,&ie);
    k = 0;
    for(ii=is;ii<ie;ii++)
      {
  i = ii/m;
  j = ii - i*m;
  if( i>0 )   { jj = ii - m; index[k] = jj; value[k++] = -1.0;}
  if( i<n-1 ) { jj = ii + m; index[k] = jj; value[k++] = -1.0;}
  if( j>0 )   { jj = ii - 1; index[k] = jj; value[k++] = -1.0;}
  if( j<m-1 ) { jj = ii + 1; index[k] = jj; value[k++] = -1.0;}
  index[k] = ii; value[k++] = 4.0;
  ptr[ii-is+1] = k;
      }
    ptr[0] = 0;
    err = lis_matrix_set_csr(ptr[ie-is],ptr,index,value,A);
    CHKERR(err);
    err = lis_matrix_assemble(A);
    CHKERR(err);

    nnz = A->nnz;
#ifdef USE_MPI
    MPI_Allreduce(&nnz,&i,1,LIS_MPI_INT,MPI_SUM,A->comm);
    nnz   = i;
#endif
#ifdef _LONGLONG
    if( my_rank==0 ) printf("matrix size = %lld x %lld (%lld nonzero entries)\n",nn,nn,nnz);
#else
    if( my_rank==0 ) printf("matrix size = %d x %d (%d nonzero entries)\n",nn,nn,nnz);
#endif

    err = lis_matrix_duplicate(A,&A0);
    CHKERR(err);
    lis_matrix_set_type(A0,mtype);
    err = lis_matrix_convert(A,A0);
    CHKERR(err);
    lis_matrix_destroy(A);
    A = A0;
    lis_vector_duplicate(A,&x);

    err = lis_esolver_create(&esolver);
    CHKERR(err);

    lis_esolver_set_option("-eprint mem",esolver);
    lis_esolver_set_optionC(esolver);
    lis_esolve(A, x, &evalue0, esolver);
    lis_esolver_get_esolver(esolver,&nesol);
    lis_esolver_get_esolvername(nesol,esolvername);
    lis_esolver_get_residualnorm(esolver, &residual);
    lis_esolver_get_iters(esolver, &iters);
    lis_esolver_get_timeex(esolver,&times,&itimes,&ptimes,&p_c_times,&p_i_times);
    if( my_rank==0 ) {
      printf("%s: mode number              = %d\n", esolvername, 0);
#ifdef _LONG__DOUBLE
      printf("%s: eigenvalue               = %Le\n", esolvername, evalue0);
#else
      printf("%s: eigenvalue               = %e\n", esolvername, evalue0);
#endif
#ifdef _LONGLONG
      printf("%s: number of iterations     = %lld\n",esolvername, iters);
#else
      printf("%s: number of iterations     = %d\n",esolvername, iters);
#endif
      printf("%s: elapsed time             = %e sec.\n", esolvername, times);
      printf("%s:   preconditioner         = %e sec.\n", esolvername, ptimes);
      printf("%s:     matrix creation      = %e sec.\n", esolvername, p_c_times);
      printf("%s:   linear solver          = %e sec.\n", esolvername, itimes);
#ifdef _LONG__DOUBLE
      printf("%s: relative residual 2-norm = %Le\n\n",esolvername, residual);
#else
      printf("%s: relative residual 2-norm = %e\n\n",esolvername, residual);
#endif
  }

    /* write solution */
    lis_output_vector(x,LIS_FMT_MM,argv[4]);

    /* write residual */
    lis_esolver_output_rhistory(esolver, argv[5]);

    CHKERR(err);

    lis_esolver_destroy(esolver);
    lis_matrix_destroy(A);
    lis_matrix_destroy(A0);
    lis_vector_destroy(x);

    lis_finalize();

    LIS_DEBUG_FUNC_OUT;

    return 0;
}
예제 #7
0
LIS_INT lis_eai_quad(LIS_ESOLVER esolver)
{
  LIS_MATRIX A;
  LIS_INT ss,ic;
  LIS_INT emaxiter,iter0,hqriter;
  LIS_REAL tol,hqrerr,D;
  LIS_INT i,j;
  LIS_INT output, niesolver;
  LIS_REAL nrm2,resid0; 
  LIS_VECTOR *v,w;
  LIS_SCALAR *h,*hq,*hr,evalue,evalue0;
  LIS_SOLVER solver;
  LIS_ESOLVER esolver2;
  char esolvername[128],solvername[128],preconname[128];
  LIS_INT nsol,precon_type;

  ss = esolver->options[LIS_EOPTIONS_SUBSPACE];
  emaxiter = esolver->options[LIS_EOPTIONS_MAXITER];
  tol = esolver->params[LIS_EPARAMS_RESID - LIS_EOPTIONS_LEN]; 
  output  = esolver->options[LIS_EOPTIONS_OUTPUT];
  niesolver = esolver->options[LIS_EOPTIONS_INNER_ESOLVER];

  h = (LIS_SCALAR *)lis_malloc(ss*ss*sizeof(LIS_SCALAR), "lis_eai_quad::h");
  hq = (LIS_SCALAR *)lis_malloc(ss*ss*sizeof(LIS_SCALAR), "lis_eai_quad::hq");
  hr = (LIS_SCALAR *)lis_malloc(ss*ss*sizeof(LIS_SCALAR), "lis_eai_quad::hr");
  
  A = esolver->A;
  w = esolver->work[0];
  v = &esolver->work[1];
  lis_vector_set_all(0.0,v[0]);
  lis_vector_set_all(1.0,w);
  lis_vector_nrm2(w, &nrm2);

  lis_solver_create(&solver);
  lis_solver_set_option("-i bicg -p none",solver);  
  lis_solver_set_optionC(solver);
  lis_solver_get_solver(solver, &nsol);
  lis_solver_get_precon(solver, &precon_type);
  lis_solver_get_solvername(nsol, solvername);
  lis_solver_get_preconname(precon_type, preconname);
  lis_esolver_get_esolvername(niesolver, esolvername);
  if( A->my_rank==0 ) printf("inner eigensolver     : %s\n", esolvername);
  if( A->my_rank==0 ) printf("linear solver         : %s\n", solvername);
  if( A->my_rank==0 ) printf("preconditioner        : %s\n", preconname);

  for (i=0;i<ss*ss;i++) h[i] = 0.0;

  j=-1;
  while (j<ss-1)
    {
      j = j+1;
      lis_vector_copy(w, v[j]);

      /* w = A * v(j) */
      lis_matvec(A, v[j], w);

      /* reorthogonalization */
      for (i=0;i<=j;i++)
	{
	  /* h(i,j) = <v(i), w> */
	  lis_vector_dot(v[i], w, &h[i+j*ss]);
	  /* w = w - h(i,j) * v(i) */
	  lis_vector_axpy(-h[i+j*ss], v[i], w); 
	}

      /* h(j+1,j) = ||w||_2 */
      lis_vector_nrm2(w, &h[j+1+j*ss]);

      /* convergence check */
      if (fabs(h[j+1+j*ss])<tol) break;

      /* v(j+1) = w / h(i+1,j) */
      lis_vector_scale(1/h[j+1+j*ss],w);
      lis_vector_copy(w,v[j+1]);
      
    }

  /* compute eigenvalues of a real upper
     Hessenberg matrix H(j) = SH'(j)S^* */
  lis_array_qr(ss,h,hq,hr,&hqriter,&hqrerr);


  if( A->my_rank==0 ) 
    {
#ifdef _LONG__LONG
      printf("size of subspace      : %lld\n\n", ss);
#else
      printf("size of subspace      : %d\n\n", ss);
#endif
      if( output ) printf("approximate eigenvalues in subspace:\n\n");


      i=0;
      while (i<ss-1) 
	{
	  i = i + 1;
	  if (fabs(h[i+(i-1)*ss])<tol)
	    {
#ifdef _LONG__LONG
	      printf("Arnoldi: mode number              = %lld\n",i-1);
#else
	      printf("Arnoldi: mode number              = %d\n",i-1);
#endif	  
#ifdef _LONG__DOUBLE
	      printf("Arnoldi: eigenvalue               = %Le\n",h[i-1+(i-1)*ss]);
#else
	      printf("Arnoldi: eigenvalue               = %e\n",h[i-1+(i-1)*ss]);
#endif
	      esolver->evalue[i-1] = h[i-1+(i-1)*ss];
	    }
	  else
	    {
	      D = (h[i-1+(i-1)*ss]-h[i+i*ss]) * (h[i-1+(i-1)*ss]-h[i+i*ss])
		+ 4 * h[i-1+i*ss] * h[i+(i-1)*ss];
	      if (D<0)
		{
#ifdef _LONG__LONG
		  printf("Arnoldi: mode number              = %lld\n",i-1);
#else
		  printf("Arnoldi: mode number              = %d\n",i-1);
#endif
#ifdef _LONG__DOUBLE	      
		  printf("Arnoldi: eigenvalue               = %Le + %Le i\n", (h[i-1+(i-1)*ss]+h[i+i*ss])/2, sqrt(-D)/2);
#else
		  printf("Arnoldi: eigenvalue               = %e + %e i\n", (h[i-1+(i-1)*ss]+h[i+i*ss])/2, sqrt(-D)/2);
#endif
#ifdef _LONG__LONG	      
		  printf("Arnoldi: mode number              = %lld\n",i);
#else
		  printf("Arnoldi: mode number              = %d\n",i);
#endif
#ifdef _LONG__DOUBLE	      	      
		  printf("Arnoldi: eigenvalue               = %Le - %Le i\n", (h[i-1+(i-1)*ss]+h[i+i*ss])/2, sqrt(-D)/2);
#else
		  printf("Arnoldi: eigenvalue               = %e - %e i\n", (h[i-1+(i-1)*ss]+h[i+i*ss])/2, sqrt(-D)/2);
#endif	      
		  esolver->evalue[i-1] = (h[i-1+(i-1)*ss]+h[i+i*ss])/2;
		  esolver->evalue[i]   = (h[i-1+(i-1)*ss]+h[i+i*ss])/2;	      
		  i=i+1;
		}
	      else
		{
#ifdef _LONG__LONG	      	      
		  printf("Arnoldi: mode number              = %lld\n",i-1);
#else
		  printf("Arnoldi: mode number              = %d\n",i-1);
#endif
#ifdef _LONG__DOUBLE	      	      	      
		  printf("Arnoldi: eigenvalue               = %Le\n",h[i-1+(i-1)*ss]);
#else
		  printf("Arnoldi: eigenvalue               = %e\n",h[i-1+(i-1)*ss]);
#endif	      
		  esolver->evalue[i-1] = h[i-1+(i-1)*ss];
		}
	    }
	}
      if (i<ss)
	{
#ifdef _LONG__LONG	            
	  printf("Arnoldi: mode number              = %lld\n",i);
#else
	  printf("Arnoldi: mode number              = %d\n",i);
#endif
#ifdef _LONG__DOUBLE	      	      	      	      
	  printf("Arnoldi: eigenvalue               = %Le\n",h[i+i*ss]);
#else
	  printf("Arnoldi: eigenvalue               = %e\n",h[i+i*ss]);
#endif	      
	}

      if( output ) printf("\n");
      if( output ) printf("compute refined (real) eigenpairs, where imaginary parts are currently neglected:\n\n");
  
    }

  lis_esolver_create(&esolver2);
  esolver2->options[LIS_EOPTIONS_ESOLVER] = niesolver;
  esolver2->options[LIS_EOPTIONS_SUBSPACE] = 1;
  esolver2->options[LIS_EOPTIONS_MAXITER] = emaxiter;
  esolver2->options[LIS_EOPTIONS_OUTPUT] = esolver->options[LIS_EOPTIONS_OUTPUT];
  esolver2->params[LIS_EPARAMS_RESID - LIS_EOPTIONS_LEN] = tol;
  esolver2->eprecision = LIS_PRECISION_QUAD;

  /* compute refined (real) eigenpairs, where imaginary parts are currently neglected */

  for (i=0;i<ss;i++)
    {
      lis_vector_duplicate(A, &esolver->evector[i]); 
      esolver2->lshift = -(esolver->evalue[i]);
      lis_esolve(A, esolver->evector[i], &evalue, esolver2);
      lis_esolver_work_destroy(esolver2); 
      esolver->evalue[i] = evalue - esolver2->lshift;
      esolver->iter[i] = esolver2->iter[0];            
      esolver->resid[i] = esolver2->resid[0];

      if (i==0) 
	{
	  evalue0 = esolver->evalue[0];
	  iter0 = esolver2->iter[0];
	  resid0 = esolver2->resid[0];
	  if( output & LIS_EPRINT_MEM ) 
	    {
	      for (ic=0;ic<iter0+1;ic++)
		{
		  esolver->rhistory[ic] = esolver2->rhistory[ic]; 
		}
	    }
	  esolver->ptime = esolver2->ptime;
	  esolver->itime = esolver2->itime;
	  esolver->p_c_time = esolver2->p_c_time;
	  esolver->p_i_time = esolver2->p_i_time;
	}

      if (A->my_rank==0) 
	{

#ifdef _LONG__LONG
	  if( output ) printf("Arnoldi: mode number          = %lld\n", i);
#else
	  if( output ) printf("Arnoldi: mode number          = %d\n", i);
#endif
#ifdef _LONG__DOUBLE
	  if( output ) printf("Arnoldi: eigenvalue           = %Le\n", esolver->evalue[i]);
#else
	  if( output ) printf("Arnoldi: eigenvalue           = %e\n", esolver->evalue[i]);
#endif
#ifdef _LONG__LONG
	  if( output ) printf("Arnoldi: number of iterations = %lld\n",esolver2->iter[0]);
#else
	  if( output ) printf("Arnoldi: number of iterations = %d\n",esolver2->iter[0]);
#endif
#ifdef _LONG__DOUBLE
	  if( output ) printf("Arnoldi: relative residual    = %Le\n\n",esolver2->resid[0]);
#else
	  if( output ) printf("Arnoldi: relative residual    = %e\n\n",esolver2->resid[0]);
#endif
	}
    }

  esolver->evalue[0] = evalue0; 
  esolver->iter[0] = iter0;
  esolver->resid[0] = resid0;

  lis_vector_copy(esolver->evector[0], esolver->x);

  lis_esolver_destroy(esolver2); 

  lis_free(h); 
  lis_free(hq);
  lis_free(hr);

  lis_solver_destroy(solver);

  LIS_DEBUG_FUNC_OUT;
  return LIS_SUCCESS;
}