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 #include <stdio.h>
36 #include <stdlib.h>
37 #include <math.h>
38 #include "lis.h"
39 
40 #undef __FUNC__
41 #define __FUNC__ "main"
main(int argc,char * argv[])42 LIS_INT main(int argc, char* argv[])
43 {
44 	LIS_Comm comm;
45 	LIS_MATRIX A0,A;
46 	LIS_VECTOR x,b,u;
47 	LIS_SOLVER solver;
48 	LIS_INT m,n,nn,nnz;
49 	LIS_INT	i,j,ii,jj,ctr;
50 	LIS_INT	is,ie;
51 	int nprocs,my_rank;
52 	LIS_INT	nsol;
53 	LIS_INT	err,iter,mtype,iter_double,iter_quad;
54 	double time,itime,ptime,p_c_time,p_i_time;
55 	LIS_REAL resid;
56 	char solvername[128];
57 	LIS_INT	*ptr,*index;
58 	LIS_SCALAR *value;
59 
60 
61 	LIS_DEBUG_FUNC_IN;
62 
63 	lis_initialize(&argc, &argv);
64 
65 	comm = LIS_COMM_WORLD;
66 
67 	#ifdef USE_MPI
68 	        MPI_Comm_size(comm,&nprocs);
69 		MPI_Comm_rank(comm,&my_rank);
70 	#else
71 		nprocs  = 1;
72 		my_rank = 0;
73 	#endif
74 
75 	if( argc < 6 )
76 	{
77 	  lis_printf(comm,"Usage: %s m n matrix_type solution_filename rhistory_filename [options]\n", argv[0]);
78 	  CHKERR(1);
79 	}
80 
81 	m  = atoi(argv[1]);
82 	n  = atoi(argv[2]);
83 	mtype = atoi(argv[3]);
84 	if( m<=0 || n<=0 )
85 	{
86 	  lis_printf(comm,"m=%D <=0 or n=%D <=0\n",m,n);
87 	  CHKERR(1);
88 	}
89 
90 	lis_printf(comm,"\n");
91 	lis_printf(comm,"number of processes = %d\n",nprocs);
92 
93 #ifdef _OPENMP
94 	lis_printf(comm,"max number of threads = %d\n",omp_get_num_procs());
95 	lis_printf(comm,"number of threads = %d\n",omp_get_max_threads());
96 #endif
97 
98 	/* create matrix and vectors */
99 	nn = m*n;
100 	err = lis_matrix_create(comm,&A);
101 	err = lis_matrix_set_size(A,0,nn);
102 	CHKERR(err);
103 
104 	ptr = (LIS_INT *)malloc((A->n+1)*sizeof(LIS_INT));
105 	if( ptr==NULL ) CHKERR(1);
106 	index = (LIS_INT *)malloc(5*A->n*sizeof(LIS_INT));
107 	if( index==NULL ) CHKERR(1);
108 	value = (LIS_SCALAR *)malloc(5*A->n*sizeof(LIS_SCALAR));
109 	if( value==NULL ) CHKERR(1);
110 
111 	lis_matrix_get_range(A,&is,&ie);
112 	ctr = 0;
113 	for(ii=is;ii<ie;ii++)
114 	{
115 		i = ii/m;
116 		j = ii - i*m;
117 		if( i>0 )   { jj = ii - m; index[ctr] = jj; value[ctr++] = -1.0;}
118 		if( i<n-1 ) { jj = ii + m; index[ctr] = jj; value[ctr++] = -1.0;}
119 		if( j>0 )   { jj = ii - 1; index[ctr] = jj; value[ctr++] = -1.0;}
120 		if( j<m-1 ) { jj = ii + 1; index[ctr] = jj; value[ctr++] = -1.0;}
121 		index[ctr] = ii; value[ctr++] = 4.0;
122 		ptr[ii-is+1] = ctr;
123 	}
124 	ptr[0] = 0;
125 	err = lis_matrix_set_csr(ptr[ie-is],ptr,index,value,A);
126 	CHKERR(err);
127 	err = lis_matrix_assemble(A);
128 	CHKERR(err);
129 
130 	nnz = A->nnz;
131 #ifdef USE_MPI
132 	MPI_Allreduce(&nnz,&i,1,LIS_MPI_INT,MPI_SUM,A->comm);
133 	nnz = i;
134 #endif
135 
136 	lis_printf(comm,"matrix size = %D x %D (%D nonzero entries)\n\n",nn,nn,nnz);
137 
138 	err = lis_matrix_duplicate(A,&A0);
139 	CHKERR(err);
140 	lis_matrix_set_type(A0,mtype);
141 	err = lis_matrix_convert(A,A0);
142 	CHKERR(err);
143 	lis_matrix_destroy(A);
144 	A = A0;
145 
146 	/*
147         err = lis_vector_create(comm,&u);
148 	CHKERR(err);
149 	err = lis_vector_set_size(u,0,nn);
150 	CHKERR(err);
151 	*/
152 
153 	err = lis_vector_duplicate(A,&u);
154 	CHKERR(err);
155 	err = lis_vector_duplicate(A,&b);
156 	CHKERR(err);
157 	err = lis_vector_duplicate(A,&x);
158 	CHKERR(err);
159 
160 	err = lis_vector_set_all(1.0,u);
161 	lis_matvec(A,u,b);
162 
163 	err = lis_solver_create(&solver); CHKERR(err);
164 	lis_solver_set_option("-print mem",solver);
165 	err = lis_solver_set_optionC(solver);
166 	CHKERR(err);
167 
168 	err = lis_solve(A,b,x,solver);
169 	CHKERR(err);
170 	lis_solver_get_iterex(solver,&iter,&iter_double,&iter_quad);
171 	lis_solver_get_timeex(solver,&time,&itime,&ptime,&p_c_time,&p_i_time);
172 	lis_solver_get_residualnorm(solver,&resid);
173 	lis_solver_get_solver(solver,&nsol);
174 	lis_solver_get_solvername(nsol,solvername);
175 
176 	lis_printf(comm,"%s: number of iterations = %D\n",solvername,iter);
177 #ifndef _LONG__DOUBLE
178 	lis_printf(comm,"%s:   double             = %D\n",solvername,iter_double);
179 	lis_printf(comm,"%s:   quad               = %D\n",solvername,iter_quad);
180 #endif
181 	lis_printf(comm,"%s: elapsed time         = %e sec.\n",solvername,time);
182 	lis_printf(comm,"%s:   preconditioner     = %e sec.\n",solvername, ptime);
183 	lis_printf(comm,"%s:     matrix creation  = %e sec.\n",solvername, p_c_time);
184 	lis_printf(comm,"%s:   linear solver      = %e sec.\n",solvername, itime);
185 	lis_printf(comm,"%s: relative residual    = %e\n\n",solvername,(double)resid);
186 
187 	/* write solution */
188 	lis_output_vector(x,LIS_FMT_MM,argv[4]);
189 
190 	/* write residual history */
191 	lis_solver_output_rhistory(solver, argv[5]);
192 
193 	lis_solver_destroy(solver);
194 	lis_matrix_destroy(A);
195 	lis_vector_destroy(b);
196 	lis_vector_destroy(x);
197 	lis_vector_destroy(u);
198 
199 	lis_finalize();
200 
201 	LIS_DEBUG_FUNC_OUT;
202 
203 	return 0;
204 }
205 
206 
207