1 /* Copyright (C) 2005 The Scalable Software Infrastructure Project. All rights reserved.
2
3 Redistribution and use in source and binary forms, with or without
4 modification, are permitted provided that the following conditions are met:
5 1. Redistributions of source code must retain the above copyright
6 notice, this list of conditions and the following disclaimer.
7 2. Redistributions in binary form must reproduce the above copyright
8 notice, this list of conditions and the following disclaimer in the
9 documentation and/or other materials provided with the distribution.
10 3. Neither the name of the project nor the names of its contributors
11 may be used to endorse or promote products derived from this software
12 without specific prior written permission.
13
14 THIS SOFTWARE IS PROVIDED BY THE SCALABLE SOFTWARE INFRASTRUCTURE PROJECT
15 ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
16 TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
17 PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE SCALABLE SOFTWARE INFRASTRUCTURE
18 PROJECT BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
19 OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
20 SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
21 INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
22 CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
23 ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
24 POSSIBILITY OF SUCH DAMAGE.
25 */
26
27 #ifdef HAVE_CONFIG_H
28 #include "lis_config.h"
29 #else
30 #ifdef HAVE_CONFIG_WIN_H
31 #include "lis_config_win.h"
32 #endif
33 #endif
34
35 #ifdef HAVE_CONFIG_H
36 #include "lis_config.h"
37 #else
38 #ifdef HAVE_CONFIG_WIN_H
39 #include "lis_config_win.h"
40 #endif
41 #endif
42
43 #include <stdio.h>
44 #include <stdlib.h>
45 #include <math.h>
46 #include "lis.h"
47
48
49 #undef __FUNC__
50 #define __FUNC__ "main"
main(int argc,char * argv[])51 LIS_INT main(int argc, char* argv[])
52 {
53 LIS_Comm comm;
54 LIS_MATRIX A;
55 LIS_VECTOR x,b,u;
56 LIS_SOLVER solver;
57 LIS_INT k,n,gn,ii,jj;
58 LIS_INT is,ie;
59 int nprocs,my_rank;
60 LIS_INT nsol;
61 LIS_INT err,iter,iter_double,iter_quad;
62 double time,itime,ptime,p_c_time,p_i_time;
63 LIS_REAL resid;
64 char solvername[128];
65 LIS_INT *ptr,*index;
66 LIS_SCALAR *value,gamma;
67
68
69 LIS_DEBUG_FUNC_IN;
70
71 lis_initialize(&argc, &argv);
72
73 comm = LIS_COMM_WORLD;
74
75 #ifdef USE_MPI
76 MPI_Comm_size(comm,&nprocs);
77 MPI_Comm_rank(comm,&my_rank);
78 #else
79 nprocs = 1;
80 my_rank = 0;
81 #endif
82
83 if( argc < 3 )
84 {
85 lis_printf(comm,"Usage: %s n gamma [options]\n", argv[0]);
86 CHKERR(1);
87 }
88
89 gn = atoi(argv[1]);
90 gamma = atof(argv[2]);
91 if( gn<=0 )
92 {
93 lis_printf(comm,"n=%D <=0 \n",gn);
94 CHKERR(1);
95 }
96
97 lis_printf(comm,"\n");
98 lis_printf(comm,"number of processes = %d\n",nprocs);
99
100 #ifdef _OPENMP
101 lis_printf(comm,"max number of threads = %d\n",omp_get_num_procs());
102 lis_printf(comm,"number of threads = %d\n",omp_get_max_threads());
103 #endif
104
105 lis_printf(comm,"n = %D, gamma = %f\n\n",gn,(double)gamma);
106
107 /* create matrix and vectors */
108 err = lis_matrix_create(comm,&A); CHKERR(err);
109 err = lis_matrix_set_size(A,0,gn); CHKERR(err);
110 err = lis_matrix_get_size(A,&n,&gn); CHKERR(err);
111 err = lis_matrix_malloc_csr(n,3*n,&ptr,&index,&value); CHKERR(err);
112 err = lis_matrix_get_range(A,&is,&ie); CHKERR(err);
113
114 k = 0;
115 ptr[0] = 0;
116 for(ii=is;ii<ie;ii++)
117 {
118 if( ii>1 ) { jj = ii - 2; index[k] = jj; value[k++] = gamma;}
119 if( ii<gn-1 ) { jj = ii + 1; index[k] = jj; value[k++] = 1.0;}
120 index[k] = ii; value[k++] = 2.0;
121 ptr[ii-is+1] = k;
122 }
123 err = lis_matrix_set_csr(ptr[ie-is],ptr,index,value,A); CHKERR(err);
124 err = lis_matrix_assemble(A); CHKERR(err);
125
126 err = lis_vector_duplicate(A,&u); CHKERR(err);
127 err = lis_vector_duplicate(u,&b); CHKERR(err);
128 err = lis_vector_duplicate(u,&x); CHKERR(err);
129
130 err = lis_vector_set_all(1.0,u);
131 lis_matvec(A,u,b);
132
133 err = lis_solver_create(&solver); CHKERR(err);
134 lis_solver_set_option("-print mem",solver);
135 err = lis_solver_set_optionC(solver);
136 CHKERR(err);
137
138 err = lis_solve(A,b,x,solver); CHKERR(err);
139 lis_solver_get_iterex(solver,&iter,&iter_double,&iter_quad);
140 lis_solver_get_timeex(solver,&time,&itime,&ptime,&p_c_time,&p_i_time);
141 lis_solver_get_residualnorm(solver,&resid);
142 lis_solver_get_solver(solver,&nsol);
143 lis_solver_get_solvername(nsol,solvername);
144
145 lis_printf(comm,"%s: number of iterations = %D\n",solvername,iter);
146 #ifndef _LONG__DOUBLE
147 lis_printf(comm,"%s: double = %D\n",solvername,iter_double);
148 lis_printf(comm,"%s: quad = %D\n",solvername,iter_quad);
149 #endif
150 lis_printf(comm,"%s: elapsed time = %e sec.\n",solvername,time);
151 lis_printf(comm,"%s: preconditioner = %e sec.\n",solvername, ptime);
152 lis_printf(comm,"%s: matrix creation = %e sec.\n",solvername, p_c_time);
153 lis_printf(comm,"%s: linear solver = %e sec.\n",solvername, itime);
154 lis_printf(comm,"%s: relative residual = %e\n\n",solvername,(double)resid);
155
156 lis_solver_destroy(solver);
157 lis_matrix_destroy(A);
158 lis_vector_destroy(b);
159 lis_vector_destroy(x);
160 lis_vector_destroy(u);
161
162 lis_finalize();
163
164 LIS_DEBUG_FUNC_OUT;
165
166 return 0;
167 }
168
169
170