Line data Source code
1 : /* This file is part of ELPA. */
2 : /* */
3 : /* The ELPA library was originally created by the ELPA consortium, */
4 : /* consisting of the following organizations: */
5 : /* */
6 : /* - Max Planck Computing and Data Facility (MPCDF), formerly known as */
7 : /* Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG), */
8 : /* - Bergische Universität Wuppertal, Lehrstuhl für angewandte */
9 : /* Informatik, */
10 : /* - Technische Universität München, Lehrstuhl für Informatik mit */
11 : /* Schwerpunkt Wissenschaftliches Rechnen , */
12 : /* - Fritz-Haber-Institut, Berlin, Abt. Theorie, */
13 : /* - Max-Plack-Institut für Mathematik in den Naturwissenschaften, */
14 : /* Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition, */
15 : /* and */
16 : /* - IBM Deutschland GmbH */
17 : /* */
18 : /* */
19 : /* More information can be found here: */
20 : /* http://elpa.mpcdf.mpg.de/ */
21 : /* */
22 : /* ELPA is free software: you can redistribute it and/or modify */
23 : /* it under the terms of the version 3 of the license of the */
24 : /* GNU Lesser General Public License as published by the Free */
25 : /* Software Foundation. */
26 : /* */
27 : /* ELPA is distributed in the hope that it will be useful, */
28 : /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
29 : /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
30 : /* GNU Lesser General Public License for more details. */
31 : /* */
32 : /* You should have received a copy of the GNU Lesser General Public License */
33 : /* along with ELPA. If not, see <http://www.gnu.org/licenses/> */
34 : /* */
35 : /* ELPA reflects a substantial effort on the part of the original */
36 : /* ELPA consortium, and we ask you to respect the spirit of the */
37 : /* license that we chose: i.e., please contribute any changes you */
38 : /* may have back to the original ELPA library distribution, and keep */
39 : /* any derivatives of ELPA under the same license that we chose for */
40 : /* the original distribution, the GNU Lesser General Public License. */
41 : /* */
42 : /* */
43 :
44 : #include "config-f90.h"
45 :
46 : #include <stdio.h>
47 : #include <stdlib.h>
48 : #ifdef WITH_MPI
49 : #include <mpi.h>
50 : #endif
51 : #include <math.h>
52 :
53 : #include <elpa/elpa_legacy.h>
54 : #include <test/shared/generated.h>
55 : #include <complex.h>
56 :
57 : #define DOUBLE_PRECISION_COMPLEX 1
58 :
59 192 : int main(int argc, char** argv) {
60 : int myid;
61 : int nprocs;
62 : #ifndef WITH_MPI
63 : int MPI_COMM_WORLD;
64 : #endif
65 : int na, nev, nblk;
66 :
67 : int status;
68 :
69 : int np_cols, np_rows, np_colsStart;
70 :
71 : int my_blacs_ctxt, my_prow, my_pcol;
72 :
73 : int mpierr;
74 :
75 : int my_mpi_comm_world;
76 : int mpi_comm_rows, mpi_comm_cols;
77 :
78 : int info, *sc_desc;
79 :
80 : int na_rows, na_cols;
81 :
82 : double startVal;
83 : #ifdef DOUBLE_PRECISION_COMPLEX
84 :
85 : complex double *a, *z, *as;
86 : double *ev;
87 : #else
88 :
89 : complex *a, *z, *as;
90 : float *ev;
91 : #endif
92 :
93 : int useGPU;
94 : int success;
95 :
96 : #ifdef WITH_MPI
97 128 : MPI_Init(&argc, &argv);
98 128 : MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
99 128 : MPI_Comm_rank(MPI_COMM_WORLD, &myid);
100 : #else
101 64 : nprocs=1;
102 64 : myid=0;
103 64 : MPI_COMM_WORLD=1;
104 : #endif
105 192 : na = 1000;
106 192 : nev = 500;
107 192 : nblk = 16;
108 :
109 192 : if (myid == 0) {
110 128 : printf("This is the c version of an ELPA test-programm\n");
111 128 : printf("\n");
112 128 : printf("It will call the 1stage ELPA complex solver for a matrix\n");
113 128 : printf("of matrix size %d. It will compute %d eigenvalues\n",na,nev);
114 128 : printf("and uses a blocksize of %d\n",nblk);
115 128 : printf("\n");
116 128 : printf("This is an example program with much less functionality\n");
117 128 : printf("as it's Fortran counterpart. It's only purpose is to show how \n");
118 128 : printf("to evoke ELPA1 from a c programm\n");
119 :
120 128 : printf("\n");
121 :
122 : #ifdef DOUBLE_PRECISION_COMPLEX
123 128 : printf("The double precision version of ELPA1 is used\n");
124 : #else
125 : printf("The single precision version of ELPA1 is used\n");
126 : #endif
127 128 : printf("\n");
128 : }
129 :
130 192 : status = 0;
131 :
132 192 : startVal = sqrt((double) nprocs);
133 192 : np_colsStart = (int) round(startVal);
134 192 : for (np_cols=np_colsStart;np_cols>1;np_cols--){
135 0 : if (nprocs %np_cols ==0){
136 0 : break;
137 : }
138 : }
139 :
140 192 : np_rows = nprocs/np_cols;
141 :
142 192 : if (myid == 0) {
143 128 : printf("\n");
144 128 : printf("Number of processor rows %d, cols %d, total %d \n",np_rows,np_cols,nprocs);
145 : }
146 :
147 : /* set up blacs */
148 : /* convert communicators before */
149 : #ifdef WITH_MPI
150 128 : my_mpi_comm_world = MPI_Comm_c2f(MPI_COMM_WORLD);
151 : #else
152 64 : my_mpi_comm_world = 1;
153 : #endif
154 192 : set_up_blacsgrid_f(my_mpi_comm_world, np_rows, np_cols, 'C', &my_blacs_ctxt, &my_prow, &my_pcol);
155 :
156 192 : if (myid == 0) {
157 128 : printf("\n");
158 128 : printf("Past BLACS_Gridinfo...\n");
159 128 : printf("\n");
160 : }
161 :
162 : /* get the ELPA row and col communicators. */
163 : /* These are NOT usable in C without calling the MPI_Comm_f2c function on them !! */
164 : #ifdef WITH_MPI
165 128 : my_mpi_comm_world = MPI_Comm_c2f(MPI_COMM_WORLD);
166 : #endif
167 192 : mpierr = elpa_get_communicators(my_mpi_comm_world, my_prow, my_pcol, &mpi_comm_rows, &mpi_comm_cols);
168 :
169 192 : if (myid == 0) {
170 128 : printf("\n");
171 128 : printf("Past split communicator setup for rows and columns...\n");
172 128 : printf("\n");
173 : }
174 :
175 192 : sc_desc = malloc(9*sizeof(int));
176 :
177 192 : set_up_blacs_descriptor_f(na, nblk, my_prow, my_pcol, np_rows, np_cols, &na_rows, &na_cols, sc_desc, my_blacs_ctxt, &info);
178 :
179 192 : if (myid == 0) {
180 128 : printf("\n");
181 128 : printf("Past scalapack descriptor setup...\n");
182 128 : printf("\n");
183 : }
184 :
185 : /* allocate the matrices needed for elpa */
186 192 : if (myid == 0) {
187 128 : printf("\n");
188 128 : printf("Allocating matrices with na_rows=%d and na_cols=%d\n",na_rows, na_cols);
189 128 : printf("\n");
190 : }
191 :
192 : #ifdef DOUBLE_PRECISION_COMPLEX
193 192 : a = malloc(na_rows*na_cols*sizeof(complex double));
194 192 : z = malloc(na_rows*na_cols*sizeof(complex double));
195 192 : as = malloc(na_rows*na_cols*sizeof(complex double));
196 192 : ev = malloc(na*sizeof(double));
197 : #else
198 : a = malloc(na_rows*na_cols*sizeof(complex));
199 : z = malloc(na_rows*na_cols*sizeof(complex));
200 : as = malloc(na_rows*na_cols*sizeof(complex));
201 : ev = malloc(na*sizeof(float));
202 : #endif
203 :
204 : #ifdef DOUBLE_PRECISION_COMPLEX
205 192 : prepare_matrix_random_complex_double_f(na, myid, na_rows, na_cols, sc_desc, a, z, as);
206 : #else
207 : prepare_matrix_random_complex_single_f(na, myid, na_rows, na_cols, sc_desc, a, z, as);
208 : #endif
209 :
210 192 : if (myid == 0) {
211 128 : printf("\n");
212 128 : printf("Entering ELPA 1stage complex solver\n");
213 128 : printf("\n");
214 : }
215 : #ifdef WITH_MPI
216 128 : mpierr = MPI_Barrier(MPI_COMM_WORLD);
217 : #endif
218 :
219 192 : useGPU = 0;
220 : #ifdef DOUBLE_PRECISION_COMPLEX
221 192 : success = elpa_solve_evp_complex_1stage_double_precision(na, nev, a, na_rows, ev, z, na_rows, nblk, na_cols, mpi_comm_rows, mpi_comm_cols, my_mpi_comm_world, useGPU);
222 : #else
223 : success = elpa_solve_evp_complex_1stage_single_precision(na, nev, a, na_rows, ev, z, na_rows, nblk, na_cols, mpi_comm_rows, mpi_comm_cols, my_mpi_comm_world, useGPU);
224 : #endif
225 :
226 192 : if (success != 1) {
227 0 : printf("error in ELPA solve \n");
228 : #ifdef WITH_MPI
229 0 : mpierr = MPI_Abort(MPI_COMM_WORLD, 99);
230 : #else
231 0 : exit(99);
232 : #endif
233 : }
234 :
235 :
236 192 : if (myid == 0) {
237 128 : printf("\n");
238 128 : printf("1stage ELPA complex solver complete\n");
239 128 : printf("\n");
240 : }
241 :
242 : /* check the results */
243 : #ifdef DOUBLE_PRECISION_COMPLEX
244 192 : status = check_correctness_evp_numeric_residuals_complex_double_f(na, nev, na_rows, na_cols, as, z, ev, sc_desc, myid);
245 : #else
246 : status = check_correctness_evp_numeric_residuals_complex_single_f(na, nev, na_rows, na_cols, as, z, ev, sc_desc, myid);
247 : #endif
248 :
249 192 : if (status !=0){
250 0 : printf("The computed EVs are not correct !\n");
251 : }
252 192 : if (status ==0){
253 192 : printf("All ok!\n");
254 : }
255 :
256 192 : free(sc_desc);
257 192 : free(a);
258 192 : free(z);
259 192 : free(as);
260 192 : free(ev);
261 :
262 : #ifdef WITH_MPI
263 128 : MPI_Finalize();
264 : #endif
265 192 : return 0;
266 : }
|